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

#include <G4PSPassageCellCurrent.hh>

Inheritance diagram for G4PSPassageCellCurrent:
G4VPrimitivePlotter G4VPrimitiveScorer G4PSPassageCellCurrent3D

Public Member Functions

virtual void clear ()
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
 G4PSPassageCellCurrent (G4String name, 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 ~G4PSPassageCellCurrent ()
 

Protected Member Functions

void CheckAndSetUnit (const G4String &unit, const G4String &category)
 
G4VSolidComputeCurrentSolid (G4Step *aStep)
 
G4VSolidComputeSolid (G4Step *aStep, G4int replicaIdx)
 
virtual G4int GetIndex (G4Step *)
 
virtual G4bool IsPassed (G4Step *)
 
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

G4THitsMap< G4double > * EvtMap
 
G4double fCurrent
 
G4int fCurrentTrkID
 
G4int HCID
 
G4bool weighted
 

Detailed Description

Definition at line 49 of file G4PSPassageCellCurrent.hh.

Constructor & Destructor Documentation

◆ G4PSPassageCellCurrent()

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

Definition at line 50 of file G4PSPassageCellCurrent.cc.

52 , HCID(-1)
53 , fCurrentTrkID(-1)
54 , fCurrent(0)
55 , EvtMap(0)
56 , weighted(true)
57{
58 SetUnit("");
59}
virtual void SetUnit(const G4String &unit)
G4THitsMap< G4double > * EvtMap
G4VPrimitivePlotter(G4String name, G4int depth=0)
const char * name(G4int ptype)

References SetUnit().

◆ ~G4PSPassageCellCurrent()

G4PSPassageCellCurrent::~G4PSPassageCellCurrent ( )
virtual

Definition at line 61 of file G4PSPassageCellCurrent.cc.

61{ ; }

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 G4PSTrackLength::SetUnit().

◆ clear()

void G4PSPassageCellCurrent::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 138 of file G4PSPassageCellCurrent.cc.

138{ 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().

◆ DrawAll()

void G4PSPassageCellCurrent::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 140 of file G4PSPassageCellCurrent.cc.

140{ ; }

◆ EndOfEvent()

void G4PSPassageCellCurrent::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 136 of file G4PSPassageCellCurrent.cc.

136{}

◆ 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(), 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(), G4PSCylinderSurfaceCurrent::ProcessHits(), G4PSCylinderSurfaceFlux::ProcessHits(), G4PSDoseDeposit::ProcessHits(), G4PSEnergyDeposit::ProcessHits(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSMinKinEAtGeneration::ProcessHits(), G4PSNofCollision::ProcessHits(), G4PSNofSecondary::ProcessHits(), G4PSNofStep::ProcessHits(), 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(), G4PSCylinderSurfaceCurrent::Initialize(), G4PSCylinderSurfaceFlux::Initialize(), G4PSDoseDeposit::Initialize(), G4PSEnergyDeposit::Initialize(), G4PSFlatSurfaceCurrent::Initialize(), G4PSFlatSurfaceFlux::Initialize(), G4PSMinKinEAtGeneration::Initialize(), G4PSNofCollision::Initialize(), G4PSNofSecondary::Initialize(), G4PSNofStep::Initialize(), 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(), G4PSCylinderSurfaceCurrent::PrintAll(), G4PSCylinderSurfaceFlux::PrintAll(), G4PSDoseDeposit::PrintAll(), G4PSEnergyDeposit::PrintAll(), G4PSFlatSurfaceCurrent::PrintAll(), G4PSFlatSurfaceFlux::PrintAll(), G4PSMinKinEAtGeneration::PrintAll(), G4PSNofCollision::PrintAll(), G4PSNofSecondary::PrintAll(), G4PSNofStep::PrintAll(), 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(), G4PSCylinderSurfaceCurrent::SetUnit(), G4PSCylinderSurfaceFlux::SetUnit(), G4PSFlatSurfaceCurrent::SetUnit(), G4PSFlatSurfaceFlux::SetUnit(), G4PSNofCollision::SetUnit(), G4PSNofSecondary::SetUnit(), G4PSNofStep::SetUnit(), 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 G4PSPassageCellCurrent::Initialize ( G4HCofThisEvent HCE)
virtual

◆ IsPassed()

G4bool G4PSPassageCellCurrent::IsPassed ( G4Step aStep)
protectedvirtual

Definition at line 93 of file G4PSPassageCellCurrent.cc.

94{
95 G4bool Passed = FALSE;
96
97 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
98 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
99
100 G4int trkid = aStep->GetTrack()->GetTrackID();
101
102 if(IsEnter && IsExit)
103 { // Passed at one step
104 Passed = TRUE;
105 }
106 else if(IsEnter)
107 { // Enter a new geometry
108 fCurrentTrkID = trkid; // Resetting the current track.
109 }
110 else if(IsExit)
111 { // Exit a current geometry
112 if(fCurrentTrkID == trkid)
113 {
114 Passed = TRUE; // if the track is same as entered.
115 }
116 }
117 else
118 { // Inside geometry
119 if(fCurrentTrkID == trkid)
120 { // Adding the track length to current one ,
121 }
122 }
123 return Passed;
124}
@ fGeomBoundary
Definition: G4StepStatus.hh:43
bool G4bool
Definition: G4Types.hh:86
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
G4StepStatus GetStepStatus() const
G4Track * GetTrack() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const

References FALSE, fCurrentTrkID, fGeomBoundary, G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetStepStatus(), G4Step::GetTrack(), G4Track::GetTrackID(), and TRUE.

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 G4PSPassageCellCurrent::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 142 of file G4PSPassageCellCurrent.cc.

143{
144 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
145 G4cout << " PrimitiveScorer " << GetName() << G4endl;
146 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
147 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
148 for(; itr != EvtMap->GetMap()->end(); itr++)
149 {
150 G4cout << " copy no.: " << itr->first
151 << " cell current : " << *(itr->second) << " [tracks] " << G4endl;
152 }
153}
G4GLOB_DLL std::ostream G4cout
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t entries() const
Definition: G4THitsMap.hh:158

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

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 63 of file G4PSPassageCellCurrent.cc.

64{
65 if(IsPassed(aStep))
66 {
67 fCurrent = 1.;
68 if(weighted)
70 G4int index = GetIndex(aStep);
71 EvtMap->add(index, fCurrent);
72
73 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
74 {
75 auto filler = G4VScoreHistFiller::Instance();
76 if(!filler)
77 {
79 "G4PSVolumeFlux::ProcessHits", "SCORER0123", JustWarning,
80 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
81 }
82 else
83 {
84 filler->FillH1(hitIDMap[index],
86 }
87 }
88 }
89
90 return TRUE;
91}
virtual G4bool IsPassed(G4Step *)
G4double GetKineticEnergy() const
G4double GetWeight() const
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

References G4VTHitsMap< T, Map_t >::add(), EvtMap, fCurrent, G4Exception(), G4VPrimitiveScorer::GetIndex(), G4StepPoint::GetKineticEnergy(), G4Step::GetPreStepPoint(), G4StepPoint::GetWeight(), G4VPrimitivePlotter::hitIDMap, G4VScoreHistFiller::Instance(), IsPassed(), JustWarning, 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 G4PSPassageCellCurrent::SetUnit ( const G4String unit)
virtual

Definition at line 155 of file G4PSPassageCellCurrent.cc.

156{
157 if(unit == "")
158 {
159 unitName = unit;
160 unitValue = 1.0;
161 }
162 else
163 {
164 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
165 GetUnit() + "] ) for " + GetName();
166 G4Exception("G4PSPassageCellCurrent::SetUnit", "DetPS0012", JustWarning,
167 msg);
168 }
169}

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

Referenced by G4PSPassageCellCurrent().

◆ SetVerboseLevel()

void G4VPrimitiveScorer::SetVerboseLevel ( G4int  vl)
inlineinherited

Definition at line 109 of file G4VPrimitiveScorer.hh.

109{ verboseLevel = vl; }

References G4VPrimitiveScorer::verboseLevel.

◆ Weighted()

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

Definition at line 55 of file G4PSPassageCellCurrent.hh.

55{ weighted = flg; }

References weighted.

Field Documentation

◆ detector

G4MultiFunctionalDetector* G4VPrimitiveScorer::detector
protectedinherited

◆ EvtMap

G4THitsMap<G4double>* G4PSPassageCellCurrent::EvtMap
private

Definition at line 78 of file G4PSPassageCellCurrent.hh.

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

◆ fCurrent

G4double G4PSPassageCellCurrent::fCurrent
private

Definition at line 77 of file G4PSPassageCellCurrent.hh.

Referenced by ProcessHits().

◆ fCurrentTrkID

G4int G4PSPassageCellCurrent::fCurrentTrkID
private

Definition at line 76 of file G4PSPassageCellCurrent.hh.

Referenced by Initialize(), and IsPassed().

◆ 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 G4PSPassageCellCurrent::HCID
private

Definition at line 75 of file G4PSPassageCellCurrent.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 G4PSPassageCellCurrent::weighted
private

Definition at line 79 of file G4PSPassageCellCurrent.hh.

Referenced by ProcessHits(), and Weighted().


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