Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4PSCylinderSurfaceCurrent Class Reference

#include <G4PSCylinderSurfaceCurrent.hh>

Inheritance diagram for G4PSCylinderSurfaceCurrent:
G4VPrimitivePlotter G4VPrimitiveScorer G4PSCylinderSurfaceCurrent3D

Public Member Functions

virtual void clear ()
 
void DivideByArea (G4bool flg=true)
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
 G4PSCylinderSurfaceCurrent (G4String name, G4int direction, const G4String &unit, G4int depth=0)
 
 G4PSCylinderSurfaceCurrent (G4String name, G4int direction, G4int depth=0)
 
G4int GetCollectionID (G4int)
 
G4VSDFilterGetFilter () const
 
G4MultiFunctionalDetectorGetMultiFunctionalDetector () const
 
G4String GetName () const
 
G4int GetNumberOfHist () const
 
const G4StringGetUnit () const
 
G4double GetUnitValue () const
 
G4int GetVerboseLevel () const
 
virtual void Initialize (G4HCofThisEvent *)
 
void Plot (G4int copyNo, G4int histID)
 
virtual void PrintAll ()
 
void SetFilter (G4VSDFilter *f)
 
void SetMultiFunctionalDetector (G4MultiFunctionalDetector *d)
 
void SetNijk (G4int i, G4int j, G4int k)
 
virtual void SetUnit (const G4String &unit)
 
void SetVerboseLevel (G4int vl)
 
void Weighted (G4bool flg=true)
 
virtual ~G4PSCylinderSurfaceCurrent ()
 

Protected Member Functions

void CheckAndSetUnit (const G4String &unit, const G4String &category)
 
G4VSolidComputeCurrentSolid (G4Step *aStep)
 
G4VSolidComputeSolid (G4Step *aStep, G4int replicaIdx)
 
virtual void DefineUnitAndCategory ()
 
virtual G4int GetIndex (G4Step *)
 
G4int IsSelectedSurface (G4Step *, G4Tubs *)
 
virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
 

Protected Attributes

G4MultiFunctionalDetectordetector
 
G4VSDFilterfilter
 
G4int fNi
 
G4int fNj
 
G4int fNk
 
std::map< G4int, G4inthitIDMap
 
G4int indexDepth
 
G4String primitiveName
 
G4String unitName
 
G4double unitValue
 
G4int verboseLevel
 

Private Member Functions

G4bool HitPrimitive (G4Step *aStep, G4TouchableHistory *ROhis)
 

Private Attributes

G4bool divideByArea
 
G4THitsMap< G4double > * EvtMap
 
G4int fDirection
 
G4int HCID
 
G4bool weighted
 

Detailed Description

Definition at line 60 of file G4PSCylinderSurfaceCurrent.hh.

Constructor & Destructor Documentation

◆ G4PSCylinderSurfaceCurrent() [1/2]

G4PSCylinderSurfaceCurrent::G4PSCylinderSurfaceCurrent ( G4String  name,
G4int  direction,
G4int  depth = 0 
)

Definition at line 59 of file G4PSCylinderSurfaceCurrent.cc.

63 , HCID(-1)
64 , fDirection(direction)
65 , EvtMap(0)
66 , weighted(true)
67 , divideByArea(true)
68{
70 SetUnit("percm2");
71}
virtual void SetUnit(const G4String &unit)
G4VPrimitivePlotter(G4String name, G4int depth=0)
const char * name(G4int ptype)

References DefineUnitAndCategory(), and SetUnit().

◆ G4PSCylinderSurfaceCurrent() [2/2]

G4PSCylinderSurfaceCurrent::G4PSCylinderSurfaceCurrent ( G4String  name,
G4int  direction,
const G4String unit,
G4int  depth = 0 
)

Definition at line 73 of file G4PSCylinderSurfaceCurrent.cc.

78 , HCID(-1)
79 , fDirection(direction)
80 , EvtMap(0)
81 , weighted(true)
82 , divideByArea(true)
83{
85 SetUnit(unit);
86}

References DefineUnitAndCategory(), and SetUnit().

◆ ~G4PSCylinderSurfaceCurrent()

G4PSCylinderSurfaceCurrent::~G4PSCylinderSurfaceCurrent ( )
virtual

Definition at line 88 of file G4PSCylinderSurfaceCurrent.cc.

88{ ; }

Member Function Documentation

◆ CheckAndSetUnit()

void G4VPrimitiveScorer::CheckAndSetUnit ( const G4String unit,
const G4String category 
)
protectedinherited

Definition at line 81 of file G4VPrimitiveScorer.cc.

