Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PSVolumeFlux.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Class G4PSVolumeFlux
27//
28// Class description:
29// Scorer that scores number of tracks coming into the associated volume.
30// Optionally number of tracks can be divided by the surface area to
31// score the volume current and divided by cos(theta) for volume flux
32// where theta is the incident angle
33//
34// - Created M. Asai, Sept. 2020
35//
36//
37#include "G4PSVolumeFlux.hh"
38
39#include "G4SystemOfUnits.hh"
40#include "G4StepStatus.hh"
41#include "G4Track.hh"
42#include "G4VSolid.hh"
43#include "G4VPhysicalVolume.hh"
45#include "G4UnitsTable.hh"
47#include "G4VScoreHistFiller.hh"
48
50 : G4VPrimitivePlotter(name, depth)
51 , fDirection(direction)
52{}
53
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}
142
144{
145 if(HCID < 0)
146 HCID = GetCollectionID(0);
148 GetName());
149 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap);
150}
151
152void G4PSVolumeFlux::clear() { EvtMap->clear(); }
153
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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double cosTheta() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4VSolid * GetSolid() const
const G4AffineTransform & GetTopTransform() const
void clear() override
void PrintAll() override
G4PSVolumeFlux(G4String name, G4int direction=1, G4int depth=0)
void Initialize(G4HCofThisEvent *) override
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
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
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
G4MultiFunctionalDetector * GetMultiFunctionalDetector() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
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
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t entries() const
Definition: G4THitsMap.hh:158
void clear()
Definition: G4THitsMap.hh:524
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:82