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

#include <G4PSTrackLength.hh>

Inheritance diagram for G4PSTrackLength:
G4VPrimitiveScorer G4PSTrackLength3D

Public Member Functions

virtual void clear ()
 
void DivideByVelocity (G4bool flg=true)
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
 G4PSTrackLength (G4String name, const G4String &unit, G4int depth=0)
 
 G4PSTrackLength (G4String name, G4int depth=0)
 
G4int GetCollectionID (G4int)
 
G4VSDFilterGetFilter () const
 
G4MultiFunctionalDetectorGetMultiFunctionalDetector () const
 
G4String GetName () const
 
const G4StringGetUnit () const
 
G4double GetUnitValue () const
 
G4int GetVerboseLevel () const
 
virtual void Initialize (G4HCofThisEvent *)
 
void MultiplyKineticEnergy (G4bool flg=true)
 
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 ~G4PSTrackLength ()
 

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 *)
 
virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
 

Protected Attributes

G4MultiFunctionalDetectordetector
 
G4VSDFilterfilter
 
G4int fNi
 
G4int fNj
 
G4int fNk
 
G4int indexDepth
 
G4String primitiveName
 
G4String unitName
 
G4double unitValue
 
G4int verboseLevel
 

Private Member Functions

G4bool HitPrimitive (G4Step *aStep, G4TouchableHistory *ROhis)
 

Private Attributes

G4bool divideByVelocity
 
G4THitsMap< G4double > * EvtMap
 
G4int HCID
 
G4bool multiplyKinE
 
G4bool weighted
 

Detailed Description

Definition at line 48 of file G4PSTrackLength.hh.

Constructor & Destructor Documentation

◆ G4PSTrackLength() [1/2]

G4PSTrackLength::G4PSTrackLength ( G4String  name,
G4int  depth = 0 
)

Definition at line 43 of file G4PSTrackLength.cc.

44 : G4VPrimitiveScorer(name, depth)
45 , HCID(-1)
46 , EvtMap(0)
47 , weighted(false)
48 , multiplyKinE(false)
49 , divideByVelocity(false)
50{
52 SetUnit("mm");
53}
virtual void SetUnit(const G4String &unit)
virtual void DefineUnitAndCategory()
G4THitsMap< G4double > * EvtMap
G4VPrimitiveScorer(G4String name, G4int depth=0)
const char * name(G4int ptype)

References DefineUnitAndCategory(), and SetUnit().

◆ G4PSTrackLength() [2/2]

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

Definition at line 55 of file G4PSTrackLength.cc.

57 : G4VPrimitiveScorer(name, depth)
58 , HCID(-1)
59 , EvtMap(0)
60 , weighted(false)
61 , multiplyKinE(false)
62 , divideByVelocity(false)
63{
65 SetUnit(unit);
66}

References DefineUnitAndCategory(), and SetUnit().

◆ ~G4PSTrackLength()

G4PSTrackLength::~G4PSTrackLength ( )
virtual

Definition at line 68 of file G4PSTrackLength.cc.

68{ ; }

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(), G4PSCylinderSurfaceCurrent::SetUnit(), G4PSCylinderSurfaceFlux::SetUnit(), G4PSDoseDeposit::SetUnit(), G4PSEnergyDeposit::SetUnit(), G4PSFlatSurfaceCurrent::SetUnit(), G4PSFlatSurfaceFlux::SetUnit(), G4PSMinKinEAtGeneration::SetUnit(), G4PSPassageCellFlux::SetUnit(), G4PSPassageTrackLength::SetUnit(), G4PSSphereSurfaceCurrent::SetUnit(), G4PSSphereSurfaceFlux::SetUnit(), and SetUnit().

◆ clear()

void G4PSTrackLength::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 112 of file G4PSTrackLength.cc.

