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 // G4PSFlatSurfaceCurrent 29 #include "G4PSFlatSurfaceCurrent.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 only Surface Flux. 44 // Current version assumes only for G4Box shape. 45 // 46 // Surface is defined at the -Z surface. 47 // Direction -Z +Z 48 // 0 IN || OUT ->|<- | 49 // 1 IN ->| | 50 // 2 OUT |<- | 51 // 52 // Created: 2005-11-14 Tsukasa ASO, Akinori Kimura. 53 // 17-Nov-2005 T.Aso, Bug fix for area definition. 54 // 31-Mar-2007 T.Aso, Add option for normalizing by the area. 55 // 2010-07-22 Introduce Unit specification. 56 // 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x) 57 // vs. Surface Current * track weight (y) (Makoto Asai) 58 // 59 /////////////////////////////////////////////////////////////////////////////// 60 61 G4PSFlatSurfaceCurrent::G4PSFlatSurfaceCurrent(const G4String& name, G4int direction, 62 G4int depth) 63 : G4PSFlatSurfaceCurrent(name, direction, "percm2", depth) 64 {} 65 66 G4PSFlatSurfaceCurrent::G4PSFlatSurfaceCurrent(const G4String& name, G4int direction, 67 const G4String& unit, 68 G4int depth) 69 : G4VPrimitivePlotter(name, depth) 70 , HCID(-1) 71 , fDirection(direction) 72 , EvtMap(nullptr) 73 , weighted(true) 74 , divideByArea(true) 75 { 76 DefineUnitAndCategory(); 77 SetUnit(unit); 78 } 79 80 G4bool G4PSFlatSurfaceCurrent::ProcessHits(G4Step* aStep, G4TouchableHistory*) 81 { 82 G4StepPoint* preStep = aStep->GetPreStepPoint(); 83 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume(); 84 G4VPVParameterisation* physParam = physVol->GetParameterisation(); 85 G4VSolid* solid = nullptr; 86 if(physParam != nullptr) 87 { // for parameterized volume 88 G4int idx = 89 ((G4TouchableHistory*) (aStep->GetPreStepPoint()->GetTouchable())) 90 ->GetReplicaNumber(indexDepth); 91 solid = physParam->ComputeSolid(idx, physVol); 92 solid->ComputeDimensions(physParam, idx, physVol); 93 } 94 else 95 { // for ordinary volume 96 solid = physVol->GetLogicalVolume()->GetSolid(); 97 } 98 99 auto boxSolid = (G4Box*) (solid); 100 101 G4int dirFlag = IsSelectedSurface(aStep, boxSolid); 102 if(dirFlag > 0) 103 { 104 if(fDirection == fCurrent_InOut || fDirection == dirFlag) 105 { 106 G4int index = GetIndex(aStep); 107 G4TouchableHandle theTouchable = preStep->GetTouchableHandle(); 108 G4double current = 1.0; 109 if(weighted) 110 current = preStep->GetWeight(); // Current (Particle Weight) 111 if(divideByArea) 112 { 113 G4double square = 114 4. * boxSolid->GetXHalfLength() * boxSolid->GetYHalfLength(); 115 current = current / square; // Normalized by Area 116 } 117 EvtMap->add(index, current); 118 119 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend()) 120 { 121 auto filler = G4VScoreHistFiller::Instance(); 122 if(filler == nullptr) 123 { 124 G4Exception("G4PSFlatSurfaceCurrent::ProcessHits", "SCORER0123", 125 JustWarning, 126 "G4TScoreHistFiller is not instantiated!! Histogram is " 127 "not filled."); 128 } 129 else 130 { 131 filler->FillH1(hitIDMap[index], preStep->GetKineticEnergy(), current); 132 } 133 } 134 } 135 } 136 137 return true; 138 } 139 140 G4int G4PSFlatSurfaceCurrent::IsSelectedSurface(G4Step* aStep, G4Box* boxSolid) 141 { 142 G4TouchableHandle theTouchable = 143 aStep->GetPreStepPoint()->GetTouchableHandle(); 144 G4double kCarTolerance = 145 G4GeometryTolerance::GetInstance()->GetSurfaceTolerance(); 146 147 if(aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary) 148 { 149 // Entering Geometry 150 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition(); 151 G4ThreeVector localpos1 = 152 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1); 153 if(std::fabs(localpos1.z() + boxSolid->GetZHalfLength()) < kCarTolerance) 154 { 155 return fCurrent_In; 156 } 157 } 158 159 if(aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary) 160 { 161 // Exiting Geometry 162 G4ThreeVector stppos2 = aStep->GetPostStepPoint()->GetPosition(); 163 G4ThreeVector localpos2 = 164 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2); 165 if(std::fabs(localpos2.z() + boxSolid->GetZHalfLength()) < kCarTolerance) 166 { 167 return fCurrent_Out; 168 } 169 } 170 171 return -1; 172 } 173 174 void G4PSFlatSurfaceCurrent::Initialize(G4HCofThisEvent* HCE) 175 { 176 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName()); 177 if(HCID < 0) 178 HCID = GetCollectionID(0); 179 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap); 180 } 181 182 void G4PSFlatSurfaceCurrent::clear() { EvtMap->clear(); } 183 184 void G4PSFlatSurfaceCurrent::PrintAll() 185 { 186 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl; 187 G4cout << " PrimitiveScorer " << GetName() << G4endl; 188 G4cout << " Number of entries " << EvtMap->entries() << G4endl; 189 for(const auto& [copy, current] : *(EvtMap->GetMap())) 190 { 191 G4cout << " copy no.: " << copy << " current : "; 192 if(divideByArea) 193 { 194 G4cout << *(current) / GetUnitValue() << " [" << GetUnit() << "]"; 195 } 196 else 197 { 198 G4cout << *(current) / GetUnitValue() << " [tracks]"; 199 } 200 G4cout << G4endl; 201 } 202 } 203 204 void G4PSFlatSurfaceCurrent::SetUnit(const G4String& unit) 205 { 206 if(divideByArea) 207 { 208 CheckAndSetUnit(unit, "Per Unit Surface"); 209 } 210 else 211 { 212 if(unit.empty()) 213 { 214 unitName = unit; 215 unitValue = 1.0; 216 } 217 else 218 { 219 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" + 220 GetUnit() + "] ) for " + GetName(); 221 G4Exception("G4PSFlatSurfaceCurrent::SetUnit", "DetPS0007", JustWarning, 222 msg); 223 } 224 } 225 } 226 227 void G4PSFlatSurfaceCurrent::DefineUnitAndCategory() 228 { 229 // Per Unit Surface 230 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface", 231 (1. / cm2)); 232 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface", 233 (1. / mm2)); 234 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2)); 235 } 236