83{
84 if(G4UnitDefinition::GetCategory(unit) == category)
85 {
86 unitName = unit;
88 }
89 else
90 {
91 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
92 GetUnit() + "] ) requested for " + GetName();
93 G4Exception("G4VPrimitiveScorer::CheckAndSetUnit", "Det0151", JustWarning,
94 msg);
95 }
96}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)
G4String GetName() const
const G4String & GetUnit() const

References G4Exception(), G4UnitDefinition::GetCategory(), G4VPrimitiveScorer::GetName(), G4VPrimitiveScorer::GetUnit(), G4UnitDefinition::GetValueOf(), JustWarning, G4VPrimitiveScorer::unitName, and G4VPrimitiveScorer::unitValue.

Referenced by G4PSCellCharge::SetUnit(), G4PSCellFlux::SetUnit(), SetUnit(), G4PSCylinderSurfaceFlux::SetUnit(), G4PSDoseDeposit::SetUnit(), G4PSEnergyDeposit::SetUnit(), G4PSFlatSurfaceCurrent::SetUnit(), G4PSFlatSurfaceFlux::SetUnit(), G4PSMinKinEAtGeneration::SetUnit(), G4PSPassageCellFlux::SetUnit(), G4PSPassageTrackLength::SetUnit(), G4PSSphereSurfaceCurrent::SetUnit(), G4PSSphereSurfaceFlux::SetUnit(), and G4PSTrackLength::SetUnit().

◆ clear()

void G4PSCylinderSurfaceCurrent::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 202 of file G4PSCylinderSurfaceCurrent.cc.

202{ EvtMap->clear(); }
void clear()
Definition: G4THitsMap.hh:524

References G4VTHitsMap< T, Map_t >::clear(), and EvtMap.

◆ ComputeCurrentSolid()

G4VSolid * G4VPrimitiveScorer::ComputeCurrentSolid ( G4Step aStep)
protectedinherited

Definition at line 128 of file G4VPrimitiveScorer.cc.

129{
130 G4StepPoint* preStep = aStep->GetPreStepPoint();
131 // The only difference: did not know the replica number
132 G4int replicaIdx =
133 (static_cast<const G4TouchableHistory*>(preStep->GetTouchable()))
134 ->GetReplicaNumber(indexDepth);
135
136 return ComputeSolid(aStep, replicaIdx);
137}
int G4int
Definition: G4Types.hh:85
const G4VTouchable * GetTouchable() const
G4StepPoint * GetPreStepPoint() const
G4VSolid * ComputeSolid(G4Step *aStep, G4int replicaIdx)

References G4VPrimitiveScorer::ComputeSolid(), G4Step::GetPreStepPoint(), G4StepPoint::GetTouchable(), and G4VPrimitiveScorer::indexDepth.

Referenced by ProcessHits(), G4PSCylinderSurfaceFlux::ProcessHits(), and G4PSSphereSurfaceCurrent::ProcessHits().

◆ ComputeSolid()

G4VSolid * G4VPrimitiveScorer::ComputeSolid ( G4Step aStep,
G4int  replicaIdx 
)
protectedinherited

Definition at line 98 of file G4VPrimitiveScorer.cc.

99{
100 G4VSolid* solid = nullptr;
101 G4StepPoint* preStep = aStep->GetPreStepPoint();
102
103 auto physVol = preStep->GetPhysicalVolume();
104 G4VPVParameterisation* physParam = physVol->GetParameterisation();
105
106 if(physParam)
107 { // for parameterized volume
108 if(replicaIdx < 0)
109 {
111 desc << "Incorrect replica number --- GetReplicaNumber : " << replicaIdx
112 << G4endl;
113 G4Exception("G4VPrimitiveScorer::ComputeSolid", "DetPS0001", JustWarning,
114 desc);
115 // replicaIdx= 0; // You must ensure that it's in range !!!
116 }
117 solid = physParam->ComputeSolid(replicaIdx, physVol);
118 solid->ComputeDimensions(physParam, replicaIdx, physVol);
119 }
120 else
121 { // for ordinary volume
122 solid = physVol->GetLogicalVolume()->GetSolid();
123 }
124
125 return solid;
126}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
G4VPhysicalVolume * GetPhysicalVolume() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137

References G4VSolid::ComputeDimensions(), G4VPVParameterisation::ComputeSolid(), G4endl, G4Exception(), G4StepPoint::GetPhysicalVolume(), G4Step::GetPreStepPoint(), and JustWarning.

Referenced by G4VPrimitiveScorer::ComputeCurrentSolid(), G4PSCellFlux::ComputeVolume(), G4PSDoseDeposit::ComputeVolume(), and G4PSPassageCellFlux::ComputeVolume().

◆ DefineUnitAndCategory()

void G4PSCylinderSurfaceCurrent::DefineUnitAndCategory ( )
protectedvirtual

Definition at line 250 of file G4PSCylinderSurfaceCurrent.cc.