112{ 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 G4PSCylinderSurfaceCurrent::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 G4PSTrackLength::DefineUnitAndCategory ( )
protectedvirtual

Definition at line 198 of file G4PSTrackLength.cc.

199{
200 // EnergyFlux
201 new G4UnitDefinition("eV_second", "eV_s", "EnergyFlux", (eV * second));
202 new G4UnitDefinition("keV_second", "keV_s", "EnergyFlux", (keV * second));
203 new G4UnitDefinition("MeV_second", "MeV_s", "EnergyFlux", (MeV * second));
204 new G4UnitDefinition("eV_millisecond", "eV_ms", "EnergyFlux", (eV * ms));
205 new G4UnitDefinition("keV_millisecond", "keV_ms", "EnergyFlux", (keV * ms));
206 new G4UnitDefinition("MeV_millisecond", "MeV_ms", "EnergyFlux", (MeV * ms));
207 // EnergyFlow
208 new G4UnitDefinition("eV_millimeter", "eV_mm", "EnergyFlow", (eV * mm));
209 new G4UnitDefinition("keV_millimeter", "keV_mm", "EnergyFlow", (keV * mm));
210 new G4UnitDefinition("MeV_millimeter", "MeV_mm", "EnergyFlow", (MeV * mm));
211 new G4UnitDefinition("eV_centimeter", "eV_cm", "EnergyFlow", (eV * cm));
212 new G4UnitDefinition("keV_centimeter", "keV_cm", "EnergyFlow", (keV * cm));
213 new G4UnitDefinition("MeV_centimeter", "MeV_cm", "EnergyFlow", (MeV * cm));
214 new G4UnitDefinition("eV_meter", "eV_m", "EnergyFlow", (eV * m));
215 new G4UnitDefinition("keV_meter", "keV_m", "EnergyFlow", (keV * m));
216 new G4UnitDefinition("MeV_meter", "MeV_m", "EnergyFlow", (MeV * m));
217}
static constexpr double ms
Definition: G4SIunits.hh:155
static constexpr double m
Definition: G4SIunits.hh:109
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double second
Definition: G4SIunits.hh:137
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double cm
Definition: G4SIunits.hh:99

References cm, eV, keV, m, MeV, mm, ms, and second.

Referenced by G4PSTrackLength().

◆ DivideByVelocity()

void G4PSTrackLength::DivideByVelocity ( G4bool  flg = true)

Definition at line 77 of file G4PSTrackLength.cc.

78{
79 divideByVelocity = flg;
80 // Default unit is set according to flags.
81 SetUnit("");
82}

References divideByVelocity, and SetUnit().

◆ DrawAll()

void G4PSTrackLength::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 114 of file G4PSTrackLength.cc.

114{ ; }

◆ EndOfEvent()

void G4PSTrackLength::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 110 of file G4PSTrackLength.cc.

110{ ; }

◆ 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(), G4PSCylinderSurfaceCurrent::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(), 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(), G4PSCylinderSurfaceCurrent::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(), 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(), G4PSCylinderSurfaceCurrent::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(), Initialize(), G4PSVolumeFlux::Initialize(), G4VScoringMesh::List(), G4PSCellCharge::PrintAll(), G4PSCellFlux::PrintAll(), G4PSCylinderSurfaceCurrent::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(), PrintAll(), G4PSVolumeFlux::PrintAll(), G4MultiFunctionalDetector::RegisterPrimitive(), G4MultiFunctionalDetector::RemovePrimitive(), G4VScoringMesh::SetFilter(), G4VScoringMesh::SetPrimitiveScorer(), G4PSCylinderSurfaceCurrent::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().

◆ 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 G4PSTrackLength::Initialize ( G4HCofThisEvent HCE)
virtual

◆ MultiplyKineticEnergy()

void G4PSTrackLength::MultiplyKineticEnergy ( G4bool  flg = true)

Definition at line 70 of file G4PSTrackLength.cc.

71{
72 multiplyKinE = flg;
73 // Default unit is set according to flags.
74 SetUnit("");
75}

References multiplyKinE, and SetUnit().

◆ PrintAll()

void G4PSTrackLength::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 116 of file G4PSTrackLength.cc.

117{
118 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
119 G4cout << " PrimitiveScorer " << GetName() << G4endl;
120 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
121 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
122 for(; itr != EvtMap->GetMap()->end(); itr++)
123 {
124 G4cout << " copy no.: " << itr->first;
125 if(multiplyKinE)
126 {
128 G4cout << " EnergyFlux: ";
129 else
130 G4cout << " EnergyFlow: ";
131 }
132 else
133 {
135 G4cout << " Time: ";
136 else
137 G4cout << " Length: ";
138 }
139 G4cout << *(itr->second) / GetUnitValue() << " [" << GetUnit() << "]";
140 G4cout << G4endl;
141 }
142}
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, divideByVelocity, G4VTHitsMap< T, Map_t >::entries(), EvtMap, G4cout, G4endl, G4VTHitsMap< T, Map_t >::GetMap(), G4VPrimitiveScorer::GetName(), G4VSensitiveDetector::GetName(), G4VPrimitiveScorer::GetUnit(), G4VPrimitiveScorer::GetUnitValue(), and multiplyKinE.

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 84 of file G4PSTrackLength.cc.

85{
86 G4double trklength = aStep->GetStepLength();
87 if(trklength == 0.)
88 return FALSE;
89 if(weighted)
90 trklength *= aStep->GetPreStepPoint()->GetWeight();
91 if(multiplyKinE)
92 trklength *= aStep->GetPreStepPoint()->GetKineticEnergy();
94 trklength /= aStep->GetPreStepPoint()->GetVelocity();
95 G4int index = GetIndex(aStep);
96 EvtMap->add(index, trklength);
97 return TRUE;
98}
double G4double
Definition: G4Types.hh:83
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
G4double GetVelocity() const
G4double GetKineticEnergy() const
G4double GetWeight() const
G4double GetStepLength() const
virtual G4int GetIndex(G4Step *)
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

References G4VTHitsMap< T, Map_t >::add(), divideByVelocity, EvtMap, FALSE, G4VPrimitiveScorer::GetIndex(), G4StepPoint::GetKineticEnergy(), G4Step::GetPreStepPoint(), G4Step::GetStepLength(), G4StepPoint::GetVelocity(), G4StepPoint::GetWeight(), multiplyKinE, 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 G4PSTrackLength::SetUnit ( const G4String unit)
virtual

Definition at line 144 of file G4PSTrackLength.cc.

145{
146 if(multiplyKinE)
147 {
149 {
150 if(unit == "")
151 {
152 CheckAndSetUnit("MeV_second", "EnergyFlux");
153 }
154 else
155 {
156 CheckAndSetUnit(unit, "EnergyFlux");
157 }
158 }
159 else
160 {
161 if(unit == "")
162 {
163 CheckAndSetUnit("MeV_mm", "EnergyFlow");
164 }
165 else
166 {
167 CheckAndSetUnit(unit, "EnergyFlow");
168 }
169 }
170 }
171 else
172 {
174 {
175 if(unit == "")
176 {
177 CheckAndSetUnit("second", "Time");
178 }
179 else
180 {
181 CheckAndSetUnit(unit, "Time");
182 }
183 }
184 else
185 {
186 if(unit == "")
187 {
188 CheckAndSetUnit("mm", "Length");
189 }
190 else
191 {
192 CheckAndSetUnit(unit, "Length");
193 }
194 }
195 }
196}
void CheckAndSetUnit(const G4String &unit, const G4String &category)

References G4VPrimitiveScorer::CheckAndSetUnit(), divideByVelocity, and multiplyKinE.

Referenced by DivideByVelocity(), G4PSTrackLength(), G4PSTrackLength3D::G4PSTrackLength3D(), and MultiplyKineticEnergy().

◆ SetVerboseLevel()

void G4VPrimitiveScorer::SetVerboseLevel ( G4int  vl)
inlineinherited

Definition at line 109 of file G4VPrimitiveScorer.hh.

109{ verboseLevel = vl; }

References G4VPrimitiveScorer::verboseLevel.

◆ Weighted()

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

Definition at line 55 of file G4PSTrackLength.hh.

55{ weighted = flg; }

References weighted.

Field Documentation

◆ detector

G4MultiFunctionalDetector* G4VPrimitiveScorer::detector
protectedinherited

◆ divideByVelocity

G4bool G4PSTrackLength::divideByVelocity
private

Definition at line 84 of file G4PSTrackLength.hh.

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

◆ EvtMap

G4THitsMap<G4double>* G4PSTrackLength::EvtMap
private

Definition at line 81 of file G4PSTrackLength.hh.

Referenced by clear(), Initialize(), PrintAll(), and 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 G4PSTrackLength::HCID
private

Definition at line 80 of file G4PSTrackLength.hh.

Referenced by Initialize().

◆ indexDepth

G4int G4VPrimitiveScorer::indexDepth
protectedinherited

◆ multiplyKinE

G4bool G4PSTrackLength::multiplyKinE
private

Definition at line 83 of file G4PSTrackLength.hh.

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

◆ primitiveName

G4String G4VPrimitiveScorer::primitiveName
protectedinherited

◆ unitName

G4String G4VPrimitiveScorer::unitName
protectedinherited

◆ unitValue

G4double G4VPrimitiveScorer::unitValue
protectedinherited

◆ verboseLevel

G4int G4VPrimitiveScorer::verboseLevel
protectedinherited

◆ weighted

G4bool G4PSTrackLength::weighted
private

Definition at line 82 of file G4PSTrackLength.hh.

Referenced by ProcessHits(), and Weighted().


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