#include <G4PSFlatSurfaceCurrent.hh>
Inheritance diagram for G4PSFlatSurfaceCurrent:
Public Member Functions | |
G4PSFlatSurfaceCurrent (G4String name, G4int direction, G4int depth=0) | |
G4PSFlatSurfaceCurrent (G4String name, G4int direction, const G4String &unit, G4int depth=0) | |
virtual | ~G4PSFlatSurfaceCurrent () |
void | Weighted (G4bool flg=true) |
void | DivideByArea (G4bool flg=true) |
virtual void | Initialize (G4HCofThisEvent *) |
virtual void | EndOfEvent (G4HCofThisEvent *) |
virtual void | clear () |
virtual void | DrawAll () |
virtual void | PrintAll () |
virtual void | SetUnit (const G4String &unit) |
Protected Member Functions | |
virtual G4bool | ProcessHits (G4Step *, G4TouchableHistory *) |
G4int | IsSelectedSurface (G4Step *, G4Box *) |
virtual void | DefineUnitAndCategory () |
Definition at line 60 of file G4PSFlatSurfaceCurrent.hh.
Definition at line 59 of file G4PSFlatSurfaceCurrent.cc.
References DefineUnitAndCategory(), and SetUnit().
00061 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction), 00062 weighted(true),divideByArea(true) 00063 { 00064 DefineUnitAndCategory(); 00065 SetUnit("percm2"); 00066 }
G4PSFlatSurfaceCurrent::G4PSFlatSurfaceCurrent | ( | G4String | name, | |
G4int | direction, | |||
const G4String & | unit, | |||
G4int | depth = 0 | |||
) |
Definition at line 68 of file G4PSFlatSurfaceCurrent.cc.
References DefineUnitAndCategory(), and SetUnit().
00072 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction), 00073 weighted(true),divideByArea(true) 00074 { 00075 DefineUnitAndCategory(); 00076 SetUnit(unit); 00077 }
G4PSFlatSurfaceCurrent::~G4PSFlatSurfaceCurrent | ( | ) | [virtual] |
void G4PSFlatSurfaceCurrent::clear | ( | ) | [virtual] |
void G4PSFlatSurfaceCurrent::DefineUnitAndCategory | ( | ) | [protected, virtual] |
Definition at line 199 of file G4PSFlatSurfaceCurrent.cc.
Referenced by G4PSFlatSurfaceCurrent().
00199 { 00200 // Per Unit Surface 00201 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2)); 00202 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2)); 00203 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2)); 00204 }
void G4PSFlatSurfaceCurrent::DivideByArea | ( | G4bool | flg = true |
) | [inline] |
Definition at line 72 of file G4PSFlatSurfaceCurrent.hh.
Referenced by G4ScoreQuantityMessenger::SetNewValue().
void G4PSFlatSurfaceCurrent::DrawAll | ( | ) | [virtual] |
void G4PSFlatSurfaceCurrent::EndOfEvent | ( | G4HCofThisEvent * | ) | [virtual] |
void G4PSFlatSurfaceCurrent::Initialize | ( | G4HCofThisEvent * | ) | [virtual] |
Reimplemented from G4VPrimitiveScorer.
Definition at line 149 of file G4PSFlatSurfaceCurrent.cc.
References G4HCofThisEvent::AddHitsCollection(), G4VPrimitiveScorer::detector, G4VPrimitiveScorer::GetCollectionID(), G4VPrimitiveScorer::GetName(), and G4VSensitiveDetector::GetName().
00150 { 00151 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName()); 00152 if ( HCID < 0 ) HCID = GetCollectionID(0); 00153 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap); 00154 }
Definition at line 120 of file G4PSFlatSurfaceCurrent.cc.
References fCurrent_In, fCurrent_Out, fGeomBoundary, G4GeometryTolerance::GetInstance(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetStepStatus(), G4GeometryTolerance::GetSurfaceTolerance(), G4StepPoint::GetTouchableHandle(), and G4Box::GetZHalfLength().
Referenced by ProcessHits().
00120 { 00121 00122 G4TouchableHandle theTouchable = 00123 aStep->GetPreStepPoint()->GetTouchableHandle(); 00124 G4double kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance(); 00125 00126 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){ 00127 // Entering Geometry 00128 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition(); 00129 G4ThreeVector localpos1 = 00130 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1); 00131 if(std::fabs( localpos1.z() + boxSolid->GetZHalfLength())<kCarTolerance ){ 00132 return fCurrent_In; 00133 } 00134 } 00135 00136 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){ 00137 // Exiting Geometry 00138 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition(); 00139 G4ThreeVector localpos2 = 00140 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2); 00141 if(std::fabs( localpos2.z() + boxSolid->GetZHalfLength())<kCarTolerance ){ 00142 return fCurrent_Out; 00143 } 00144 } 00145 00146 return -1; 00147 }
void G4PSFlatSurfaceCurrent::PrintAll | ( | ) | [virtual] |
Reimplemented from G4VPrimitiveScorer.
Definition at line 166 of file G4PSFlatSurfaceCurrent.cc.
References G4VPrimitiveScorer::detector, G4cout, G4endl, G4VPrimitiveScorer::GetName(), G4VSensitiveDetector::GetName(), G4VPrimitiveScorer::GetUnit(), and G4VPrimitiveScorer::GetUnitValue().
00167 { 00168 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl; 00169 G4cout << " PrimitiveScorer " << GetName() <<G4endl; 00170 G4cout << " Number of entries " << EvtMap->entries() << G4endl; 00171 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin(); 00172 for(; itr != EvtMap->GetMap()->end(); itr++) { 00173 G4cout << " copy no.: " << itr->first << " current : " ; 00174 if ( divideByArea ) { 00175 G4cout << *(itr->second)/GetUnitValue() 00176 << " ["<<GetUnit()<<"]"; 00177 }else { 00178 G4cout << *(itr->second)/GetUnitValue() << " [tracks]"; 00179 } 00180 G4cout << G4endl; 00181 } 00182 }
G4bool G4PSFlatSurfaceCurrent::ProcessHits | ( | G4Step * | , | |
G4TouchableHistory * | ||||
) | [protected, virtual] |
Implements G4VPrimitiveScorer.
Definition at line 82 of file G4PSFlatSurfaceCurrent.cc.
References G4VSolid::ComputeDimensions(), G4VPVParameterisation::ComputeSolid(), fCurrent_InOut, G4VPrimitiveScorer::GetIndex(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetParameterisation(), G4StepPoint::GetPhysicalVolume(), G4Step::GetPreStepPoint(), G4LogicalVolume::GetSolid(), G4StepPoint::GetTouchable(), G4StepPoint::GetTouchableHandle(), G4StepPoint::GetWeight(), G4Box::GetXHalfLength(), G4Box::GetYHalfLength(), G4VPrimitiveScorer::indexDepth, IsSelectedSurface(), and TRUE.
00083 { 00084 G4StepPoint* preStep = aStep->GetPreStepPoint(); 00085 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume(); 00086 G4VPVParameterisation* physParam = physVol->GetParameterisation(); 00087 G4VSolid * solid = 0; 00088 if(physParam) 00089 { // for parameterized volume 00090 G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable())) 00091 ->GetReplicaNumber(indexDepth); 00092 solid = physParam->ComputeSolid(idx, physVol); 00093 solid->ComputeDimensions(physParam,idx,physVol); 00094 } 00095 else 00096 { // for ordinary volume 00097 solid = physVol->GetLogicalVolume()->GetSolid(); 00098 } 00099 00100 G4Box* boxSolid = (G4Box*)(solid); 00101 00102 G4int dirFlag =IsSelectedSurface(aStep,boxSolid); 00103 if ( dirFlag > 0 ) { 00104 if ( fDirection == fCurrent_InOut || fDirection == dirFlag ){ 00105 G4int index = GetIndex(aStep); 00106 G4TouchableHandle theTouchable = preStep->GetTouchableHandle(); 00107 G4double current = 1.0; 00108 if ( weighted ) current=preStep->GetWeight(); // Current (Particle Weight) 00109 if ( divideByArea ){ 00110 G4double square = 4.*boxSolid->GetXHalfLength()*boxSolid->GetYHalfLength(); 00111 current = current/square; // Normalized by Area 00112 } 00113 EvtMap->add(index,current); 00114 } 00115 } 00116 00117 return TRUE; 00118 }
void G4PSFlatSurfaceCurrent::SetUnit | ( | const G4String & | unit | ) | [virtual] |
Reimplemented from G4VPrimitiveScorer.
Definition at line 184 of file G4PSFlatSurfaceCurrent.cc.
References G4VPrimitiveScorer::CheckAndSetUnit(), G4Exception(), G4VPrimitiveScorer::GetName(), G4VPrimitiveScorer::GetUnit(), JustWarning, G4VPrimitiveScorer::unitName, and G4VPrimitiveScorer::unitValue.
Referenced by G4PSFlatSurfaceCurrent(), G4PSFlatSurfaceCurrent3D::G4PSFlatSurfaceCurrent3D(), and G4ScoreQuantityMessenger::SetNewValue().
00185 { 00186 if ( divideByArea ) { 00187 CheckAndSetUnit(unit,"Per Unit Surface"); 00188 } else { 00189 if (unit == "" ){ 00190 unitName = unit; 00191 unitValue = 1.0; 00192 }else{ 00193 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName(); 00194 G4Exception("G4PSFlatSurfaceCurrent::SetUnit","DetPS0007",JustWarning,msg); 00195 } 00196 } 00197 }
void G4PSFlatSurfaceCurrent::Weighted | ( | G4bool | flg = true |
) | [inline] |
Definition at line 70 of file G4PSFlatSurfaceCurrent.hh.
Referenced by G4ScoreQuantityMessenger::SetNewValue().