251{
252 // Per Unit Surface
253 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
254 (1. / cm2));
255 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
256 (1. / mm2));
257 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
258}
static constexpr double mm2
Definition: G4SIunits.hh:96
static constexpr double cm2
Definition: G4SIunits.hh:100
static constexpr double m2
Definition: G4SIunits.hh:110

References cm2, m2, and mm2.

Referenced by G4PSCylinderSurfaceCurrent().

◆ DivideByArea()

void G4PSCylinderSurfaceCurrent::DivideByArea ( G4bool  flg = true)
inline

Definition at line 71 of file G4PSCylinderSurfaceCurrent.hh.

71{ divideByArea = flg; }

References divideByArea.

◆ DrawAll()

void G4PSCylinderSurfaceCurrent::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 204 of file G4PSCylinderSurfaceCurrent.cc.

204{ ; }

◆ EndOfEvent()

void G4PSCylinderSurfaceCurrent::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 200 of file G4PSCylinderSurfaceCurrent.cc.

200{ ; }

◆ GetCollectionID()

G4int G4VPrimitiveScorer::GetCollectionID ( G4int  )
inherited

Definition at line 55 of file G4VPrimitiveScorer.cc.

56{
57 if(detector)
59 "/" + primitiveName);
60 else
61 return -1;
62}
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:38
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:142
G4MultiFunctionalDetector * detector

References G4VPrimitiveScorer::detector, G4SDManager::GetCollectionID(), G4VSensitiveDetector::GetName(), G4SDManager::GetSDMpointer(), and G4VPrimitiveScorer::primitiveName.

Referenced by G4PSCellCharge::Initialize(), G4PSCellFlux::Initialize(), Initialize(), G4PSCylinderSurfaceFlux::Initialize(), G4PSDoseDeposit::Initialize(), G4PSEnergyDeposit::Initialize(), G4PSFlatSurfaceCurrent::Initialize(), G4PSFlatSurfaceFlux::Initialize(), G4PSMinKinEAtGeneration::Initialize(), G4PSNofCollision::Initialize(), G4PSNofSecondary::Initialize(), G4PSNofStep::Initialize(), G4PSPassageCellCurrent::Initialize(), G4PSPassageCellFlux::Initialize(), G4PSPassageTrackLength::Initialize(), G4PSPopulation::Initialize(), G4PSSphereSurfaceCurrent::Initialize(), G4PSSphereSurfaceFlux::Initialize(), G4PSTermination::Initialize(), G4PSTrackCounter::Initialize(), G4PSTrackLength::Initialize(), and G4PSVolumeFlux::Initialize().

◆ GetFilter()

G4VSDFilter * G4VPrimitiveScorer::GetFilter ( ) const
inlineinherited

Definition at line 108 of file G4VPrimitiveScorer.hh.

108{ return filter; }

References G4VPrimitiveScorer::filter.

Referenced by G4VScoringMesh::List(), and G4VScoringMesh::SetFilter().

◆ GetIndex()

G4int G4VPrimitiveScorer::GetIndex ( G4Step aStep)
protectedvirtualinherited

Reimplemented in G4PSCellCharge3D, G4PSCellFlux3D, G4PSCylinderSurfaceCurrent3D, G4PSCylinderSurfaceFlux3D, G4PSDoseDeposit3D, G4PSEnergyDeposit3D, G4PSFlatSurfaceCurrent3D, G4PSFlatSurfaceFlux3D, G4PSMinKinEAtGeneration3D, G4PSNofCollision3D, G4PSNofSecondary3D, G4PSNofStep3D, G4PSPassageCellCurrent3D, G4PSPassageCellFlux3D, G4PSPassageTrackLength3D, G4PSPopulation3D, G4PSSphereSurfaceCurrent3D, G4PSSphereSurfaceFlux3D, G4PSStepChecker3D, G4PSTermination3D, G4PSTrackCounter3D, G4PSTrackLength3D, and G4PSVolumeFlux3D.

Definition at line 74 of file G4VPrimitiveScorer.cc.

75{
76 G4StepPoint* preStep = aStep->GetPreStepPoint();
78 return th->GetReplicaNumber(indexDepth);
79}
G4int GetReplicaNumber(G4int depth=0) const

References G4Step::GetPreStepPoint(), G4TouchableHistory::GetReplicaNumber(), G4StepPoint::GetTouchable(), and G4VPrimitiveScorer::indexDepth.

