BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSPSCellFluxScaled3D.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef BDSPSCELLFLUXSCALED3D_H
20#define BDSPSCELLFLUXSCALED3D_H
21
22#include "globals.hh"
23#include "G4THitsMap.hh"
24#include "G4VPrimitiveScorer.hh"
25
27class G4PhysicsVector;
28
42class BDSPSCellFluxScaled3D: public G4VPrimitiveScorer
43{
44public:
47 BDSPSCellFluxScaled3D(const G4String& scorerName,
48 const BDSHistBinMapper* mapperIn,
49 const G4String& unitIn = "percm2",
50 G4int ni=1, G4int nj=1, G4int nk=1,
51 G4int depi=2, G4int depj=1, G4int depk=0);
52
55 BDSPSCellFluxScaled3D(const G4String& scorerName,
56 const BDSHistBinMapper* mapperIn,
57 const G4String& filename,
58 const G4String& unitIn = "percm2",
59 G4int ni=1, G4int nj=1, G4int nk=1,
60 G4int depi=2, G4int depj=1, G4int depk=0);
61
62 virtual ~BDSPSCellFluxScaled3D() override;
63
64public:
65 void Initialize(G4HCofThisEvent* HCE) override;
66 void EndOfEvent(G4HCofThisEvent* HCE) override;
67 void clear() override;
68 G4bool ProcessHits(G4Step* aStep, G4TouchableHistory*) override;
69 G4int GetIndex(G4Step* aStep) override;
70
71 virtual G4double GetConversionFactor(G4int particleID,
72 G4double kineticEnergy) const;
73
74private:
76 void DefineUnitAndCategory() const;
77
78 G4int HCID3D;
80
82 G4int fDepthi;
83 G4int fDepthj;
84 G4int fDepthk;
86
88 G4PhysicsVector* conversionFactor;
89
92
93 std::map<G4VSolid*, G4double> volumeCache;
94};
95
96#endif
Mapping from axis indices to 1D index.
Primitive scorer for cell flux in a 3D mesh with a conversion factor.
G4THitsMap< G4double > * evtMap3D
Hits map.
const BDSHistBinMapper * mapper
Mapping from coordinate systems in mesh to global replica number.
G4int fDepthk
Depth in replica to look for each dimension.
G4PhysicsVector * conversionFactor
Conversion factor interpolator object.
G4int fDepthi
Depth in replica to look for each dimension.
G4int fDepthj
Depth in replica to look for each dimension.
void DefineUnitAndCategory() const
Define units -> taken from G4PSCellFlux.
G4int HCID3D
Collection ID.