Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PSPassageCellFlux Class Reference

#include <G4PSPassageCellFlux.hh>

+ Inheritance diagram for G4PSPassageCellFlux:

Public Member Functions

 G4PSPassageCellFlux (G4String name, G4int depth=0)
 
 G4PSPassageCellFlux (G4String name, const G4String &unit, G4int depth=0)
 
 ~G4PSPassageCellFlux () override=default
 
void Weighted (G4bool flg=true)
 
void Initialize (G4HCofThisEvent *) override
 
void clear () override
 
void PrintAll () override
 
virtual void SetUnit (const G4String &unit)
 
- Public Member Functions inherited from G4VPrimitivePlotter
 ~G4VPrimitivePlotter () override=default
 
void Plot (G4int copyNo, G4int histID)
 
G4int GetNumberOfHist () const
 
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
- Public Member Functions inherited from G4VPrimitiveScorer
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
virtual ~G4VPrimitiveScorer ()=default
 
G4int GetCollectionID (G4int)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void DrawAll ()
 
void SetUnit (const G4String &unit)
 
const G4StringGetUnit () const
 
G4double GetUnitValue () const
 
void SetMultiFunctionalDetector (G4MultiFunctionalDetector *d)
 
G4MultiFunctionalDetectorGetMultiFunctionalDetector () const
 
G4String GetName () const
 
void SetFilter (G4VSDFilter *f)
 
G4VSDFilterGetFilter () const
 
void SetVerboseLevel (G4int vl)
 
G4int GetVerboseLevel () const
 
void SetNijk (G4int i, G4int j, G4int k)
 

Protected Member Functions

G4bool ProcessHits (G4Step *, G4TouchableHistory *) override
 
virtual G4bool IsPassed (G4Step *)
 
virtual G4double ComputeVolume (G4Step *, G4int idx)
 
virtual void DefineUnitAndCategory ()
 
- Protected Member Functions inherited from G4VPrimitiveScorer
G4VSolidComputeSolid (G4Step *aStep, G4int replicaIdx)
 
G4VSolidComputeCurrentSolid (G4Step *aStep)
 
virtual G4int GetIndex (G4Step *)
 
void CheckAndSetUnit (const G4String &unit, const G4String &category)
 

Additional Inherited Members

- Protected Attributes inherited from G4VPrimitivePlotter
std::map< G4int, G4inthitIDMap
 
- Protected Attributes inherited from G4VPrimitiveScorer
G4String primitiveName
 
G4MultiFunctionalDetectordetector {nullptr}
 
G4VSDFilterfilter {nullptr}
 
G4int verboseLevel {0}
 
G4int indexDepth
 
G4String unitName {"NoUnit"}
 
G4double unitValue {1.0}
 
G4int fNi {0}
 
G4int fNj {0}
 
G4int fNk {0}
 

Detailed Description

Definition at line 53 of file G4PSPassageCellFlux.hh.

Constructor & Destructor Documentation

◆ G4PSPassageCellFlux() [1/2]

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

Definition at line 58 of file G4PSPassageCellFlux.cc.

59 : G4PSPassageCellFlux(name, "percm2", depth)
60{}
G4PSPassageCellFlux(G4String name, G4int depth=0)

◆ G4PSPassageCellFlux() [2/2]

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

Definition at line 62 of file G4PSPassageCellFlux.cc.

64 : G4VPrimitivePlotter(name, depth)
65 , HCID(-1)
66 , fCurrentTrkID(-1)
67 , fCellFlux(0)
68 , EvtMap(nullptr)
69 , weighted(true)
70{
72 SetUnit(unit);
73}
virtual void SetUnit(const G4String &unit)
virtual void DefineUnitAndCategory()

◆ ~G4PSPassageCellFlux()

G4PSPassageCellFlux::~G4PSPassageCellFlux ( )
overridedefault

Member Function Documentation

◆ clear()

void G4PSPassageCellFlux::clear ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 159 of file G4PSPassageCellFlux.cc.

159{ EvtMap->clear(); }

◆ ComputeVolume()

G4double G4PSPassageCellFlux::ComputeVolume ( G4Step * aStep,
G4int idx )
protectedvirtual

Reimplemented in G4PSPassageCellFluxForCylinder3D.

Definition at line 189 of file G4PSPassageCellFlux.cc.

190{
191 return ComputeSolid(aStep, idx)->GetCubicVolume();
192}
G4VSolid * ComputeSolid(G4Step *aStep, G4int replicaIdx)
virtual G4double GetCubicVolume()
Definition G4VSolid.cc:188

Referenced by ProcessHits().

◆ DefineUnitAndCategory()

