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

#include <G4PSSphereSurfaceFlux.hh>

+ Inheritance diagram for G4PSSphereSurfaceFlux:

Public Member Functions

 G4PSSphereSurfaceFlux (G4String name, G4int direction, G4int depth=0)
 
 G4PSSphereSurfaceFlux (G4String name, G4int direction, const G4String &unit, G4int depth=0)
 
 ~G4PSSphereSurfaceFlux () override=default
 
void Weighted (G4bool flg=true)
 
void DivideByArea (G4bool flg=true)
 
void Initialize (G4HCofThisEvent *) override
 
void clear () override
 
void PrintAll () override
 
virtual void SetUnit (const G4String &unit)
 
- 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
 
G4int IsSelectedSurface (G4Step *, G4Sphere *)
 
virtual void DefineUnitAndCategory ()
 
- 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 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 59 of file G4PSSphereSurfaceFlux.hh.

Constructor & Destructor Documentation

◆ G4PSSphereSurfaceFlux() [1/2]

G4PSSphereSurfaceFlux::G4PSSphereSurfaceFlux ( G4String name,
G4int direction,
G4int depth = 0 )

Definition at line 59 of file G4PSSphereSurfaceFlux.cc.

61 : G4PSSphereSurfaceFlux(name, direction, "percm2", depth)
62{}
G4PSSphereSurfaceFlux(G4String name, G4int direction, G4int depth=0)

◆ G4PSSphereSurfaceFlux() [2/2]

G4PSSphereSurfaceFlux::G4PSSphereSurfaceFlux ( G4String name,
G4int direction,
const G4String & unit,
G4int depth = 0 )

Definition at line 64 of file G4PSSphereSurfaceFlux.cc.

66 : G4VPrimitiveScorer(name, depth)
67 , HCID(-1)
68 , fDirection(direction)
69 , EvtMap(nullptr)
70 , weighted(true)
71 , divideByArea(true)
72{
74 SetUnit(unit);
75}
virtual void SetUnit(const G4String &unit)
G4VPrimitiveScorer(G4String name, G4int depth=0)

◆ ~G4PSSphereSurfaceFlux()

G4PSSphereSurfaceFlux::~G4PSSphereSurfaceFlux ( )
overridedefault

Member Function Documentation

◆ clear()

void G4PSSphereSurfaceFlux::clear ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 221 of file G4PSSphereSurfaceFlux.cc.

221{ EvtMap->clear(); }

◆ DefineUnitAndCategory()

void G4PSSphereSurfaceFlux::DefineUnitAndCategory ( )
protectedvirtual

Definition at line 259 of file G4PSSphereSurfaceFlux.cc.

260{
261 // Per Unit Surface
262 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
263 (1. / cm2));
264 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
265 (1. / mm2));
266 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
267}

Referenced by G4PSSphereSurfaceFlux().

◆ DivideByArea()

void G4PSSphereSurfaceFlux::DivideByArea ( G4bool flg = true)
inline

Definition at line 70 of file G4PSSphereSurfaceFlux.hh.

70{ divideByArea = flg; }

◆ Initialize()

void G4PSSphereSurfaceFlux::Initialize ( G4HCofThisEvent * HCE)
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 213 of file G4PSSphereSurfaceFlux.cc.

214{
215 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
216 if(HCID < 0)
217 HCID = GetCollectionID(0);
218 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap);
219}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector

◆ IsSelectedSurface()

G4int G4PSSphereSurfaceFlux::IsSelectedSurface ( G4Step * aStep,
G4Sphere * sphereSolid )
protected

Definition at line 158 of file G4PSSphereSurfaceFlux.cc.

160{
161 G4TouchableHandle theTouchable =
165
167 {
168 // Entering Geometry
169 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition();
170 G4ThreeVector localpos1 =
171 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
172 G4double localR2 = localpos1.x() * localpos1.x() +
173 localpos1.y() * localpos1.y() +
174 localpos1.z() * localpos1.z();
175 // G4double InsideRadius2 =
176 // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
177 // if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
178 G4double InsideRadius = sphereSolid->GetInnerRadius();
179 if(localR2 >
180 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
181 localR2 <
182 (InsideRadius + kCarTolerance) * (InsideRadius + kCarTolerance))
183 {
184 return fFlux_In;
185 }
186 }
187
189 {
190 // Exiting Geometry
191 G4ThreeVector stppos2 = aStep->GetPostStepPoint()->GetPosition();
192 G4ThreeVector localpos2 =
193 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
194 G4double localR2 = localpos2.x() * localpos2.x() +
195 localpos2.y() * localpos2.y() +
196 localpos2.z() * localpos2.z();
197 // G4double InsideRadius2 =
198 // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
199 // if(std::facb(localR2 - InsideRadius2) ) < kCarTolerance ){
200 G4double InsideRadius = sphereSolid->GetInnerRadius();
201 if(localR2 >
202 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
203 localR2 <
204 (InsideRadius + kCarTolerance) * (InsideRadius + kCarTolerance))
205 {
206 return fFlux_Out;
207 }
208 }
209
210 return -1;
211}
const G4double kCarTolerance
@ fFlux_Out
@ fFlux_In
@ fGeomBoundary
double G4double
Definition G4Types.hh:83
double z() const
double x() const
double y() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
const G4AffineTransform & GetTopTransform() const
G4double GetInnerRadius() const
G4StepStatus GetStepStatus() const
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
virtual const G4NavigationHistory * GetHistory() const