Referenced by G4PSCellCharge::ProcessHits(), G4PSCellFlux::ProcessHits(), ProcessHits(), G4PSCylinderSurfaceFlux::ProcessHits(), G4PSDoseDeposit::ProcessHits(), G4PSEnergyDeposit::ProcessHits(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSMinKinEAtGeneration::ProcessHits(), G4PSNofCollision::ProcessHits(), G4PSNofSecondary::ProcessHits(), G4PSNofStep::ProcessHits(), G4PSPassageCellCurrent::ProcessHits(), G4PSPassageCellFlux::ProcessHits(), G4PSPassageTrackLength::ProcessHits(), G4PSPopulation::ProcessHits(), G4PSSphereSurfaceCurrent::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSStepChecker::ProcessHits(), G4PSTermination::ProcessHits(), G4PSTrackCounter::ProcessHits(), G4PSTrackLength::ProcessHits(), and G4PSVolumeFlux::ProcessHits().

◆ GetMultiFunctionalDetector()

G4MultiFunctionalDetector * G4VPrimitiveScorer::GetMultiFunctionalDetector ( ) const
inlineinherited

◆ GetName()

G4String G4VPrimitiveScorer::GetName ( ) const
inlineinherited

Definition at line 106 of file G4VPrimitiveScorer.hh.

106{ return primitiveName; }

References G4VPrimitiveScorer::primitiveName.

Referenced by G4VPrimitiveScorer::CheckAndSetUnit(), G4VScoringMesh::GetPrimitiveScorer(), G4PSCellCharge::Initialize(), G4PSCellFlux::Initialize(), Initialize(), G4PSCylinderSurfaceFlux::Initialize(), G4PSDoseDeposit::Initialize(), G4PSEnergyDeposit::Initialize(), G4PSFlatSurfaceCurrent::Initialize(), G4PSFlatSurfaceFlux::Initialize(), G4PSMinKinEAtGeneration::Initialize(), G4PSNofCollision::Initialize(), G4PSNofSecondary::Initialize(), G4PSNofStep::Initialize(), G4PSPassageCellCurrent::Initialize(), G4PSPassageCellFlux::Initialize(), G4PSPassageTrackLength::Initialize(), G4PSPopulation::Initialize(), G4PSSphereSurfaceCurrent::Initialize(), G4PSSphereSurfaceFlux::Initialize(), G4PSTermination::Initialize(), G4PSTrackCounter::Initialize(), G4PSTrackLength::Initialize(), G4PSVolumeFlux::Initialize(), G4VScoringMesh::List(), G4PSCellCharge::PrintAll(), G4PSCellFlux::PrintAll(), PrintAll(), G4PSCylinderSurfaceFlux::PrintAll(), G4PSDoseDeposit::PrintAll(), G4PSEnergyDeposit::PrintAll(), G4PSFlatSurfaceCurrent::PrintAll(), G4PSFlatSurfaceFlux::PrintAll(), G4PSMinKinEAtGeneration::PrintAll(), G4PSNofCollision::PrintAll(), G4PSNofSecondary::PrintAll(), G4PSNofStep::PrintAll(), G4PSPassageCellCurrent::PrintAll(), G4PSPassageCellFlux::PrintAll(), G4PSPassageTrackLength::PrintAll(), G4PSPopulation::PrintAll(), G4PSSphereSurfaceCurrent::PrintAll(), G4PSSphereSurfaceFlux::PrintAll(), G4PSTermination::PrintAll(), G4PSTrackCounter::PrintAll(), G4PSTrackLength::PrintAll(), G4PSVolumeFlux::PrintAll(), G4MultiFunctionalDetector::RegisterPrimitive(), G4MultiFunctionalDetector::RemovePrimitive(), G4VScoringMesh::SetFilter(), G4VScoringMesh::SetPrimitiveScorer(), SetUnit(), G4PSCylinderSurfaceFlux::SetUnit(), G4PSFlatSurfaceCurrent::SetUnit(), G4PSFlatSurfaceFlux::SetUnit(), G4PSNofCollision::SetUnit(), G4PSNofSecondary::SetUnit(), G4PSNofStep::SetUnit(), G4PSPassageCellCurrent::SetUnit(), G4PSPopulation::SetUnit(), G4PSSphereSurfaceCurrent::SetUnit(), G4PSSphereSurfaceFlux::SetUnit(), G4PSTermination::SetUnit(), and G4PSTrackCounter::SetUnit().

◆ GetNumberOfHist()

G4int G4VPrimitivePlotter::GetNumberOfHist ( ) const
inlineinherited

Definition at line 51 of file G4VPrimitivePlotter.hh.

51{ return hitIDMap.size(); }
std::map< G4int, G4int > hitIDMap

References G4VPrimitivePlotter::hitIDMap.

◆ GetUnit()

const G4String & G4VPrimitiveScorer::GetUnit ( ) const
inlineinherited

◆ GetUnitValue()

G4double G4VPrimitiveScorer::GetUnitValue ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4VPrimitiveScorer::GetVerboseLevel ( ) const
inlineinherited

Definition at line 110 of file G4VPrimitiveScorer.hh.

110{ return verboseLevel; }

References G4VPrimitiveScorer::verboseLevel.

◆ HitPrimitive()

G4bool G4VPrimitiveScorer::HitPrimitive ( G4Step aStep,
G4TouchableHistory ROhis 
)
inlineprivateinherited

Definition at line 113 of file G4VPrimitiveScorer.hh.

114 {
115 if(filter)
116 {
117 if(!(filter->Accept(aStep)))
118 return false;
119 }
120 return ProcessHits(aStep, ROhis);
121 }
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)=0
virtual G4bool Accept(const G4Step *) const =0

