BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSParallelWorldImportance.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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 BDSPARALLELWORLDIMPORTANCE_H
20#define BDSPARALLELWORLDIMPORTANCE_H
21
22#include "BDSExtent.hh"
23#include "BDSImportanceVolumeStore.hh"
24
25#include "globals.hh" // geant4 types / globals
26#include "G4GeometryCell.hh"
27#include "G4VUserParallelWorld.hh"
28
29#include <map>
30
31class G4UserLimits;
32class G4VisAttributes;
33class G4VPhysicalVolume;
34
35namespace GMAD {
36 struct Element;
37 template<typename T> class FastList;
38 class Placement;
39}
40
47class BDSParallelWorldImportance: public G4VUserParallelWorld
48{
49public:
50 BDSParallelWorldImportance(G4String name,
51 G4String importanceWorldGeometryFile,
52 G4String importanceValuesFile);
54
57 void Construct();
58
61
64
66 G4GeometryCell GetGeometryCell(G4int i) const;
67
69 void AddIStore();
70
71 virtual void ConstructSD();
72
74 inline G4VPhysicalVolume* GetWorldVolume() {return imWorldPV;}
75
76private:
80
83 void BuildWorld();
84
86 G4VPhysicalVolume* imWorldPV;
87
89
92
94 std::map<G4String, G4double> imVolumesAndValues;
95
96 G4String imGeomFile;
97 G4String imVolMap;
98 const G4String componentName;
99
102 G4UserLimits* userLimits;
103 G4VisAttributes* visAttr;
105
107 G4double GetCellImportanceValue(const G4String& cellName);
108};
109
110#endif
111
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
Registry of importance values.
Class that constructs a parallel importance world.
const G4String componentName
String preprended to geometry with preprocessGDML.
BDSImportanceVolumeStore imVolumeStore
Store of geometry cells for importance volumes.
G4double GetCellImportanceValue(const G4String &cellName)
Get importance value of a given physical volume name.
BDSParallelWorldImportance & operator=(const BDSParallelWorldImportance &)
assignment and copy constructor not implemented nor used
void AddIStore()
Create IStore for all importance sampling geometry cells.
std::map< G4String, G4double > imVolumesAndValues
Container for all user placed physical volumes and corresponding importance values.
G4VPhysicalVolume * GetWorldVolume()
World volume getter required in parallel world utilities.
BDSExtent worldExtent
Record of the world extent.
G4VPhysicalVolume * imWorldPV
Importance sampling world volume.
G4int verbosity
Cached global constants values.
BDSExtent WorldExtent() const
Public access to the world extent.
G4GeometryCell GetGeometryCell(G4int i) const
Get geometry cell from store.
void BuildPhysicsBias()
Create biasing operations.
G4VisAttributes * visAttr
Cached global constants values.
G4UserLimits * userLimits
Cached global constants values.
Parser namespace for GMAD language. Combination of Geant4 and MAD.