Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // 27 // 28 // G4PSFlatSurfaceFlux 29 #include "G4PSFlatSurfaceFlux.hh" 30 31 #include "G4SystemOfUnits.hh" 32 #include "G4StepStatus.hh" 33 #include "G4Track.hh" 34 #include "G4VSolid.hh" 35 #include "G4VPhysicalVolume.hh" 36 #include "G4VPVParameterisation.hh" 37 #include "G4UnitsTable.hh" 38 #include "G4GeometryTolerance.hh" 39 #include "G4VScoreHistFiller.hh" 40 41 //////////////////////////////////////////////////////////////////////////////// 42 // (Description) 43 // This is a primitive scorer class for scoring Surface Flux. 44 // Current version assumes only for G4Box shape, and the surface 45 // is fixed on -Z plane. 46 // 47 // Surface is defined at the -Z surface. 48 // Direction -Z +Z 49 // 0 IN || OUT ->|<- | 50 // 1 IN ->| | 51 // 2 OUT |<- | 52 // 53 // Created: 2005-11-14 Tsukasa ASO, Akinori Kimura. 54 // 55 // 18-Nov-2005 T.Aso, To use always positive value for anglefactor. 56 // 29-Mar-2007 T.Aso, Bug fix for momentum direction at outgoing flux. 57 // 2010-07-22 Introduce Unit specification. 58 // 2010-07-22 Add weighted and divideByAre options 59 // 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x) 60 // vs. cell flux * track weight (Makoto Asai) 61 /////////////////////////////////////////////////////////////////////////////// 62 63 G4PSFlatSurfaceFlux::G4PSFlatSurfaceFlux(const G4String& name, G4int direction, 64 G4int depth) 65 : G4PSFlatSurfaceFlux(name, direction, "percm2", depth) 66 {} 67 68 G4PSFlatSurfaceFlux::G4PSFlatSurfaceFlux(const G4String& name, G4int direction, 69 const G4String& unit, G4int depth) 70 : G4VPrimitivePlotter(name, depth) 71 , HCID(-1) 72 , fDirection(direction) 73 , EvtMap(nullptr) 74 , weighted(true) 75 , divideByArea(true) 76 { 77 DefineUnitAndCategory(); 78 SetUnit(unit); 79 } 80 81 G4bool G4PSFlatSurfaceFlux::ProcessHits(G4Step* aStep, G4TouchableHistory*) 82 { 83 G4StepPoint* preStep = aStep->GetPreStepPoint(); 84 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume(); 85 G4VPVParameterisation* physParam = physVol->GetParameterisation(); 86 G4VSolid* solid = nullptr; 87 if(physParam != nullptr) 88 { // for parameterized volume 89 G4int idx = 90 ((G4TouchableHistory*) (aStep->GetPreStepPoint()->GetTouchable())) 91 ->GetReplicaNumber(indexDepth); 92 solid = physParam->ComputeSolid(idx, physVol); 93 solid->ComputeDimensions(physParam, idx, physVol); 94 } 95 else 96 { // for ordinary volume 97 solid = physVol->GetLogicalVolume()->GetSolid(); 98 } 99 100 auto boxSolid = (G4Box*) (solid); 101 102 G4int dirFlag = IsSelectedSurface(aStep, boxSolid); 103 if(dirFlag > 0) 104 { 105 if(fDirection == fFlux_InOut || fDirection == dirFlag) 106 { 107 G4StepPoint* thisStep = nullptr; 108 if(dirFlag == fFlux_In) 109 { 110 thisStep = preStep; 111 } 112 else if(dirFlag == fFlux_Out) 113 { 114 thisStep = aStep->GetPostStepPoint(); 115 } 116 else 117 { 118 return false; 119 } 120 121 G4TouchableHandle theTouchable = thisStep->GetTouchableHandle(); 122 G4ThreeVector pdirection = thisStep->GetMomentumDirection(); 123 G4ThreeVector localdir = 124 theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection); 125 // 126 G4double angleFactor = localdir.z(); 127 if(angleFactor < 0) 128 angleFactor *= -1.; 129 G4double flux = 1.0; 130 if(weighted) 131 flux *= preStep->GetWeight(); // Current (Particle Weight) 132 // 133 G4double square = 134 4. * boxSolid->GetXHalfLength() * boxSolid->GetYHalfLength(); 135 // 136 flux = flux / angleFactor; // Flux with angle. 137 if(divideByArea) 138 flux /= square; 139 // 140 G4int index = GetIndex(aStep); 141 EvtMap->add(index, flux); 142 143 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend()) 144 { 145 auto filler = G4VScoreHistFiller::Instance(); 146 if(filler == nullptr) 147 { 148 G4Exception("G4PSFlatSurfaceFlux::ProcessHits", "SCORER0123", 149 JustWarning, 150 "G4TScoreHistFiller is not instantiated!! Histogram is " 151 "not filled."); 152 } 153 else 154 { 155 filler->FillH1(hitIDMap[index], preStep->GetKineticEnergy(), flux); 156 } 157 } 158 } 159 } 160 #ifdef debug 161 G4cout << " PASSED vol " << index << " trk " << trkid << " len " 162 << fFlatSurfaceFlux << G4endl; 163 #endif 164 165 return true; 166 } 167 168 G4int G4PSFlatSurfaceFlux::IsSelectedSurface(G4Step* aStep, G4Box* boxSolid) 169 { 170 G4TouchableHandle theTouchable = 171 aStep->GetPreStepPoint()->GetTouchableHandle(); 172 G4double kCarTolerance = 173 G4GeometryTolerance::GetInstance()->GetSurfaceTolerance(); 174 175 if(aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary) 176 { 177 // Entering Geometry 178 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition(); 179 G4ThreeVector localpos1 = 180 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1); 181 if(std::fabs(localpos1.z() + boxSolid->GetZHalfLength()) < kCarTolerance) 182 { 183 return fFlux_In; 184 } 185 } 186 187 if(aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary) 188 { 189 // Exiting Geometry 190 G4ThreeVector stppos2 = aStep->GetPostStepPoint()->GetPosition(); 191 G4ThreeVector localpos2 = 192 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2); 193 if(std::fabs(localpos2.z() + boxSolid->GetZHalfLength()) < kCarTolerance) 194 { 195 return fFlux_Out; 196 } 197 } 198 199 return -1; 200 } 201 202 void G4PSFlatSurfaceFlux::Initialize(G4HCofThisEvent* HCE) 203 { 204 EvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(), 205 GetName()); 206 if(HCID < 0) 207 HCID = GetCollectionID(0); 208 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap); 209 } 210 211 void G4PSFlatSurfaceFlux::clear() { EvtMap->clear(); } 212 213 void G4PSFlatSurfaceFlux::PrintAll() 214 { 215 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl; 216 G4cout << " PrimitiveScorer" << GetName() << G4endl; 217 G4cout << " Number of entries " << EvtMap->entries() << G4endl; 218 for(const auto& [copy, flux] : *(EvtMap->GetMap())) 219 { 220 G4cout << " copy no.: " << copy 221 << " flux : " << *(flux) / GetUnitValue() << " [" 222 << GetUnit() << "]" << G4endl; 223 } 224 } 225 226 void G4PSFlatSurfaceFlux::SetUnit(const G4String& unit) 227 { 228 if(divideByArea) 229 { 230 CheckAndSetUnit(unit, "Per Unit Surface"); 231 } 232 else 233 { 234 if(unit.empty()) 235 { 236 unitName = unit; 237 unitValue = 1.0; 238 } 239 else 240 { 241 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" + 242 GetUnit() + "] ) for " + GetName(); 243 G4Exception("G4PSFlatSurfaceFlux::SetUnit", "DetPS0008", JustWarning, 244 msg); 245 } 246 } 247 } 248 249 void G4PSFlatSurfaceFlux::DefineUnitAndCategory() 250 { 251 // Per Unit Surface 252 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface", 253 (1. / cm2)); 254 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface", 255 (1. / mm2)); 256 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2)); 257 } 258