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

#include <G4PSVolumeFlux.hh>

+ Inheritance diagram for G4PSVolumeFlux:

Public Member Functions

 G4PSVolumeFlux (G4String name, G4int direction=1, G4int depth=0)
 
 ~G4PSVolumeFlux () override=default
 
void Initialize (G4HCofThisEvent *) override
 
void clear () override
 
void PrintAll () override
 
void SetDivAre (G4bool val)
 
void SetDivCos (G4bool val)
 
- 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
 
- 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 43 of file G4PSVolumeFlux.hh.

Constructor & Destructor Documentation

◆ G4PSVolumeFlux()

G4PSVolumeFlux::G4PSVolumeFlux ( G4String name,
G4int direction = 1,
G4int depth = 0 )

Definition at line 49 of file G4PSVolumeFlux.cc.

50 : G4VPrimitivePlotter(name, depth)
51 , fDirection(direction)
52{}

◆ ~G4PSVolumeFlux()

G4PSVolumeFlux::~G4PSVolumeFlux ( )
overridedefault

Member Function Documentation

◆ clear()

void G4PSVolumeFlux::clear ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 152 of file G4PSVolumeFlux.cc.

152{ EvtMap->clear(); }

◆ Initialize()

void G4PSVolumeFlux::Initialize ( G4HCofThisEvent * HCE)
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 143 of file G4PSVolumeFlux.cc.

144{
145 if(HCID < 0)
146 HCID = GetCollectionID(0);
148 GetName());
149 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap);
150}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * GetMultiFunctionalDetector() const

◆ PrintAll()

void G4PSVolumeFlux::PrintAll ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 154 of file G4PSVolumeFlux.cc.

155{
156 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
157 G4cout << " PrimitiveScorer" << GetName() << G4endl;
158 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
159 for(const auto& [copy, flux] : *(EvtMap->GetMap()))
160 {
161 G4cout << " copy no.: " << copy << " flux : " << *(flux)
162 << G4endl;
163 }
164}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4MultiFunctionalDetector * detector
Map_t * GetMap() const
size_t entries() const
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 54 of file G4PSVolumeFlux.cc.

55{
56 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
57 G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
58 G4StepPoint* thisStepPoint = nullptr;
59 if(fDirection == 1) // Score only the inward particle
60 {
61 if(preStepPoint->GetStepStatus() == fGeomBoundary)
62 {
63 thisStepPoint = preStepPoint;
64 }
65 else
66 {
67 return false;
68 }
69 }
70 else if(fDirection == 2) // Score only the outward particle
71 {
72 if(postStepPoint->GetStepStatus() == fGeomBoundary)
73 {
74 thisStepPoint = postStepPoint;
75 }
76 else
77 {
78 return false;
79 }
80 }
81
82 G4double flux = preStepPoint->GetWeight();
83
84 if(divare || divcos)
85 {
86 G4VPhysicalVolume* physVol = preStepPoint->GetPhysicalVolume();
87 G4VPVParameterisation* physParam = physVol->GetParameterisation();
88 G4VSolid* solid = nullptr;
89 if(physParam != nullptr)
90 { // for parameterized volume
91 auto idx = ((G4TouchableHistory*) (preStepPoint->GetTouchable()))
92 ->GetReplicaNumber(indexDepth);
93 solid = physParam->ComputeSolid(idx, physVol);
94 solid->ComputeDimensions(physParam, idx, physVol);
95 }
96 else
97 { // for ordinary volume
98 solid = physVol->GetLogicalVolume()->GetSolid();
99 }
100
101 if(divare)
102 {
103 flux /= solid->GetSurfaceArea();
104 }
105
106 if(divcos)
107 {
108 G4TouchableHandle theTouchable = thisStepPoint->GetTouchableHandle();
109 G4ThreeVector pdirection = thisStepPoint->GetMomentumDirection();
110 G4ThreeVector localdir =
111 theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
112 G4ThreeVector globalPos = thisStepPoint->GetPosition();
113 G4ThreeVector localPos =
114 theTouchable->GetHistory()->GetTopTransform().TransformPoint(globalPos);
115 G4ThreeVector surfNormal = solid->SurfaceNormal(localPos);
116 G4double cosT = surfNormal.cosTheta(localdir);
117 if(cosT != 0.)
118 flux /= std::abs(cosT);
119 }
120 }
121
122 G4int index = GetIndex(aStep);
123 EvtMap->add(index, flux);
124
125 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend())
126 {
127 auto filler = G4VScoreHistFiller::Instance();
128 if(filler == nullptr)
129 {
131 "G4PSVolumeFlux::ProcessHits", "SCORER0123", JustWarning,
132 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
133 }
134 else
135 {
136 filler->FillH1(hitIDMap[index], thisStepPoint->GetKineticEnergy(), flux);
137 }
138 }
139
140 return true;
141}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ fGeomBoundary
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
double cosTheta() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4VSolid * GetSolid() const
const G4AffineTransform & GetTopTransform() const
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
const G4TouchableHandle & GetTouchableHandle() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetKineticEnergy() const
G4double GetWeight() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
virtual const G4NavigationHistory * GetHistory() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition G4VSolid.cc:137
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual G4double GetSurfaceArea()
Definition G4VSolid.cc:250
size_t add(const G4int &key, U *&aHit) const

◆ SetDivAre()

void G4PSVolumeFlux::SetDivAre ( G4bool val)
inline

Definition at line 54 of file G4PSVolumeFlux.hh.

54{ divare = val; }

◆ SetDivCos()

void G4PSVolumeFlux::SetDivCos ( G4bool val)
inline

Definition at line 55 of file G4PSVolumeFlux.hh.

55{ divcos = val; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


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