Geant4 9.6.0
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)
 
virtual ~G4PSPassageCellFlux ()
 
void Weighted (G4bool flg=true)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
virtual void SetUnit (const G4String &unit)
 
- Public Member Functions inherited from G4VPrimitiveScorer
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
virtual ~G4VPrimitiveScorer ()
 
G4int GetCollectionID (G4int)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
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

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

Additional Inherited Members

- Protected Attributes inherited from G4VPrimitiveScorer
G4String primitiveName
 
G4MultiFunctionalDetectordetector
 
G4VSDFilterfilter
 
G4int verboseLevel
 
G4int indexDepth
 
G4String unitName
 
G4double unitValue
 
G4int fNi
 
G4int fNj
 
G4int fNk
 

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 56 of file G4PSPassageCellFlux.cc.

57 : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
58 weighted(true)
59{
61 SetUnit("percm2");
62}
virtual void SetUnit(const G4String &unit)
virtual void DefineUnitAndCategory()

◆ G4PSPassageCellFlux() [2/2]

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

Definition at line 64 of file G4PSPassageCellFlux.cc.

66 : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
67 weighted(true)
68{
70 SetUnit(unit);
71}

◆ ~G4PSPassageCellFlux()

G4PSPassageCellFlux::~G4PSPassageCellFlux ( )
virtual

Definition at line 73 of file G4PSPassageCellFlux.cc.

74{;}

Member Function Documentation

◆ clear()

void G4PSPassageCellFlux::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 137 of file G4PSPassageCellFlux.cc.

137 {
138 EvtMap->clear();
139}
void clear()
Definition: G4THitsMap.hh:209

◆ ComputeVolume()

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

Reimplemented in G4PSPassageCellFluxForCylinder3D.

Definition at line 171 of file G4PSPassageCellFlux.cc.

171 {
172
174 G4VPVParameterisation* physParam = physVol->GetParameterisation();
175 G4VSolid* solid = 0;
176 if(physParam)
177 { // for parameterized volume
178 if(idx<0)
179 {
181 ED << "Incorrect replica number --- GetReplicaNumber : " << idx << G4endl;
182 G4Exception("G4PSPassageCellFlux::ComputeVolume","DetPS0013",JustWarning,ED);
183 }
184 solid = physParam->ComputeSolid(idx, physVol);
185 solid->ComputeDimensions(physParam,idx,physVol);
186 }
187 else
188 { // for ordinary volume
189 solid = physVol->GetLogicalVolume()->GetSolid();
190 }
191
192 return solid->GetCubicVolume();
193}
@ JustWarning
#define G4endl
Definition: G4ios.hh:52
G4VSolid * GetSolid() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Referenced by ProcessHits().

◆ DefineUnitAndCategory()

void G4PSPassageCellFlux::DefineUnitAndCategory ( )
protectedvirtual

Definition at line 163 of file G4PSPassageCellFlux.cc.

163 {
164 // Per Unit Surface
165 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
166 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
167 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
168}

Referenced by G4PSPassageCellFlux().

◆ DrawAll()

void G4PSPassageCellFlux::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 141 of file G4PSPassageCellFlux.cc.

142{;}

◆ EndOfEvent()

void G4PSPassageCellFlux::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 134 of file G4PSPassageCellFlux.cc.

135{;}

◆ Initialize()

void G4PSPassageCellFlux::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 123 of file G4PSPassageCellFlux.cc.

124{
125 fCurrentTrkID = -1;
126
127 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
128 GetName());
129 if ( HCID < 0 ) HCID = GetCollectionID(0);
130 HCE->AddHitsCollection(HCID,EvtMap);
131
132}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ IsPassed()

G4bool G4PSPassageCellFlux::IsPassed ( G4Step aStep)
protectedvirtual

Definition at line 93 of file G4PSPassageCellFlux.cc.

93 {
94 G4bool Passed = FALSE;
95
96 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
97 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
98
99 G4int trkid = aStep->GetTrack()->GetTrackID();
100 G4double trklength = aStep->GetStepLength();
101 if ( weighted ) trklength *= aStep->GetPreStepPoint()->GetWeight();
102
103 if ( IsEnter &&IsExit ){ // Passed at one step
104 fCellFlux = trklength; // Track length is absolutely given.
105 Passed = TRUE;
106 }else if ( IsEnter ){ // Enter a new geometry
107 fCurrentTrkID = trkid; // Resetting the current track.
108 fCellFlux = trklength;
109 }else if ( IsExit ){ // Exit a current geometry
110 if ( fCurrentTrkID == trkid ) {// if the track is same as entered,
111 fCellFlux += trklength; // add the track length to current one.
112 Passed = TRUE;
113 }
114 }else{ // Inside geometry
115 if ( fCurrentTrkID == trkid ){ // if the track is same as entered,
116 fCellFlux += trklength; // adding the track length to current one.
117 }
118 }
119
120 return Passed;
121}
@ fGeomBoundary
Definition: G4StepStatus.hh:54
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4StepStatus GetStepStatus() const
G4double GetWeight() const
G4Track * GetTrack() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52

Referenced by ProcessHits().

◆ PrintAll()

void G4PSPassageCellFlux::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 144 of file G4PSPassageCellFlux.cc.

145{
146 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
147 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
148 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
149 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
150 for(; itr != EvtMap->GetMap()->end(); itr++) {
151 G4cout << " copy no.: " << itr->first
152 << " cell flux : " << *(itr->second)/GetUnitValue()
153 << " [" << GetUnit()
154 << G4endl;
155 }
156}
G4DLLIMPORT std::ostream G4cout
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
G4int entries() const
Definition: G4THitsMap.hh:79
const G4String & GetUnit() const
G4double GetUnitValue() const

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 76 of file G4PSPassageCellFlux.cc.

77{
78
79 if ( IsPassed(aStep) ) {
81 (aStep->GetPreStepPoint()->GetTouchable()))
82 ->GetReplicaNumber(indexDepth);
83 G4double cubicVolume = ComputeVolume(aStep, idx);
84
85 fCellFlux /= cubicVolume;
86 G4int index = GetIndex(aStep);
87 EvtMap->add(index,fCellFlux);
88 }
89
90 return TRUE;
91}
virtual G4double ComputeVolume(G4Step *, G4int idx)
virtual G4bool IsPassed(G4Step *)
const G4VTouchable * GetTouchable() const
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:138
virtual G4int GetIndex(G4Step *)

◆ SetUnit()

void G4PSPassageCellFlux::SetUnit ( const G4String unit)
virtual

Definition at line 158 of file G4PSPassageCellFlux.cc.

159{
160 CheckAndSetUnit(unit,"Per Unit Surface");
161}
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 62 of file G4PSPassageCellFlux.hh.

62{ weighted = flg; }

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