Geant4 9.6.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 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)
 

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 58 of file G4PSCellFlux.hh.

Constructor & Destructor Documentation

◆ G4PSCellFlux() [1/2]

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

Definition at line 56 of file G4PSCellFlux.cc.

57 :G4VPrimitiveScorer(name,depth),HCID(-1),weighted(true)
58{
60 SetUnit("percm2");
61 //verboseLevel = 10;
62}
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 64 of file G4PSCellFlux.cc.

65 :G4VPrimitiveScorer(name,depth),HCID(-1),weighted(true)
66{
68 SetUnit(unit);
69}

◆ ~G4PSCellFlux()

G4PSCellFlux::~G4PSCellFlux ( )
virtual

Definition at line 71 of file G4PSCellFlux.cc.

72{;}

Member Function Documentation

◆ clear()

void G4PSCellFlux::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 103 of file G4PSCellFlux.cc.

103 {
104 EvtMap->clear();
105}
void clear()
Definition: G4THitsMap.hh:209

◆ ComputeVolume()

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

Reimplemented in G4PSCellFluxForCylinder3D.

Definition at line 136 of file G4PSCellFlux.cc.

136 {
137
139 G4VPVParameterisation* physParam = physVol->GetParameterisation();
140 G4VSolid* solid = 0;
141 if(physParam)
142 { // for parameterized volume
143 if(idx<0)
144 {
146 ED << "Incorrect replica number --- GetReplicaNumber : " << idx << G4endl;
147 G4Exception("G4PSCellFlux::ComputeVolume","DetPS0001",JustWarning,ED);
148 }
149 solid = physParam->ComputeSolid(idx, physVol);
150 solid->ComputeDimensions(physParam,idx,physVol);
151 }
152 else
153 { // for ordinary volume
154 solid = physVol->GetLogicalVolume()->GetSolid();
155 }
156
157 return solid->GetCubicVolume();
158}
@ 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 G4PSCellFlux::DefineUnitAndCategory ( )
protectedvirtual

Definition at line 129 of file G4PSCellFlux.cc.

129 {
130 // Per Unit Surface
131 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
132 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
133 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
134}

Referenced by G4PSCellFlux().

◆ DrawAll()

void G4PSCellFlux::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 107 of file G4PSCellFlux.cc.

108{;}

◆ EndOfEvent()

void G4PSCellFlux::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 100 of file G4PSCellFlux.cc.

101{;}

◆ Initialize()

void G4PSCellFlux::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 92 of file G4PSCellFlux.cc.

93{
95 GetName());
96 if ( HCID < 0 ) HCID = GetCollectionID(0);
97 HCE->AddHitsCollection(HCID,EvtMap);
98}
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 110 of file G4PSCellFlux.cc.

111{
112 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
113 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
114 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
115 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
116 for(; itr != EvtMap->GetMap()->end(); itr++) {
117 G4cout << " copy no.: " << itr->first
118 << " cell flux : " << *(itr->second)/GetUnitValue()
119 << " [" << GetUnit() << "]"
120 << G4endl;
121 }
122}
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 G4PSCellFlux::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 74 of file G4PSCellFlux.cc.

75{
76 G4double stepLength = aStep->GetStepLength();
77 if ( stepLength == 0. ) return FALSE;
78
80 (aStep->GetPreStepPoint()->GetTouchable()))
81 ->GetReplicaNumber(indexDepth);
82 G4double cubicVolume = ComputeVolume(aStep, idx);
83
84 G4double CellFlux = stepLength / cubicVolume;
85 if (weighted) CellFlux *= aStep->GetPreStepPoint()->GetWeight();
86 G4int index = GetIndex(aStep);
87 EvtMap->add(index,CellFlux);
88
89 return TRUE;
90}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
virtual G4double ComputeVolume(G4Step *, G4int idx)
const G4VTouchable * GetTouchable() const
G4double GetWeight() const
G4double GetStepLength() const
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:138
virtual G4int GetIndex(G4Step *)
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52

◆ SetUnit()

void G4PSCellFlux::SetUnit ( const G4String unit)
virtual

Definition at line 124 of file G4PSCellFlux.cc.

125{
126 CheckAndSetUnit(unit,"Per Unit Surface");
127}
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: