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

#include <G4PSDoseDeposit3D.hh>

Inheritance diagram for G4PSDoseDeposit3D:
G4PSDoseDeposit G4VPrimitivePlotter G4VPrimitiveScorer G4PSDoseDepositForCylinder3D

Public Member Functions

virtual void clear ()
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
 G4PSDoseDeposit3D (G4String name, const G4String &unit, G4int ni=1, G4int nj=1, G4int nk=1, G4int depi=2, G4int depj=1, G4int depk=0)
 
 G4PSDoseDeposit3D (G4String name, G4int ni=1, G4int nj=1, G4int nk=1, G4int depi=2, G4int depj=1, G4int depk=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)
 
virtual ~G4PSDoseDeposit3D ()
 

Protected Member Functions

void CheckAndSetUnit (const G4String &unit, const G4String &category)
 
G4VSolidComputeCurrentSolid (G4Step *aStep)
 
G4VSolidComputeSolid (G4Step *aStep, G4int replicaIdx)
 
virtual G4double ComputeVolume (G4Step *, G4int idx)
 
virtual G4int GetIndex (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
 
G4int fDepthi
 
G4int fDepthj
 
G4int fDepthk
 
G4int HCID
 

Detailed Description

Definition at line 43 of file G4PSDoseDeposit3D.hh.

Constructor & Destructor Documentation

◆ G4PSDoseDeposit3D() [1/2]

G4PSDoseDeposit3D::G4PSDoseDeposit3D ( G4String  name,
G4int  ni = 1,
G4int  nj = 1,
G4int  nk = 1,
G4int  depi = 2,
G4int  depj = 1,
G4int  depk = 0 
)

Definition at line 40 of file G4PSDoseDeposit3D.cc.

43 , fDepthi(di)
44 , fDepthj(dj)
45 , fDepthk(dk)
46{
47 fNi = ni;
48 fNj = nj;
49 fNk = nk;
50}
G4PSDoseDeposit(G4String name, G4int depth=0)
const char * name(G4int ptype)

References G4VPrimitiveScorer::fNi, G4VPrimitiveScorer::fNj, and G4VPrimitiveScorer::fNk.

◆ G4PSDoseDeposit3D() [2/2]

G4PSDoseDeposit3D::G4PSDoseDeposit3D ( G4String  name,
const G4String unit,
G4int  ni = 1,
G4int  nj = 1,
G4int  nk = 1,
G4int  depi = 2,
G4int  depj = 1,
G4int  depk = 0 
)

Definition at line 52 of file G4PSDoseDeposit3D.cc.

56 , fDepthi(di)
57 , fDepthj(dj)
58 , fDepthk(dk)
59{
60 fNi = ni;
61 fNj = nj;
62 fNk = nk;
63 SetUnit(unit);
64}
virtual void SetUnit(const G4String &unit)

References G4VPrimitiveScorer::fNi, G4VPrimitiveScorer::fNj, G4VPrimitiveScorer::fNk, and G4PSDoseDeposit::SetUnit().

◆ ~G4PSDoseDeposit3D()

G4PSDoseDeposit3D::~G4PSDoseDeposit3D ( )
virtual

Definition at line 66 of file G4PSDoseDeposit3D.cc.

66{ ; }

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 G4PSDoseDeposit::clear ( )
virtualinherited

Reimplemented from G4VPrimitiveScorer.

Definition at line 119 of file G4PSDoseDeposit.cc.

119{ EvtMap->clear(); }
G4THitsMap< G4double > * EvtMap
void clear()
Definition: G4THitsMap.hh:524

References G4VTHitsMap< T, Map_t >::clear(), and G4PSDoseDeposit::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().

◆ ComputeVolume()

G4double G4PSDoseDeposit::ComputeVolume ( G4Step aStep,
G4int  idx 
)
protectedvirtualinherited

Reimplemented in G4PSDoseDepositForCylinder3D.

Definition at line 142 of file G4PSDoseDeposit.cc.

143{
144 G4VSolid* solid = ComputeSolid(aStep, idx);
145 return solid->GetCubicVolume();
146}
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188

References G4VPrimitiveScorer::ComputeSolid(), and G4VSolid::GetCubicVolume().

Referenced by G4PSDoseDeposit::ProcessHits().

◆ DrawAll()

void G4PSDoseDeposit::DrawAll ( )
virtualinherited

Reimplemented from G4VPrimitiveScorer.

Definition at line 121 of file G4PSDoseDeposit.cc.

121{ ; }

◆ EndOfEvent()

void G4PSDoseDeposit::EndOfEvent ( G4HCofThisEvent )
virtualinherited

Reimplemented from G4VPrimitiveScorer.

Definition at line 117 of file G4PSDoseDeposit.cc.

117{ ; }

◆ 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(), 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 G4PSDoseDeposit3D::GetIndex ( G4Step aStep)
protectedvirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 68 of file G4PSDoseDeposit3D.cc.

69{
70 const G4VTouchable* touchable = aStep->GetPreStepPoint()->GetTouchable();
71 G4int i = touchable->GetReplicaNumber(fDepthi);
72 G4int j = touchable->GetReplicaNumber(fDepthj);
73 G4int k = touchable->GetReplicaNumber(fDepthk);
74
75 if(i < 0 || j < 0 || k < 0)
76 {
78 ED << "GetReplicaNumber is negative" << G4endl
79 << "touchable->GetReplicaNumber(fDepthi) returns i,j,k = " << i << ","
80 << j << "," << k << " for volume "
81 << touchable->GetVolume(fDepthi)->GetName() << ","
82 << touchable->GetVolume(fDepthj)->GetName() << ","
83 << touchable->GetVolume(fDepthk)->GetName() << G4endl;
84 G4Exception("G4PSDoseDeposit3D::GetIndex", "DetPS0005", JustWarning, ED);
85 }
86
87 return i * fNj * fNk + j * fNk + k;
88}
const G4String & GetName() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:50

References fDepthi, fDepthj, fDepthk, G4VPrimitiveScorer::fNj, G4VPrimitiveScorer::fNk, G4endl, G4Exception(), G4VPhysicalVolume::GetName(), G4Step::GetPreStepPoint(), G4VTouchable::GetReplicaNumber(), G4StepPoint::GetTouchable(), G4VTouchable::GetVolume(), and JustWarning.

◆ 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(), 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(), 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(), 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().

◆ 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 G4PSDoseDeposit::Initialize ( G4HCofThisEvent HCE)
virtualinherited

◆ 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 G4PSDoseDeposit::PrintAll ( )
virtualinherited

Reimplemented from G4VPrimitiveScorer.

Definition at line 123 of file G4PSDoseDeposit.cc.

124{
125 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
126 G4cout << " PrimitiveScorer " << GetName() << G4endl;
127 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
128 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
129 for(; itr != EvtMap->GetMap()->end(); itr++)
130 {
131 G4cout << " copy no.: " << itr->first
132 << " dose deposit: " << *(itr->second) / GetUnitValue() << " ["
133 << GetUnit() << "]" << G4endl;
134 }
135}
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, G4VTHitsMap< T, Map_t >::entries(), G4PSDoseDeposit::EvtMap, G4cout, G4endl, G4VTHitsMap< T, Map_t >::GetMap(), G4VPrimitiveScorer::GetName(), G4VSensitiveDetector::GetName(), G4VPrimitiveScorer::GetUnit(), and G4VPrimitiveScorer::GetUnitValue().

◆ ProcessHits()

G4bool G4PSDoseDeposit::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtualinherited

Implements G4VPrimitiveScorer.

Definition at line 67 of file G4PSDoseDeposit.cc.

68{
69 G4double edep = aStep->GetTotalEnergyDeposit();
70 if(edep == 0.)
71 return FALSE;
72
74 ->GetReplicaNumber(indexDepth);
75 G4double cubicVolume = ComputeVolume(aStep, idx);
76
77 G4double density = aStep->GetTrack()
78 ->GetStep()
80 ->GetMaterial()
81 ->GetDensity();
82 G4double dose = edep / (density * cubicVolume);
83 G4double wei = aStep->GetPreStepPoint()->GetWeight();
84 G4int index = GetIndex(aStep);
85 G4double dosew = dose * wei;
86 EvtMap->add(index, dosew);
87
88 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
89 {
90 auto filler = G4VScoreHistFiller::Instance();
91 if(!filler)
92 {
94 "G4PSDoseDeposit::ProcessHits", "SCORER0123", JustWarning,
95 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
96 }
97 else
98 {
99 filler->FillH1(hitIDMap[index], dose, wei);
100 }
101 }
102
103 return TRUE;
104}
double G4double
Definition: G4Types.hh:83
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
G4double GetDensity() const
Definition: G4Material.hh:176
virtual G4double ComputeVolume(G4Step *, G4int idx)
G4Material * GetMaterial() const
G4double GetWeight() const
G4Track * GetTrack() const
G4double GetTotalEnergyDeposit() const
const G4Step * GetStep() 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(), G4PSDoseDeposit::ComputeVolume(), G4PSDoseDeposit::EvtMap, FALSE, G4Exception(), G4Material::GetDensity(), G4VPrimitiveScorer::GetIndex(), G4StepPoint::GetMaterial(), G4Step::GetPreStepPoint(), G4Track::GetStep(), G4Step::GetTotalEnergyDeposit(), G4StepPoint::GetTouchable(), G4Step::GetTrack(), G4StepPoint::GetWeight(), G4VPrimitivePlotter::hitIDMap, G4VPrimitiveScorer::indexDepth, G4VScoreHistFiller::Instance(), JustWarning, and TRUE.

◆ 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

Definition at line 133 of file G4VPrimitiveScorer.hh.

134 {
135 fNi = i;
136 fNj = j;
137 fNk = k;
138 }

References G4VPrimitiveScorer::fNi, G4VPrimitiveScorer::fNj, and G4VPrimitiveScorer::fNk.

Referenced by G4VScoringMesh::SetPrimitiveScorer().

◆ SetUnit()

void G4PSDoseDeposit::SetUnit ( const G4String unit)
virtualinherited

Definition at line 137 of file G4PSDoseDeposit.cc.

138{
139 CheckAndSetUnit(unit, "Dose");
140}
void CheckAndSetUnit(const G4String &unit, const G4String &category)

References G4VPrimitiveScorer::CheckAndSetUnit().

Referenced by G4PSDoseDeposit::G4PSDoseDeposit(), G4PSDoseDeposit3D(), and G4ScoreQuantityMessenger::SetNewValue().

◆ SetVerboseLevel()

void G4VPrimitiveScorer::SetVerboseLevel ( G4int  vl)
inlineinherited

Definition at line 109 of file G4VPrimitiveScorer.hh.

109{ verboseLevel = vl; }

References G4VPrimitiveScorer::verboseLevel.

Field Documentation

◆ detector

G4MultiFunctionalDetector* G4VPrimitiveScorer::detector
protectedinherited

◆ EvtMap

G4THitsMap<G4double>* G4PSDoseDeposit::EvtMap
privateinherited

◆ fDepthi

G4int G4PSDoseDeposit3D::fDepthi
private

Definition at line 58 of file G4PSDoseDeposit3D.hh.

Referenced by GetIndex().

◆ fDepthj

G4int G4PSDoseDeposit3D::fDepthj
private

Definition at line 58 of file G4PSDoseDeposit3D.hh.

Referenced by GetIndex().

◆ fDepthk

G4int G4PSDoseDeposit3D::fDepthk
private

Definition at line 58 of file G4PSDoseDeposit3D.hh.

Referenced by GetIndex().

◆ 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(), 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(), 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(), 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(), 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 G4PSDoseDeposit::HCID
privateinherited

Definition at line 68 of file G4PSDoseDeposit.hh.

Referenced by G4PSDoseDeposit::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

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