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

#include <G4PSCellFlux.hh>

+ Inheritance diagram for G4PSCellFlux:

Public Member Functions

 G4PSCellFlux (G4String name, G4int depth=0)
 
 G4PSCellFlux (G4String name, const G4String &unit, G4int depth=0)
 
virtual ~G4PSCellFlux ()
 
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 G4VPrimitivePlotter
 G4VPrimitivePlotter (G4String name, G4int depth=0)
 
virtual ~G4VPrimitivePlotter ()
 
void Plot (G4int copyNo, G4int histID)
 
G4int GetNumberOfHist () const
 
- 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 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)
 
G4VSolidComputeSolid (G4Step *aStep, G4int replicaIdx)
 
G4VSolidComputeCurrentSolid (G4Step *aStep)
 

Additional Inherited Members

- Protected Attributes inherited from G4VPrimitivePlotter
std::map< G4int, G4inthitIDMap
 
- 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 58 of file G4PSCellFlux.hh.

Constructor & Destructor Documentation

◆ G4PSCellFlux() [1/2]

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

Definition at line 58 of file G4PSCellFlux.cc.

59 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0),weighted(true)
60{
62 SetUnit("percm2");
63 //verboseLevel = 10;
64}
virtual void SetUnit(const G4String &unit)
virtual void DefineUnitAndCategory()

◆ G4PSCellFlux() [2/2]

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

Definition at line 66 of file G4PSCellFlux.cc.

67 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0),weighted(true)
68{
70 SetUnit(unit);
71}

◆ ~G4PSCellFlux()

G4PSCellFlux::~G4PSCellFlux ( )
virtual

Definition at line 73 of file G4PSCellFlux.cc.

74{;}

Member Function Documentation

◆ clear()

void G4PSCellFlux::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 119 of file G4PSCellFlux.cc.

119 {
120 EvtMap->clear();
121}
void clear()
Definition: G4THitsMap.hh:537

◆ ComputeVolume()

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

Reimplemented in G4PSCellFluxForCylinder3D.

Definition at line 152 of file G4PSCellFlux.cc.

153{
154 G4VSolid* solid = ComputeSolid( aStep, idx );
155 assert(solid);
156 return solid->GetCubicVolume();
157}
G4VSolid * ComputeSolid(G4Step *aStep, G4int replicaIdx)
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:176

Referenced by ProcessHits().

◆ DefineUnitAndCategory()

void G4PSCellFlux::DefineUnitAndCategory ( )
protectedvirtual

Definition at line 145 of file G4PSCellFlux.cc.

145 {
146 // Per Unit Surface
147 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
148 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
149 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
150}

Referenced by G4PSCellFlux().

◆ DrawAll()

void G4PSCellFlux::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 123 of file G4PSCellFlux.cc.

124{;}

◆ EndOfEvent()

void G4PSCellFlux::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 116 of file G4PSCellFlux.cc.

117{;}

◆ Initialize()

void G4PSCellFlux::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 108 of file G4PSCellFlux.cc.

109{
110 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
111 GetName());
112 if ( HCID < 0 ) HCID = GetCollectionID(0);
113 HCE->AddHitsCollection(HCID,EvtMap);
114}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ PrintAll()

void G4PSCellFlux::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 126 of file G4PSCellFlux.cc.

127{
128 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
129 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
130 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
131 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
132 for(; itr != EvtMap->GetMap()->end(); itr++) {
133 G4cout << " copy no.: " << itr->first
134 << " cell flux : " << *(itr->second)/GetUnitValue()
135 << " [" << GetUnit() << "]"
136 << G4endl;
137 }
138}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetUnit() const
G4double GetUnitValue() const
Map_t * GetMap() const
Definition: G4THitsMap.hh:151
size_t entries() const
Definition: G4THitsMap.hh:155

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 76 of file G4PSCellFlux.cc.

77{
78 G4double stepLength = aStep->GetStepLength();
79 if ( stepLength == 0. ) return FALSE;
80
82 (aStep->GetPreStepPoint()->GetTouchable()))
83 ->GetReplicaNumber(indexDepth);
84 G4double cubicVolume = ComputeVolume(aStep, idx);
85
86 G4double CellFlux = stepLength / cubicVolume;
87 if (weighted) CellFlux *= aStep->GetPreStepPoint()->GetWeight();
88 G4int index = GetIndex(aStep);
89 EvtMap->add(index,CellFlux);
90
91 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
92 {
93 auto filler = G4VScoreHistFiller::Instance();
94 if(!filler)
95 {
96 G4Exception("G4PSCellFlux::ProcessHits","SCORER0123",JustWarning,
97 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
98 }
99 else
100 {
101 filler->FillH1(hitIDMap[index],aStep->GetPreStepPoint()->GetKineticEnergy(),CellFlux);
102 }
103 }
104
105 return TRUE;
106}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
virtual G4double ComputeVolume(G4Step *, G4int idx)
const G4VTouchable * GetTouchable() const
G4double GetKineticEnergy() const
G4double GetWeight() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

◆ SetUnit()

void G4PSCellFlux::SetUnit ( const G4String unit)
virtual

Definition at line 140 of file G4PSCellFlux.cc.

141{
142 CheckAndSetUnit(unit,"Per Unit Surface");
143}
void CheckAndSetUnit(const G4String &unit, const G4String &category)

Referenced by G4PSCellFlux(), G4PSCellFlux3D::G4PSCellFlux3D(), and G4ScoreQuantityMessenger::SetNewValue().

◆ Weighted()

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

Definition at line 65 of file G4PSCellFlux.hh.

65{ weighted = flg; }

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