BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes
BDSPSCellFluxScaled3D Class Reference

Primitive scorer for cell flux in a 3D mesh with a conversion factor. More...

#include <BDSPSCellFluxScaled3D.hh>

Inheritance diagram for BDSPSCellFluxScaled3D:
Inheritance graph
Collaboration diagram for BDSPSCellFluxScaled3D:
Collaboration graph

Public Member Functions

 BDSPSCellFluxScaled3D (const G4String &scorerName, const BDSHistBinMapper *mapperIn, const G4String &unitIn="percm2", G4int ni=1, G4int nj=1, G4int nk=1, G4int depi=2, G4int depj=1, G4int depk=0)
 
 BDSPSCellFluxScaled3D (const G4String &scorerName, const BDSHistBinMapper *mapperIn, const G4String &filename, const G4String &unitIn="percm2", G4int ni=1, G4int nj=1, G4int nk=1, G4int depi=2, G4int depj=1, G4int depk=0)
 
void Initialize (G4HCofThisEvent *HCE) override
 
void EndOfEvent (G4HCofThisEvent *HCE) override
 
void clear () override
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *) override
 
G4int GetIndex (G4Step *aStep) override
 
virtual G4double GetConversionFactor (G4int particleID, G4double kineticEnergy) const
 

Private Member Functions

void DefineUnitAndCategory () const
 Define units -> taken from G4PSCellFlux.
 

Private Attributes

G4int HCID3D
 Collection ID.
 
G4THitsMap< G4double > * evtMap3D
 Hits map.
 
G4PhysicsVector * conversionFactor
 Conversion factor interpolator object.
 
const BDSHistBinMappermapper
 Mapping from coordinate systems in mesh to global replica number.
 
std::map< G4VSolid *, G4double > volumeCache
 
G4int fDepthi
 Depth in replica to look for each dimension.
 
G4int fDepthj
 Depth in replica to look for each dimension.
 
G4int fDepthk
 Depth in replica to look for each dimension.
 

Detailed Description

Primitive scorer for cell flux in a 3D mesh with a conversion factor.

The cell flux (step length / cell volume) * weight is also multiplied by a conversion factor (just a number) as listed vs energy in a supplied file. The default is none and just a factor of 1.

The implementation also differs from G4PSCellFlux3D as we cache the volume to avoid repeated calculation.

Author
Robin Tesse

Definition at line 42 of file BDSPSCellFluxScaled3D.hh.

Constructor & Destructor Documentation

◆ BDSPSCellFluxScaled3D() [1/2]

BDSPSCellFluxScaled3D::BDSPSCellFluxScaled3D ( const G4String &  scorerName,
const BDSHistBinMapper mapperIn,
const G4String &  unitIn = "percm2",
G4int  ni = 1,
G4int  nj = 1,
G4int  nk = 1,
G4int  depi = 2,
G4int  depj = 1,
G4int  depk = 0 
)

Constructor where no conversion factor file is provided and all cell fluxes just use conversion factor 1.0.

Definition at line 42 of file BDSPSCellFluxScaled3D.cc.

References DefineUnitAndCategory().

Here is the call graph for this function:

◆ BDSPSCellFluxScaled3D() [2/2]

BDSPSCellFluxScaled3D::BDSPSCellFluxScaled3D ( const G4String &  scorerName,
const BDSHistBinMapper mapperIn,
const G4String &  filename,
const G4String &  unitIn = "percm2",
G4int  ni = 1,
G4int  nj = 1,
G4int  nk = 1,
G4int  depi = 2,
G4int  depj = 1,
G4int  depk = 0 
)

Constructor where conversion factor file is provided and loaded into a physics vector. Cell fluxes are multiplied by the factor as a function of the particle kinetic energy.

Definition at line 65 of file BDSPSCellFluxScaled3D.cc.

References conversionFactor, BDS::GetFullPath(), and BDSScorerConversionLoader< T >::Load().

Here is the call graph for this function:

◆ ~BDSPSCellFluxScaled3D()

BDSPSCellFluxScaled3D::~BDSPSCellFluxScaled3D ( )
overridevirtual

Definition at line 97 of file BDSPSCellFluxScaled3D.cc.

Member Function Documentation

◆ clear()

void BDSPSCellFluxScaled3D::clear ( )
override

Definition at line 153 of file BDSPSCellFluxScaled3D.cc.

◆ DefineUnitAndCategory()

void BDSPSCellFluxScaled3D::DefineUnitAndCategory ( ) const
private

Define units -> taken from G4PSCellFlux.

Definition at line 182 of file BDSPSCellFluxScaled3D.cc.

Referenced by BDSPSCellFluxScaled3D().

Here is the caller graph for this function:

◆ EndOfEvent()

void BDSPSCellFluxScaled3D::EndOfEvent ( G4HCofThisEvent *  HCE)
override

Definition at line 150 of file BDSPSCellFluxScaled3D.cc.

◆ GetConversionFactor()

G4double BDSPSCellFluxScaled3D::GetConversionFactor ( G4int  particleID,
G4double  kineticEnergy 
) const
virtual

Definition at line 135 of file BDSPSCellFluxScaled3D.cc.

◆ GetIndex()

G4int BDSPSCellFluxScaled3D::GetIndex ( G4Step *  aStep)
override

Definition at line 158 of file BDSPSCellFluxScaled3D.cc.

◆ Initialize()

void BDSPSCellFluxScaled3D::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 141 of file BDSPSCellFluxScaled3D.cc.

◆ ProcessHits()

G4bool BDSPSCellFluxScaled3D::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
override

Definition at line 102 of file BDSPSCellFluxScaled3D.cc.

Field Documentation

◆ conversionFactor

G4PhysicsVector* BDSPSCellFluxScaled3D::conversionFactor
private

Conversion factor interpolator object.

Definition at line 88 of file BDSPSCellFluxScaled3D.hh.

Referenced by BDSPSCellFluxScaled3D().

◆ evtMap3D

G4THitsMap<G4double>* BDSPSCellFluxScaled3D::evtMap3D
private

Hits map.

Definition at line 79 of file BDSPSCellFluxScaled3D.hh.

◆ fDepthi

G4int BDSPSCellFluxScaled3D::fDepthi
private

Depth in replica to look for each dimension.

Definition at line 82 of file BDSPSCellFluxScaled3D.hh.

◆ fDepthj

G4int BDSPSCellFluxScaled3D::fDepthj
private

Depth in replica to look for each dimension.

Definition at line 83 of file BDSPSCellFluxScaled3D.hh.

◆ fDepthk

G4int BDSPSCellFluxScaled3D::fDepthk
private

Depth in replica to look for each dimension.

Definition at line 84 of file BDSPSCellFluxScaled3D.hh.

◆ HCID3D

G4int BDSPSCellFluxScaled3D::HCID3D
private

Collection ID.

Definition at line 78 of file BDSPSCellFluxScaled3D.hh.

◆ mapper

const BDSHistBinMapper* BDSPSCellFluxScaled3D::mapper
private

Mapping from coordinate systems in mesh to global replica number.

We don't own this.

Definition at line 91 of file BDSPSCellFluxScaled3D.hh.

◆ volumeCache

std::map<G4VSolid*, G4double> BDSPSCellFluxScaled3D::volumeCache
private

Definition at line 93 of file BDSPSCellFluxScaled3D.hh.


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