BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
BDSFieldQuery Class Reference

Class for querying the Geant4 model for field at any point. More...

#include <BDSFieldQuery.hh>

Inheritance diagram for BDSFieldQuery:
Inheritance graph
Collaboration diagram for BDSFieldQuery:
Collaboration graph

Public Member Functions

virtual void QueryField (const BDSFieldQueryInfo *query)
 Query the field in the Geant4 model according to information in query. More...
 
void QueryFields (const std::vector< BDSFieldQueryInfo * > &fieldQueries)
 Vector version of above function. Unique output files for each query. More...
 
virtual void CleanUp ()
 Reset any member variables used during a query. Closes any files if open. More...
 

Static Public Member Functions

static void AttachWorldVolumeToNavigator (G4VPhysicalVolume *worldPVIn)
 Setup the navigator w.r.t. to a world volume. More...
 

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. More...
 

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. More...
 
void QuerySpecificPoints (const BDSFieldQueryInfo *query)
 Different algorithm where we query a specific list of points defined in the query info object. More...
 
virtual void PrintBAndEInfo (const BDSFieldQueryInfo *query) const
 Print out whether B and E files are being generated and which ones. More...
 
virtual void OpenFiles (const BDSFieldQueryInfo *query)
 Open potentially both electric and magnetic field map files. More...
 
virtual void WriteHeader (std::ofstream &out, const BDSFieldQueryInfo *query) const
 Utility to write required header information. More...
 
virtual void CloseFiles ()
 Close any files if open. More...
 
G4ThreeVector LocalToGlobalPoint (const G4AffineTransform &localToGlobalTransform, G4double xLocal, G4double yLocal, G4double zLocal) const
 Apply a transform to the coordinates. Does not apply to time. More...
 
void GlobalToLocalAxisField (const G4AffineTransform &globalToLocalTransform, const G4double globalBEField[6], G4double localBEField[6])
 Convert a global field axis to a local one. More...
 
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. More...
 

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. More...
 

Detailed Description

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.

Author
Laurie Nevay

Definition at line 42 of file BDSFieldQuery.hh.

Constructor & Destructor Documentation

◆ BDSFieldQuery()

BDSFieldQuery::BDSFieldQuery ( )

Definition at line 46 of file BDSFieldQuery.cc.

◆ ~BDSFieldQuery()

BDSFieldQuery::~BDSFieldQuery ( )
virtual

Definition at line 57 of file BDSFieldQuery.cc.

Member Function Documentation

◆ AttachWorldVolumeToNavigator()

void BDSFieldQuery::AttachWorldVolumeToNavigator ( G4VPhysicalVolume *  worldPVIn)
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().

Here is the caller graph for this function:

◆ CheckIfFieldObjectSpecified()

void BDSFieldQuery::CheckIfFieldObjectSpecified ( const BDSFieldQueryInfo query) const
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().

Here is the caller graph for this function:

◆ CheckNStepsAndRange()

void BDSFieldQuery::CheckNStepsAndRange ( const BDSFieldQueryInfo::QueryDimensionInfo dimensionInfo,
const G4String &  dimensionName,
const G4String &  queryName 
) const
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CleanUp()

void BDSFieldQuery::CleanUp ( )
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CloseFiles()

void BDSFieldQuery::CloseFiles ( )
privatevirtual

Close any files if open.

Reimplemented in BDSFieldQueryForVis.

Definition at line 328 of file BDSFieldQuery.cc.

Referenced by CleanUp(), QueryField(), and QuerySpecificPoints().

Here is the caller graph for this function:

◆ GetFieldValue()

void BDSFieldQuery::GetFieldValue ( const G4ThreeVector &  globalXYZ,
const G4ThreeVector &  globalDirection,
G4double  tGlobal,
G4double  fieldValue[6] 
)
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().

Here is the caller graph for this function:

◆ GlobalToLocalAxisField()

void BDSFieldQuery::GlobalToLocalAxisField ( const G4AffineTransform &  globalToLocalTransform,
const G4double  globalBEField[6],
G4double  localBEField[6] 
)
private

Convert a global field axis to a local one.

Definition at line 346 of file BDSFieldQuery.cc.

Referenced by QueryField().

Here is the caller graph for this function:

◆ LocalToGlobalPoint()

G4ThreeVector BDSFieldQuery::LocalToGlobalPoint ( const G4AffineTransform &  localToGlobalTransform,
G4double  xLocal,
G4double  yLocal,
G4double  zLocal 
) const
private

Apply a transform to the coordinates. Does not apply to time.

Definition at line 336 of file BDSFieldQuery.cc.

Referenced by QueryField().

Here is the caller graph for this function:

◆ OpenFiles()

void BDSFieldQuery::OpenFiles ( const BDSFieldQueryInfo query)
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ PrintBAndEInfo()

void BDSFieldQuery::PrintBAndEInfo ( const BDSFieldQueryInfo query) const
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().

Here is the caller graph for this function:

◆ QueryField()

void BDSFieldQuery::QueryField ( const BDSFieldQueryInfo query)
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ QuerySpecificPoints()

void BDSFieldQuery::QuerySpecificPoints ( const BDSFieldQueryInfo query)
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ WriteFieldValue()

void BDSFieldQuery::WriteFieldValue ( const G4ThreeVector &  xyzGlobal,
G4double  tGlobal,
const G4double  fieldValue[6] 
)
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().

Here is the caller graph for this function:

◆ WriteHeader()

void BDSFieldQuery::WriteHeader ( std::ofstream &  out,
const BDSFieldQueryInfo query 
) const
privatevirtual

Utility to write required header information.

Reimplemented in BDSFieldQueryForVis.

Definition at line 292 of file BDSFieldQuery.cc.

Referenced by OpenFiles().

Here is the caller graph for this function:

Field Documentation

◆ navigator

G4Navigator * BDSFieldQuery::navigator = new G4Navigator()
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().

◆ oFileElectric

std::ofstream BDSFieldQuery::oFileElectric
private

Definition at line 109 of file BDSFieldQuery.hh.

◆ oFileMagnetic

std::ofstream BDSFieldQuery::oFileMagnetic
private

Definition at line 108 of file BDSFieldQuery.hh.

◆ queryElectric

G4bool BDSFieldQuery::queryElectric
private

Definition at line 111 of file BDSFieldQuery.hh.

◆ queryMagnetic

G4bool BDSFieldQuery::queryMagnetic
private

Definition at line 110 of file BDSFieldQuery.hh.

◆ writeT

G4bool BDSFieldQuery::writeT
private

Definition at line 115 of file BDSFieldQuery.hh.

◆ writeX

G4bool BDSFieldQuery::writeX
private

Definition at line 112 of file BDSFieldQuery.hh.

◆ writeY

G4bool BDSFieldQuery::writeY
private

Definition at line 113 of file BDSFieldQuery.hh.

◆ writeZ

G4bool BDSFieldQuery::writeZ
private

Definition at line 114 of file BDSFieldQuery.hh.


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