BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldQuery.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 BDSFIELDQUERY_H
20#define BDSFIELDQUERY_H
21#include "BDSFieldQueryInfo.hh"
22
23#include "G4AffineTransform.hh"
24#include "G4ThreeVector.hh"
25#include "G4Types.hh"
26
27#include <fstream>
28#include <vector>
29
30class G4Navigator;
31class G4VPhysicalVolume;
32
33
43{
44public:
46 virtual ~BDSFieldQuery();
47
49 virtual void QueryField(const BDSFieldQueryInfo* query);
51 void QueryFields(const std::vector<BDSFieldQueryInfo*>& fieldQueries);
52
54 static void AttachWorldVolumeToNavigator(G4VPhysicalVolume* worldPVIn);
55
57 virtual void CleanUp();
58
59protected:
63 virtual void GetFieldValue(const G4ThreeVector& globalXYZ,
64 const G4ThreeVector& globalDirection,
65 G4double tGlobal,
66 G4double fieldValue[6]);
67
69 virtual void CheckIfFieldObjectSpecified(const BDSFieldQueryInfo* query) const;
70
71private:
74 const G4String& dimensionName,
75 const G4String& queryName) const;
76
78 void QuerySpecificPoints(const BDSFieldQueryInfo* query);
79
81 virtual void PrintBAndEInfo(const BDSFieldQueryInfo* query) const;
82
84 virtual void OpenFiles(const BDSFieldQueryInfo* query);
85
87 virtual void WriteHeader(std::ofstream& out, const BDSFieldQueryInfo* query) const;
88
90 virtual void CloseFiles();
91
93 G4ThreeVector LocalToGlobalPoint(const G4AffineTransform& localToGlobalTransform,
94 G4double xLocal,
95 G4double yLocal,
96 G4double zLocal) const;
97
99 void GlobalToLocalAxisField(const G4AffineTransform& globalToLocalTransform,
100 const G4double globalBEField[6],
101 G4double localBEField[6]);
102
104 virtual void WriteFieldValue(const G4ThreeVector& xyzGlobal,
105 G4double tGlobal,
106 const G4double fieldValue[6]);
107
108 std::ofstream oFileMagnetic;
109 std::ofstream oFileElectric;
110 G4bool queryMagnetic;
111 G4bool queryElectric;
112 G4bool writeX;
113 G4bool writeY;
114 G4bool writeZ;
115 G4bool writeT;
116
118 static G4Navigator* navigator;
119};
120
121#endif
Holder class for all information required for a field query.
Class for querying the Geant4 model for field at any point.
static G4Navigator * navigator
One navigator used - static so we can externally link it up to the world PV.
static void AttachWorldVolumeToNavigator(G4VPhysicalVolume *worldPVIn)
Setup the navigator w.r.t. to a world volume.
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.
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.
G4ThreeVector LocalToGlobalPoint(const G4AffineTransform &localToGlobalTransform, G4double xLocal, G4double yLocal, G4double zLocal) const
Apply a transform to the coordinates. Does not apply to time.
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.
virtual void CheckIfFieldObjectSpecified(const BDSFieldQueryInfo *query) const
Warn the user if the fieldObject variable is use when it shouldn't be.
void QuerySpecificPoints(const BDSFieldQueryInfo *query)
Different algorithm where we query a specific list of points defined in the query info object.
virtual void CloseFiles()
Close any files if open.
virtual void PrintBAndEInfo(const BDSFieldQueryInfo *query) const
Print out whether B and E files are being generated and which ones.
virtual void QueryField(const BDSFieldQueryInfo *query)
Query the field in the Geant4 model according to information in query.
virtual void GetFieldValue(const G4ThreeVector &globalXYZ, const G4ThreeVector &globalDirection, G4double tGlobal, G4double fieldValue[6])
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,...
void GlobalToLocalAxisField(const G4AffineTransform &globalToLocalTransform, const G4double globalBEField[6], G4double localBEField[6])
Convert a global field axis to a local one.