References G4VSDFilter::Accept(), G4VPrimitiveScorer::filter, and G4VPrimitiveScorer::ProcessHits().

◆ Initialize()

void G4PSCylinderSurfaceCurrent::Initialize ( G4HCofThisEvent HCE)
virtual

◆ IsSelectedSurface()

G4int G4PSCylinderSurfaceCurrent::IsSelectedSurface ( G4Step aStep,
G4Tubs tubsSolid 
)
protected

Definition at line 141 of file G4PSCylinderSurfaceCurrent.cc.

143{
144 G4TouchableHandle theTouchable =
148
150 {
151 // Entering Geometry
152 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition();
153 G4ThreeVector localpos1 =
154 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
155 if(std::fabs(localpos1.z()) > tubsSolid->GetZHalfLength())
156 return -1;
157 G4double localR2 =
158 localpos1.x() * localpos1.x() + localpos1.y() * localpos1.y();
159 G4double InsideRadius = tubsSolid->GetInnerRadius();
160 if(localR2 >
161 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
162 localR2 <
163 (InsideRadius + kCarTolerance) * (InsideRadius + kCarTolerance))
164 {
165 return fCurrent_In;
166 }
167 }
168
170 {
171 // Exiting Geometry
172 G4ThreeVector stppos2 = aStep->GetPostStepPoint()->GetPosition();
173 G4ThreeVector localpos2 =
174 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
175 if(std::fabs(localpos2.z()) > tubsSolid->GetZHalfLength())
176 return -1;
177 G4double localR2 =
178 localpos2.x() * localpos2.x() + localpos2.y() * localpos2.y();
179 G4double InsideRadius = tubsSolid->GetInnerRadius();
180 if(localR2 >
181 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
182 localR2 <
183 (InsideRadius + kCarTolerance) * (InsideRadius + kCarTolerance))
184 {
185 return fCurrent_Out;
186 }
187 }
188
189 return -1;
190}
const G4double kCarTolerance
@ fCurrent_Out
@ fCurrent_In
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
double z() const
double x() const
double y() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
const G4AffineTransform & GetTopTransform() const
G4StepStatus GetStepStatus() const
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4StepPoint * GetPostStepPoint() const
G4double GetZHalfLength() const
G4double GetInnerRadius() const
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:82