Referenced by ProcessHits().

◆ PrintAll()

void G4PSSphereSurfaceFlux::PrintAll ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 223 of file G4PSSphereSurfaceFlux.cc.

224{
225 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
226 G4cout << " PrimitiveScorer " << GetName() << G4endl;
227 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
228 for(const auto& [copy, flux] : *(EvtMap->GetMap()))
229 {
230 G4cout << " copy no.: " << copy
231 << " Flux : " << *(flux) / GetUnitValue() << " ["
232 << GetUnit() << "]" << G4endl;
233 }
234}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetUnit() const
G4double GetUnitValue() const
Map_t * GetMap() const
size_t entries() const
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 77 of file G4PSSphereSurfaceFlux.cc.

78{
79 G4StepPoint* preStep = aStep->GetPreStepPoint();
80
81 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
82 G4VPVParameterisation* physParam = physVol->GetParameterisation();
83 G4VSolid* solid = nullptr;
84 if(physParam != nullptr)
85 { // for parameterized volume
86 G4int idx =
88 ->GetReplicaNumber(indexDepth);
89 solid = physParam->ComputeSolid(idx, physVol);
90 solid->ComputeDimensions(physParam, idx, physVol);
91 }
92 else
93 { // for ordinary volume
94 solid = physVol->GetLogicalVolume()->GetSolid();
95 }
96
97 auto sphereSolid = (G4Sphere*) (solid);
98
99 G4int dirFlag = IsSelectedSurface(aStep, sphereSolid);
100 if(dirFlag > 0)
101 {
102 if(fDirection == fFlux_InOut || fDirection == dirFlag)
103 {
104 G4StepPoint* thisStep = nullptr;
105 if(dirFlag == fFlux_In)
106 {
107 thisStep = preStep;
108 }
109 else if(dirFlag == fFlux_Out)
110 {
111 thisStep = aStep->GetPostStepPoint();
112 }
113 else
114 {
115 return false;
116 }
117
118 G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
119 G4ThreeVector pdirection = thisStep->GetMomentumDirection();
120 G4ThreeVector localdir =
121 theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
122 G4double localdirL2 = localdir.x() * localdir.x() +
123 localdir.y() * localdir.y() +
124 localdir.z() * localdir.z();
125 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition();
126 G4ThreeVector localpos1 =
127 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
128 G4double localR2 = localpos1.x() * localpos1.x() +
129 localpos1.y() * localpos1.y() +
130 localpos1.z() * localpos1.z();
131 G4double anglefactor =
132 (localdir.x() * localpos1.x() + localdir.y() * localpos1.y() +
133 localdir.z() * localpos1.z()) /
134 std::sqrt(localdirL2) / std::sqrt(localR2);
135 if(anglefactor < 0.0)
136 anglefactor *= -1.0;
137
138 G4double current = 1.0 / anglefactor;
139 if(weighted)
140 current *= thisStep->GetWeight(); // Flux (Particle Weight)
141 if(divideByArea) // Flux with angle.
142 {
143 G4double radi = sphereSolid->GetInnerRadius();
144 G4double dph = sphereSolid->GetDeltaPhiAngle() / radian;
145 G4double stth = sphereSolid->GetStartThetaAngle() / radian;
146 G4double enth = stth + sphereSolid->GetDeltaThetaAngle() / radian;
147 current /= radi * radi * dph * (-std::cos(enth) + std::cos(stth));
148 }
149
150 G4int index = GetIndex(aStep);
151 EvtMap->add(index, current);
152 }
153 }
154
155 return true;
156}
@ fFlux_InOut
int G4int
Definition G4Types.hh:85
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4VSolid * GetSolid() const
G4int IsSelectedSurface(G4Step *, G4Sphere *)
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetMomentumDirection() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetWeight() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4int GetIndex(G4Step *)
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition G4VSolid.cc:137
size_t add(const G4int &key, U *&aHit) const

◆ SetUnit()

void G4PSSphereSurfaceFlux::SetUnit ( const G4String & unit)
virtual

Definition at line 236 of file G4PSSphereSurfaceFlux.cc.

237{
238 if(divideByArea)
239 {
240 CheckAndSetUnit(unit, "Per Unit Surface");
241 }
242 else
243 {
244 if(unit.empty())
245 {
246 unitName = unit;
247 unitValue = 1.0;
248 }
249 else
250 {
251 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
252 GetUnit() + "] ) for " + GetName();
253 G4Exception("G4PSSphereSurfaceFlux::SetUnit", "DetPS0016", JustWarning,
254 msg);
255 }
256 }
257}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void CheckAndSetUnit(const G4String &unit, const G4String &category)

Referenced by G4PSSphereSurfaceFlux(), and G4PSSphereSurfaceFlux3D::G4PSSphereSurfaceFlux3D().

◆ Weighted()

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

Definition at line 67 of file G4PSSphereSurfaceFlux.hh.

67{ weighted = flg; }

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