void G4PSPassageCellFlux::DefineUnitAndCategory ( )
protectedvirtual

Definition at line 179 of file G4PSPassageCellFlux.cc.

180{
181 // Per Unit Surface
182 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
183 (1. / cm2));
184 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
185 (1. / mm2));
186 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
187}

Referenced by G4PSPassageCellFlux().

◆ Initialize()

void G4PSPassageCellFlux::Initialize ( G4HCofThisEvent * HCE)
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 149 of file G4PSPassageCellFlux.cc.

150{
151 fCurrentTrkID = -1;
152
153 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
154 if(HCID < 0)
155 HCID = GetCollectionID(0);
156 HCE->AddHitsCollection(HCID, EvtMap);
157}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector

◆ IsPassed()

G4bool G4PSPassageCellFlux::IsPassed ( G4Step * aStep)
protectedvirtual

Definition at line 108 of file G4PSPassageCellFlux.cc.

109{
110 G4bool Passed = false;
111
112 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
113 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
114
115 G4int trkid = aStep->GetTrack()->GetTrackID();
116 G4double trklength = aStep->GetStepLength();
117 if(weighted)
118 trklength *= aStep->GetPreStepPoint()->GetWeight();
119
120 if(IsEnter && IsExit)
121 { // Passed at one step
122 fCellFlux = trklength; // Track length is absolutely given.
123 Passed = true;
124 }
125 else if(IsEnter)
126 { // Enter a new geometry
127 fCurrentTrkID = trkid; // Resetting the current track.
128 fCellFlux = trklength;
129 }
130 else if(IsExit)
131 { // Exit a current geometry
132 if(fCurrentTrkID == trkid)
133 { // if the track is same as entered,
134 fCellFlux += trklength; // add the track length to current one.
135 Passed = true;
136 }
137 }
138 else
139 { // Inside geometry
140 if(fCurrentTrkID == trkid)
141 { // if the track is same as entered,
142 fCellFlux += trklength; // adding the track length to current one.
143 }
144 }
145
146 return Passed;
147}
@ fGeomBoundary
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4StepStatus GetStepStatus() const
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const

Referenced by ProcessHits().

◆ PrintAll()

void G4PSPassageCellFlux::PrintAll ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 161 of file G4PSPassageCellFlux.cc.

162{
163 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
164 G4cout << " PrimitiveScorer " << GetName() << G4endl;
165 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
166 for(const auto& [copy, flux] : *(EvtMap->GetMap()))
167 {
168 G4cout << " copy no.: " << copy
169 << " cell flux : " << *(flux) / GetUnitValue() << " ["
170 << GetUnit() << G4endl;
171 }
172}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetUnit() const
G4double GetUnitValue() const
Map_t * GetMap() const
size_t entries() const
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)

◆ ProcessHits()

G4bool G4PSPassageCellFlux::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )
overrideprotectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 75 of file G4PSPassageCellFlux.cc.

76{
77 if(IsPassed(aStep))
78 {
79 G4int idx =
81 ->GetReplicaNumber(indexDepth);
82 G4double cubicVolume = ComputeVolume(aStep, idx);
83
84 fCellFlux /= cubicVolume;
85 G4int index = GetIndex(aStep);
86 EvtMap->add(index, fCellFlux);
87
88 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend())
89 {
90 auto filler = G4VScoreHistFiller::Instance();
91 if(filler == nullptr)
92 {
94 "G4PSPassageCellFlux::ProcessHits", "SCORER0123", JustWarning,
95 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
96 }
97 else
98 {
99 filler->FillH1(hitIDMap[index],
100 aStep->GetPreStepPoint()->GetKineticEnergy(), fCellFlux);
101 }
102 }
103 }
104
105 return true;
106}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
virtual G4double ComputeVolume(G4Step *, G4int idx)
virtual G4bool IsPassed(G4Step *)
const G4VTouchable * GetTouchable() const
G4double GetKineticEnergy() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
size_t add(const G4int &key, U *&aHit) const

◆ SetUnit()

void G4PSPassageCellFlux::SetUnit ( const G4String & unit)
virtual

Definition at line 174 of file G4PSPassageCellFlux.cc.

175{
176 CheckAndSetUnit(unit, "Per Unit Surface");
177}
void CheckAndSetUnit(const G4String &unit, const G4String &category)

Referenced by G4PSPassageCellFlux(), G4PSPassageCellFlux3D::G4PSPassageCellFlux3D(), and G4ScoreQuantityMessenger::SetNewValue().

◆ Weighted()

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

Definition at line 60 of file G4PSPassageCellFlux.hh.

60{ weighted = flg; }

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