References fCurrent_In, fCurrent_Out, fGeomBoundary, G4VTouchable::GetHistory(), G4Tubs::GetInnerRadius(), G4GeometryTolerance::GetInstance(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetStepStatus(), G4GeometryTolerance::GetSurfaceTolerance(), G4NavigationHistory::GetTopTransform(), G4StepPoint::GetTouchableHandle(), G4Tubs::GetZHalfLength(), kCarTolerance, G4AffineTransform::TransformPoint(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by ProcessHits().

◆ Plot()

void G4VPrimitivePlotter::Plot ( G4int  copyNo,
G4int  histID 
)
inlineinherited

Definition at line 45 of file G4VPrimitivePlotter.hh.

45{ hitIDMap[copyNo] = histID; }

References G4VPrimitivePlotter::hitIDMap.

◆ PrintAll()

void G4PSCylinderSurfaceCurrent::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 206 of file G4PSCylinderSurfaceCurrent.cc.

207{
208 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
209 G4cout << " PrimitiveScorer " << GetName() << G4endl;
210 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
211 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
212 for(; itr != EvtMap->GetMap()->end(); itr++)
213 {
214 G4cout << " copy no.: " << itr->first << " current : ";
215 if(divideByArea)
216 {
217 G4cout << *(itr->second) / GetUnitValue() << " [" << GetUnit() << "]";
218 }
219 else
220 {
221 G4cout << *(itr->second) << " [tracks]";
222 }
223 G4cout << G4endl;
224 }
225}
G4GLOB_DLL std::ostream G4cout
G4double GetUnitValue() const
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t entries() const
Definition: G4THitsMap.hh:158

References G4VPrimitiveScorer::detector, divideByArea, G4VTHitsMap< T, Map_t >::entries(), EvtMap, G4cout, G4endl, G4VTHitsMap< T, Map_t >::GetMap(), G4VPrimitiveScorer::GetName(), G4VSensitiveDetector::GetName(), G4VPrimitiveScorer::GetUnit(), and G4VPrimitiveScorer::GetUnitValue().

◆ ProcessHits()

G4bool G4PSCylinderSurfaceCurrent::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 90 of file G4PSCylinderSurfaceCurrent.cc.

92{
93 G4StepPoint* preStep = aStep->GetPreStepPoint();
94 G4VSolid* solid = ComputeCurrentSolid(aStep);
95 G4Tubs* tubsSolid = static_cast<G4Tubs*>(solid);
96
97 G4int dirFlag = IsSelectedSurface(aStep, tubsSolid);
98 // G4cout << " pos " << preStep->GetPosition() <<" dirFlag " << G4endl;
99 if(dirFlag > 0)
100 {
101 if(fDirection == fCurrent_InOut || fDirection == dirFlag)
102 {
103 G4TouchableHandle theTouchable = preStep->GetTouchableHandle();
104 //
105 G4double current = 1.0;
106 if(weighted)
107 current = preStep->GetWeight(); // Current (Particle Weight)
108 //
109 if(divideByArea)
110 {
111 G4double square = 2. * tubsSolid->GetZHalfLength() *
112 tubsSolid->GetInnerRadius() *
113 tubsSolid->GetDeltaPhiAngle() / radian;
114 current = current / square; // Current normalized by Area
115 }
116
117 G4int index = GetIndex(aStep);
118 EvtMap->add(index, current);
119
120 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
121 {
122 auto filler = G4VScoreHistFiller::Instance();
123 if(!filler)
124 {
125 G4Exception("G4PSCylinderSurfaceCurrent::ProcessHits", "SCORER0123",
127 "G4TScoreHistFiller is not instantiated!! Histogram is "
128 "not filled.");
129 }
130 else
131 {
132 filler->FillH1(hitIDMap[index], preStep->GetKineticEnergy(), current);
133 }
134 }
135 }
136 }
137
138 return TRUE;
139}
@ fCurrent_InOut
static constexpr double radian
Definition: G4SIunits.hh:122
#define TRUE
Definition: Globals.hh:27
G4int IsSelectedSurface(G4Step *, G4Tubs *)
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Tubs.hh:75
G4double GetDeltaPhiAngle() const
virtual G4int GetIndex(G4Step *)
G4VSolid * ComputeCurrentSolid(G4Step *aStep)
static G4VScoreHistFiller * Instance()
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

References G4VTHitsMap< T, Map_t >::add(), G4VPrimitiveScorer::ComputeCurrentSolid(), divideByArea, EvtMap, fCurrent_InOut, fDirection, G4Exception(), G4Tubs::GetDeltaPhiAngle(), G4VPrimitiveScorer::GetIndex(), G4Tubs::GetInnerRadius(), G4StepPoint::GetKineticEnergy(), G4Step::GetPreStepPoint(), G4StepPoint::GetTouchableHandle(), G4StepPoint::GetWeight(), G4Tubs::GetZHalfLength(), G4VPrimitivePlotter::hitIDMap, G4VScoreHistFiller::Instance(), IsSelectedSurface(), JustWarning, radian, TRUE, and weighted.

◆ SetFilter()

void G4VPrimitiveScorer::SetFilter ( G4VSDFilter f)
inlineinherited

Definition at line 107 of file G4VPrimitiveScorer.hh.

107{ filter = f; }

References G4VPrimitiveScorer::filter.

Referenced by G4VScoringMesh::SetFilter().

◆ SetMultiFunctionalDetector()

void G4VPrimitiveScorer::SetMultiFunctionalDetector ( G4MultiFunctionalDetector d)
inlineinherited

◆ SetNijk()

void G4VPrimitiveScorer::SetNijk ( G4int  i,
G4int  j,
G4int  k 
)
inlineinherited

◆ SetUnit()

void G4PSCylinderSurfaceCurrent::SetUnit ( const G4String unit)
virtual

Definition at line 227 of file G4PSCylinderSurfaceCurrent.cc.

228{
229 if(divideByArea)
230 {
231 CheckAndSetUnit(unit, "Per Unit Surface");
232 }
233 else
234 {
235 if(unit == "")
236 {
237 unitName = unit;
238 unitValue = 1.0;
239 }
240 else
241 {
242 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
243 GetUnit() + "] ) for " + GetName();
244 G4Exception("G4PSCylinderSurfaceCurrent::SetUnit", "DetPS0002",
245 JustWarning, msg);
246 }
247 }
248}
void CheckAndSetUnit(const G4String &unit, const G4String &category)

References G4VPrimitiveScorer::CheckAndSetUnit(), divideByArea, G4Exception(), G4VPrimitiveScorer::GetName(), G4VPrimitiveScorer::GetUnit(), JustWarning, G4VPrimitiveScorer::unitName, and G4VPrimitiveScorer::unitValue.

Referenced by G4PSCylinderSurfaceCurrent(), and G4PSCylinderSurfaceCurrent3D::G4PSCylinderSurfaceCurrent3D().

◆ SetVerboseLevel()

void G4VPrimitiveScorer::SetVerboseLevel ( G4int  vl)
inlineinherited

Definition at line 109 of file G4VPrimitiveScorer.hh.

109{ verboseLevel = vl; }

References G4VPrimitiveScorer::verboseLevel.

◆ Weighted()

void G4PSCylinderSurfaceCurrent::Weighted ( G4bool  flg = true)
inline

Definition at line 68 of file G4PSCylinderSurfaceCurrent.hh.

68{ weighted = flg; }

References weighted.

Field Documentation

◆ detector

G4MultiFunctionalDetector* G4VPrimitiveScorer::detector
protectedinherited

◆ divideByArea

G4bool G4PSCylinderSurfaceCurrent::divideByArea
private

Definition at line 95 of file G4PSCylinderSurfaceCurrent.hh.

Referenced by DivideByArea(), PrintAll(), ProcessHits(), and SetUnit().

◆ EvtMap

G4THitsMap<G4double>* G4PSCylinderSurfaceCurrent::EvtMap
private

Definition at line 93 of file G4PSCylinderSurfaceCurrent.hh.

Referenced by clear(), Initialize(), PrintAll(), and ProcessHits().

◆ fDirection

G4int G4PSCylinderSurfaceCurrent::fDirection
private

Definition at line 92 of file G4PSCylinderSurfaceCurrent.hh.

Referenced by ProcessHits().

◆ filter

G4VSDFilter* G4VPrimitiveScorer::filter
protectedinherited

◆ fNi

G4int G4VPrimitiveScorer::fNi
protectedinherited

◆ fNj

G4int G4VPrimitiveScorer::fNj
protectedinherited

Definition at line 131 of file G4VPrimitiveScorer.hh.

Referenced by G4PSCellCharge3D::G4PSCellCharge3D(), G4PSCellFlux3D::G4PSCellFlux3D(), G4PSCylinderSurfaceCurrent3D::G4PSCylinderSurfaceCurrent3D(), G4PSCylinderSurfaceFlux3D::G4PSCylinderSurfaceFlux3D(), G4PSDoseDeposit3D::G4PSDoseDeposit3D(), G4PSEnergyDeposit3D::G4PSEnergyDeposit3D(), G4PSFlatSurfaceCurrent3D::G4PSFlatSurfaceCurrent3D(), G4PSFlatSurfaceFlux3D::G4PSFlatSurfaceFlux3D(), G4PSMinKinEAtGeneration3D::G4PSMinKinEAtGeneration3D(), G4PSNofCollision3D::G4PSNofCollision3D(), G4PSNofSecondary3D::G4PSNofSecondary3D(), G4PSNofStep3D::G4PSNofStep3D(), G4PSPassageCellCurrent3D::G4PSPassageCellCurrent3D(), G4PSPassageCellFlux3D::G4PSPassageCellFlux3D(), G4PSPassageTrackLength3D::G4PSPassageTrackLength3D(), G4PSPopulation3D::G4PSPopulation3D(), G4PSSphereSurfaceCurrent3D::G4PSSphereSurfaceCurrent3D(), G4PSSphereSurfaceFlux3D::G4PSSphereSurfaceFlux3D(), G4PSStepChecker3D::G4PSStepChecker3D(), G4PSTermination3D::G4PSTermination3D(), G4PSTrackCounter3D::G4PSTrackCounter3D(), G4PSTrackLength3D::G4PSTrackLength3D(), G4PSVolumeFlux3D::G4PSVolumeFlux3D(), G4PSCellCharge3D::GetIndex(), G4PSCellFlux3D::GetIndex(), G4PSCylinderSurfaceCurrent3D::GetIndex(), G4PSCylinderSurfaceFlux3D::GetIndex(), G4PSDoseDeposit3D::GetIndex(), G4PSEnergyDeposit3D::GetIndex(), G4PSFlatSurfaceCurrent3D::GetIndex(), G4PSFlatSurfaceFlux3D::GetIndex(), G4PSMinKinEAtGeneration3D::GetIndex(), G4PSNofCollision3D::GetIndex(), G4PSNofSecondary3D::GetIndex(), G4PSNofStep3D::GetIndex(), G4PSPassageCellCurrent3D::GetIndex(), G4PSPassageCellFlux3D::GetIndex(), G4PSPassageTrackLength3D::GetIndex(), G4PSPopulation3D::GetIndex(), G4PSSphereSurfaceCurrent3D::GetIndex(), G4PSSphereSurfaceFlux3D::GetIndex(), G4PSStepChecker3D::GetIndex(), G4PSTermination3D::GetIndex(), G4PSTrackCounter3D::GetIndex(), G4PSTrackLength3D::GetIndex(), G4PSVolumeFlux3D::GetIndex(), and G4VPrimitiveScorer::SetNijk().

◆ fNk

G4int G4VPrimitiveScorer::fNk
protectedinherited

Definition at line 131 of file G4VPrimitiveScorer.hh.

Referenced by G4PSCellCharge3D::G4PSCellCharge3D(), G4PSCellFlux3D::G4PSCellFlux3D(), G4PSCylinderSurfaceCurrent3D::G4PSCylinderSurfaceCurrent3D(), G4PSCylinderSurfaceFlux3D::G4PSCylinderSurfaceFlux3D(), G4PSDoseDeposit3D::G4PSDoseDeposit3D(), G4PSEnergyDeposit3D::G4PSEnergyDeposit3D(), G4PSFlatSurfaceCurrent3D::G4PSFlatSurfaceCurrent3D(), G4PSFlatSurfaceFlux3D::G4PSFlatSurfaceFlux3D(), G4PSMinKinEAtGeneration3D::G4PSMinKinEAtGeneration3D(), G4PSNofCollision3D::G4PSNofCollision3D(), G4PSNofSecondary3D::G4PSNofSecondary3D(), G4PSNofStep3D::G4PSNofStep3D(), G4PSPassageCellCurrent3D::G4PSPassageCellCurrent3D(), G4PSPassageCellFlux3D::G4PSPassageCellFlux3D(), G4PSPassageTrackLength3D::G4PSPassageTrackLength3D(), G4PSPopulation3D::G4PSPopulation3D(), G4PSSphereSurfaceCurrent3D::G4PSSphereSurfaceCurrent3D(), G4PSSphereSurfaceFlux3D::G4PSSphereSurfaceFlux3D(), G4PSStepChecker3D::G4PSStepChecker3D(), G4PSTermination3D::G4PSTermination3D(), G4PSTrackCounter3D::G4PSTrackCounter3D(), G4PSTrackLength3D::G4PSTrackLength3D(), G4PSVolumeFlux3D::G4PSVolumeFlux3D(), G4PSCellCharge3D::GetIndex(), G4PSCellFlux3D::GetIndex(), G4PSCylinderSurfaceCurrent3D::GetIndex(), G4PSCylinderSurfaceFlux3D::GetIndex(), G4PSDoseDeposit3D::GetIndex(), G4PSEnergyDeposit3D::GetIndex(), G4PSFlatSurfaceCurrent3D::GetIndex(), G4PSFlatSurfaceFlux3D::GetIndex(), G4PSMinKinEAtGeneration3D::GetIndex(), G4PSNofCollision3D::GetIndex(), G4PSNofSecondary3D::GetIndex(), G4PSNofStep3D::GetIndex(), G4PSPassageCellCurrent3D::GetIndex(), G4PSPassageCellFlux3D::GetIndex(), G4PSPassageTrackLength3D::GetIndex(), G4PSPopulation3D::GetIndex(), G4PSSphereSurfaceCurrent3D::GetIndex(), G4PSSphereSurfaceFlux3D::GetIndex(), G4PSStepChecker3D::GetIndex(), G4PSTermination3D::GetIndex(), G4PSTrackCounter3D::GetIndex(), G4PSTrackLength3D::GetIndex(), G4PSVolumeFlux3D::GetIndex(), and G4VPrimitiveScorer::SetNijk().

◆ HCID

G4int G4PSCylinderSurfaceCurrent::HCID
private

Definition at line 91 of file G4PSCylinderSurfaceCurrent.hh.

Referenced by Initialize().

◆ hitIDMap

std::map<G4int, G4int> G4VPrimitivePlotter::hitIDMap
protectedinherited

◆ indexDepth

G4int G4VPrimitiveScorer::indexDepth
protectedinherited

◆ primitiveName

G4String G4VPrimitiveScorer::primitiveName
protectedinherited

◆ unitName

G4String G4VPrimitiveScorer::unitName
protectedinherited

◆ unitValue

G4double G4VPrimitiveScorer::unitValue
protectedinherited

◆ verboseLevel

G4int G4VPrimitiveScorer::verboseLevel
protectedinherited

◆ weighted

G4bool G4PSCylinderSurfaceCurrent::weighted
private

Definition at line 94 of file G4PSCylinderSurfaceCurrent.hh.

Referenced by ProcessHits(), and Weighted().


The documentation for this class was generated from the following files: