BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Class that constructs a parallel importance world. More...
#include <BDSParallelWorldImportance.hh>
Public Member Functions | |
BDSParallelWorldImportance (G4String name, G4String importanceWorldGeometryFile, G4String importanceValuesFile) | |
void | Construct () |
void | BuildPhysicsBias () |
Create biasing operations. | |
BDSExtent | WorldExtent () const |
Public access to the world extent. | |
G4GeometryCell | GetGeometryCell (G4int i) const |
Get geometry cell from store. | |
void | AddIStore () |
Create IStore for all importance sampling geometry cells. | |
virtual void | ConstructSD () |
G4VPhysicalVolume * | GetWorldVolume () |
World volume getter required in parallel world utilities. | |
Private Member Functions | |
BDSParallelWorldImportance & | operator= (const BDSParallelWorldImportance &) |
assignment and copy constructor not implemented nor used | |
BDSParallelWorldImportance (BDSParallelWorldImportance &) | |
void | BuildWorld () |
G4double | GetCellImportanceValue (const G4String &cellName) |
Get importance value of a given physical volume name. | |
Private Attributes | |
G4VPhysicalVolume * | imWorldPV |
Importance sampling world volume. | |
BDSExtent | worldExtent |
Record of the world extent. | |
BDSImportanceVolumeStore | imVolumeStore |
Store of geometry cells for importance volumes. | |
std::map< G4String, G4double > | imVolumesAndValues |
Container for all user placed physical volumes and corresponding importance values. | |
G4String | imGeomFile |
G4String | imVolMap |
const G4String | componentName |
String preprended to geometry with preprocessGDML. | |
G4int | verbosity |
Cached global constants values. | |
G4UserLimits * | userLimits |
Cached global constants values. | |
G4VisAttributes * | visAttr |
Cached global constants values. | |
Class that constructs a parallel importance world.
Definition at line 47 of file BDSParallelWorldImportance.hh.
BDSParallelWorldImportance::BDSParallelWorldImportance | ( | G4String | name, |
G4String | importanceWorldGeometryFile, | ||
G4String | importanceValuesFile | ||
) |
Definition at line 45 of file BDSParallelWorldImportance.cc.
|
virtual |
Definition at line 85 of file BDSParallelWorldImportance.cc.
void BDSParallelWorldImportance::AddIStore | ( | ) |
Create IStore for all importance sampling geometry cells.
Definition at line 147 of file BDSParallelWorldImportance.cc.
References GetCellImportanceValue(), imVolumesAndValues, imVolumeStore, imWorldPV, and verbosity.
Referenced by BDS::AddIStore().
|
private |
Build the world volume using the extent of the BDSBeamline instance created in BuildBeamline()
Definition at line 88 of file BDSParallelWorldImportance.cc.
References BDSImportanceVolumeStore::AddPVolume(), BDSGeometryFactory::BuildGeometry(), componentName, BDSGeometryComponent::GetContainerLogicalVolume(), imVolumeStore, imWorldPV, BDSGeometryFactory::Instance(), userLimits, and visAttr.
Referenced by Construct().
void BDSParallelWorldImportance::Construct | ( | ) |
Overridden Geant4 method that must be implemented. Constructs the Geant4 geometry and returns the finished world physical volume.
Definition at line 62 of file BDSParallelWorldImportance.cc.
References BuildWorld(), BDS::GetFullPath(), and imVolumesAndValues.
|
virtual |
Definition at line 236 of file BDSParallelWorldImportance.cc.
|
private |
Get importance value of a given physical volume name.
Definition at line 184 of file BDSParallelWorldImportance.cc.
References componentName, imVolumesAndValues, BDS::IsFinite(), and BDS::StrContains().
Referenced by AddIStore().
G4GeometryCell BDSParallelWorldImportance::GetGeometryCell | ( | G4int | i | ) | const |
Get geometry cell from store.
Definition at line 135 of file BDSParallelWorldImportance.cc.
References BDSImportanceVolumeStore::GetPVolume(), imVolumeStore, and imWorldPV.
|
inline |
World volume getter required in parallel world utilities.
Definition at line 74 of file BDSParallelWorldImportance.hh.
References imWorldPV.
Referenced by BDS::RegisterImportanceBiasing().
|
inline |
Public access to the world extent.
Definition at line 63 of file BDSParallelWorldImportance.hh.
References worldExtent.
|
private |
String preprended to geometry with preprocessGDML.
Definition at line 98 of file BDSParallelWorldImportance.hh.
Referenced by BuildWorld(), and GetCellImportanceValue().
|
private |
Definition at line 96 of file BDSParallelWorldImportance.hh.
|
private |
Definition at line 97 of file BDSParallelWorldImportance.hh.
|
private |
Container for all user placed physical volumes and corresponding importance values.
Definition at line 94 of file BDSParallelWorldImportance.hh.
Referenced by AddIStore(), Construct(), and GetCellImportanceValue().
|
private |
Store of geometry cells for importance volumes.
Definition at line 91 of file BDSParallelWorldImportance.hh.
Referenced by AddIStore(), BuildWorld(), and GetGeometryCell().
|
private |
Importance sampling world volume.
Definition at line 86 of file BDSParallelWorldImportance.hh.
Referenced by AddIStore(), BuildWorld(), GetGeometryCell(), and GetWorldVolume().
|
private |
Cached global constants values.
Definition at line 102 of file BDSParallelWorldImportance.hh.
Referenced by BuildWorld().
|
private |
Cached global constants values.
Definition at line 101 of file BDSParallelWorldImportance.hh.
Referenced by AddIStore().
|
private |
Cached global constants values.
Definition at line 103 of file BDSParallelWorldImportance.hh.
Referenced by BuildWorld().
|
private |
Record of the world extent.
Definition at line 88 of file BDSParallelWorldImportance.hh.
Referenced by WorldExtent().