BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Class for querying the Geant4 model for field at any point. More...
#include <BDSFieldQuery.hh>
Public Member Functions | |
virtual void | QueryField (const BDSFieldQueryInfo *query) |
Query the field in the Geant4 model according to information in query. | |
void | QueryFields (const std::vector< BDSFieldQueryInfo * > &fieldQueries) |
Vector version of above function. Unique output files for each query. | |
virtual void | CleanUp () |
Reset any member variables used during a query. Closes any files if open. | |
Static Public Member Functions | |
static void | AttachWorldVolumeToNavigator (G4VPhysicalVolume *worldPVIn) |
Setup the navigator w.r.t. to a world volume. | |
Protected Member Functions | |
virtual void | GetFieldValue (const G4ThreeVector &globalXYZ, const G4ThreeVector &globalDirection, G4double tGlobal, G4double fieldValue[6]) |
virtual void | CheckIfFieldObjectSpecified (const BDSFieldQueryInfo *query) const |
Warn the user if the fieldObject variable is use when it shouldn't be. | |
Private Member Functions | |
void | CheckNStepsAndRange (const BDSFieldQueryInfo::QueryDimensionInfo &dimensionInfo, const G4String &dimensionName, const G4String &queryName) const |
Throw an exception if the number of steps is >1 and the difference between max and min is 0. | |
void | QuerySpecificPoints (const BDSFieldQueryInfo *query) |
Different algorithm where we query a specific list of points defined in the query info object. | |
virtual void | PrintBAndEInfo (const BDSFieldQueryInfo *query) const |
Print out whether B and E files are being generated and which ones. | |
virtual void | OpenFiles (const BDSFieldQueryInfo *query) |
Open potentially both electric and magnetic field map files. | |
virtual void | WriteHeader (std::ofstream &out, const BDSFieldQueryInfo *query) const |
Utility to write required header information. | |
virtual void | CloseFiles () |
Close any files if open. | |
G4ThreeVector | LocalToGlobalPoint (const G4AffineTransform &localToGlobalTransform, G4double xLocal, G4double yLocal, G4double zLocal) const |
Apply a transform to the coordinates. Does not apply to time. | |
void | GlobalToLocalAxisField (const G4AffineTransform &globalToLocalTransform, const G4double globalBEField[6], G4double localBEField[6]) |
Convert a global field axis to a local one. | |
virtual void | WriteFieldValue (const G4ThreeVector &xyzGlobal, G4double tGlobal, const G4double fieldValue[6]) |
Write an entry ta line of the output file(s). The array is assumed to be Bx,By,Bz,Ex,Ey,Ez as in Geant4. | |
Private Attributes | |
std::ofstream | oFileMagnetic |
std::ofstream | oFileElectric |
G4bool | queryMagnetic |
G4bool | queryElectric |
G4bool | writeX |
G4bool | writeY |
G4bool | writeZ |
G4bool | writeT |
Static Private Attributes | |
static G4Navigator * | navigator = new G4Navigator() |
One navigator used - static so we can externally link it up to the world PV. | |
Class for querying the Geant4 model for field at any point.
Output is a BDSIM-format field map. Unique files for electric and magnetic field maps as would be required to read the field maps back into BDSIM.
Definition at line 42 of file BDSFieldQuery.hh.
BDSFieldQuery::BDSFieldQuery | ( | ) |
Definition at line 46 of file BDSFieldQuery.cc.
|
virtual |
Definition at line 57 of file BDSFieldQuery.cc.
|
static |
Setup the navigator w.r.t. to a world volume.
Definition at line 60 of file BDSFieldQuery.cc.
References navigator.
Referenced by BDSDetectorConstruction::BuildWorld().
|
protectedvirtual |
Warn the user if the fieldObject variable is use when it shouldn't be.
Reimplemented in BDSFieldQueryRaw.
Definition at line 200 of file BDSFieldQuery.cc.
References BDSFieldQueryInfo::fieldObject.
Referenced by QueryField().
|
private |
Throw an exception if the number of steps is >1 and the difference between max and min is 0.
Definition at line 184 of file BDSFieldQuery.cc.
References BDS::IsFinite().
Referenced by QueryField().
|
virtual |
Reset any member variables used during a query. Closes any files if open.
Reimplemented in BDSFieldQueryForVis.
Definition at line 65 of file BDSFieldQuery.cc.
References CloseFiles(), and navigator.
Referenced by BDSFieldQueryForVis::CleanUp(), and QueryField().
|
privatevirtual |
Close any files if open.
Reimplemented in BDSFieldQueryForVis.
Definition at line 328 of file BDSFieldQuery.cc.
Referenced by CleanUp(), QueryField(), and QuerySpecificPoints().
|
protectedvirtual |
Get the electric and magnetic field at the specified coordinates. The navigator requires a direction for safe hierarchy searching. The output is written to array argument, where the values are Bx,By,Bz,Ex,Ey,Ez as in Geant4.
Reimplemented in BDSFieldQueryForVis, and BDSFieldQueryRaw.
Definition at line 361 of file BDSFieldQuery.cc.
References navigator.
Referenced by BDSFieldQueryForVis::GetFieldValue(), QueryField(), and QuerySpecificPoints().
|
private |
Convert a global field axis to a local one.
Definition at line 346 of file BDSFieldQuery.cc.
Referenced by QueryField().
|
private |
Apply a transform to the coordinates. Does not apply to time.
Definition at line 336 of file BDSFieldQuery.cc.
Referenced by QueryField().
|
privatevirtual |
Open potentially both electric and magnetic field map files.
Reimplemented in BDSFieldQueryForVis.
Definition at line 245 of file BDSFieldQuery.cc.
References BDS::FileExists(), and WriteHeader().
Referenced by QueryField(), and QuerySpecificPoints().
|
privatevirtual |
Print out whether B and E files are being generated and which ones.
Reimplemented in BDSFieldQueryForVis.
Definition at line 234 of file BDSFieldQuery.cc.
Referenced by QueryField(), and QuerySpecificPoints().
|
virtual |
Query the field in the Geant4 model according to information in query.
Ensure field transform navigator is in a completely clean state.
Reimplemented in BDSFieldQueryForVis, and BDSFieldQueryRaw.
Definition at line 83 of file BDSFieldQuery.cc.
References CheckIfFieldObjectSpecified(), CheckNStepsAndRange(), BDSFieldQueryInfo::checkParameters, CleanUp(), CloseFiles(), GetFieldValue(), GlobalToLocalAxisField(), LocalToGlobalPoint(), OpenFiles(), PrintBAndEInfo(), QuerySpecificPoints(), BDSFieldQueryInfo::SpecificPoints(), and WriteFieldValue().
Referenced by BDSFieldQueryForVis::QueryField(), and QueryFields().
void BDSFieldQuery::QueryFields | ( | const std::vector< BDSFieldQueryInfo * > & | fieldQueries | ) |
Vector version of above function. Unique output files for each query.
Definition at line 77 of file BDSFieldQuery.cc.
References QueryField().
Referenced by BDSRunManager::Initialize().
|
private |
Different algorithm where we query a specific list of points defined in the query info object.
Definition at line 210 of file BDSFieldQuery.cc.
References CloseFiles(), GetFieldValue(), OpenFiles(), PrintBAndEInfo(), and WriteFieldValue().
Referenced by QueryField().
|
privatevirtual |
Write an entry ta line of the output file(s). The array is assumed to be Bx,By,Bz,Ex,Ey,Ez as in Geant4.
Reimplemented in BDSFieldQueryForVis.
Definition at line 383 of file BDSFieldQuery.cc.
Referenced by QueryField(), and QuerySpecificPoints().
|
privatevirtual |
Utility to write required header information.
Reimplemented in BDSFieldQueryForVis.
Definition at line 292 of file BDSFieldQuery.cc.
Referenced by OpenFiles().
|
staticprivate |
One navigator used - static so we can externally link it up to the world PV.
Definition at line 118 of file BDSFieldQuery.hh.
Referenced by AttachWorldVolumeToNavigator(), CleanUp(), and GetFieldValue().
|
private |
Definition at line 109 of file BDSFieldQuery.hh.
|
private |
Definition at line 108 of file BDSFieldQuery.hh.
|
private |
Definition at line 111 of file BDSFieldQuery.hh.
|
private |
Definition at line 110 of file BDSFieldQuery.hh.
|
private |
Definition at line 115 of file BDSFieldQuery.hh.
|
private |
Definition at line 112 of file BDSFieldQuery.hh.
|
private |
Definition at line 113 of file BDSFieldQuery.hh.
|
private |
Definition at line 114 of file BDSFieldQuery.hh.