BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSAuxiliaryNavigator.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 BDSAUXILIARYNAVIGATOR_H
20#define BDSAUXILIARYNAVIGATOR_H
21
22#include "BDSMagnetStrength.hh"
23
24#include "globals.hh" // geant4 types / globals
25#include "G4AffineTransform.hh"
26#include "G4Navigator.hh"
27#include "G4ThreeVector.hh"
28#include "G4Transform3D.hh"
29
30class BDSStep;
31class G4Step;
32class G4VPhysicalVolume;
33
65{
66public:
69
71 static void AttachWorldVolumeToNavigator(G4VPhysicalVolume* worldPVIn)
72 {auxNavigator->SetWorldVolume(worldPVIn); worldPV = worldPVIn;}
73
76 static void AttachWorldVolumeToNavigatorCL(G4VPhysicalVolume* curvilinearWorldPVIn)
77 {auxNavigatorCL->SetWorldVolume(curvilinearWorldPVIn); curvilinearWorldPV = curvilinearWorldPVIn;}
78
79 static void RegisterCurvilinearBridgeWorld(G4VPhysicalVolume* curvilinearBridgeWorldPVIn)
80 {auxNavigatorCLB->SetWorldVolume(curvilinearBridgeWorldPVIn); curvilinearBridgeWorldPV = curvilinearBridgeWorldPVIn;}
81
82 static void ResetNavigatorStates();
83
85 G4VPhysicalVolume* LocateGlobalPointAndSetup(const G4ThreeVector& point,
86 const G4ThreeVector* direction = nullptr,
87 const G4bool pRelativeSearch = true,
88 const G4bool ignoreDirection = true,
89 G4bool useCurvilinear = true) const;
90
94 G4VPhysicalVolume* LocateGlobalPointAndSetup(G4Step const* const step,
95 G4bool useCurvilinear = true) const;
96
101 BDSStep ConvertToLocal(G4Step const* const step,
102 G4bool useCurvilinear = true) const;
103
111 BDSStep ConvertToLocal(const G4ThreeVector& globalPosition,
112 const G4ThreeVector& globalDirection,
113 const G4double stepLength = 0,
114 const G4bool useCurvilinear = true,
115 const G4double marginLength = 1) const;
116
121 BDSStep ConvertToGlobalStep(const G4ThreeVector& localPosition,
122 const G4ThreeVector& localDirection,
123 const G4bool useCurvilinear = true) const;
124
128 G4ThreeVector ConvertToLocal(const G4double globalPoint[3],
129 const G4bool useCurvilinear = true) const;
130
132 G4ThreeVector ConvertToLocal(const G4ThreeVector& globalPosition,
133 const G4bool useCurvilinear = true) const;
134
136 G4ThreeVector ConvertToLocalNoSetup(const G4ThreeVector& globalPosition,
137 const G4bool useCurvilinear = true) const;
138
144 G4ThreeVector ConvertAxisToLocal(const G4ThreeVector& globalAxis,
145 const G4bool useCurvilinear = true) const;
146
150 G4ThreeVector ConvertAxisToLocal(const G4double globalPoint[3],
151 const G4double globalAxis[3],
152 const G4bool useCurvilinear = true) const;
153
155 G4ThreeVector ConvertAxisToLocal(const G4ThreeVector& globalPoint,
156 const G4ThreeVector& globalAxis,
157 const G4bool useCurvilinear = true) const;
158
163 G4ThreeVector ConvertAxisToGlobal(const G4ThreeVector& localAxis,
164 const G4bool useCurvilinear = true) const;
165
170 std::pair<G4ThreeVector, G4ThreeVector> ConvertAxisToGlobal(const std::pair<G4ThreeVector, G4ThreeVector>& localAxis,
171 const G4bool useCurvilinear = true) const;
172
175 G4ThreeVector ConvertToGlobal(const G4ThreeVector& localPosition,
176 const G4bool useCurvilinear = true) const;
177
180 G4ThreeVector ConvertAxisToGlobal(const G4ThreeVector& globalPosition,
181 const G4ThreeVector& localAxis,
182 const G4bool useCurvilinear = true) const;
183
186 G4ThreeVector ConvertToGlobal(const G4ThreeVector& globalPosition,
187 const G4ThreeVector& localPosition,
188 const G4bool useCurvilinear = true) const;
189
192 BDSStep GlobalToCurvilinear(const G4double fieldArcLength,
193 const G4ThreeVector& unitField,
194 const G4double angle,
195 const G4ThreeVector& position,
196 const G4ThreeVector& unitMomentum,
197 const G4double h,
198 const G4bool useCurvilinearWorld,
199 const G4double FCof,
200 const G4double tilt = 0);
201
202 BDSStep GlobalToCurvilinear(const G4ThreeVector& position,
203 const G4ThreeVector& unitMomentum,
204 const G4double h,
205 const G4bool useCurvilinearWorld);
206
207 BDSStep CurvilinearToGlobal(const G4ThreeVector& localPosition,
208 const G4ThreeVector& localMomentum,
209 const G4bool useCurvilinearWorld);
210
211 BDSStep CurvilinearToGlobal(const G4double fieldArcLength,
212 const G4ThreeVector& unitField,
213 const G4double angle,
214 const G4ThreeVector& CLPosition,
215 const G4ThreeVector& CLMomentum,
216 const G4bool useCurvilinearWorld,
217 const G4double FCof,
218 const G4double tilt = 0);
219
220protected:
221 mutable G4AffineTransform globalToLocal;
222 mutable G4AffineTransform localToGlobal;
223 mutable G4AffineTransform globalToLocalCL;
224 mutable G4AffineTransform localToGlobalCL;
225 mutable G4bool bridgeVolumeWasUsed;
226
229 static G4Navigator* auxNavigator;
230
234 static G4Navigator* auxNavigatorCL;
235
239 static G4Navigator* auxNavigatorCLB;
240
241private:
243 G4Navigator* Navigator(G4bool curvilinear) const;
244
246 const G4AffineTransform& GlobalToLocal(G4bool curvilinear) const;
247 const G4AffineTransform& LocalToGlobal(G4bool curvilinear) const;
249
250 void InitialiseTransform(const G4bool massworld = true,
251 const G4bool curvilinearWorld = true) const;
252
259 void InitialiseTransform(const G4ThreeVector& globalPosition) const;
260
265 void InitialiseTransform(const G4ThreeVector& globalPosition,
266 const G4ThreeVector& globalMomentum,
267 const G4double stepLength);
268
271 static G4int numberOfInstances;
272
274 static G4VPhysicalVolume* worldPV;
275 static G4VPhysicalVolume* curvilinearWorldPV;
276 static G4VPhysicalVolume* curvilinearBridgeWorldPV;
278
283 G4double volumeMargin;
284};
285
286
287#endif
Extra G4Navigator to get coordinate transforms.
G4ThreeVector ConvertToGlobal(const G4ThreeVector &localPosition, const G4bool useCurvilinear=true) const
static G4VPhysicalVolume * curvilinearBridgeWorldPV
Cache of world PV to test if we're getting the wrong volume for the transform.
G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true, G4bool useCurvilinear=true) const
A wrapper for the underlying static navigator instance located within this class.
BDSStep GlobalToCurvilinear(const G4double fieldArcLength, const G4ThreeVector &unitField, const G4double angle, const G4ThreeVector &position, const G4ThreeVector &unitMomentum, const G4double h, const G4bool useCurvilinearWorld, const G4double FCof, const G4double tilt=0)
G4Navigator * Navigator(G4bool curvilinear) const
Utility function to select appropriate navigator.
static G4VPhysicalVolume * worldPV
Cache of world PV to test if we're getting the wrong volume for the transform.
G4ThreeVector ConvertAxisToLocal(const G4ThreeVector &globalAxis, const G4bool useCurvilinear=true) const
BDSStep ConvertToGlobalStep(const G4ThreeVector &localPosition, const G4ThreeVector &localDirection, const G4bool useCurvilinear=true) const
static G4Navigator * auxNavigator
G4ThreeVector ConvertToLocalNoSetup(const G4ThreeVector &globalPosition, const G4bool useCurvilinear=true) const
Similar to above function but does NOT initialise the transforms.
static void AttachWorldVolumeToNavigatorCL(G4VPhysicalVolume *curvilinearWorldPVIn)
G4ThreeVector ConvertAxisToGlobal(const G4ThreeVector &localAxis, const G4bool useCurvilinear=true) const
static G4Navigator * auxNavigatorCLB
static G4Navigator * auxNavigatorCL
const G4AffineTransform & LocalToGlobal(G4bool curvilinear) const
Utility function to select appropriate transform.
static G4VPhysicalVolume * curvilinearWorldPV
Cache of world PV to test if we're getting the wrong volume for the transform.
BDSStep ConvertToLocal(G4Step const *const step, G4bool useCurvilinear=true) const
static void AttachWorldVolumeToNavigator(G4VPhysicalVolume *worldPVIn)
Setup the navigator w.r.t. to a world volume - typically real world.
const G4AffineTransform & GlobalToLocal(G4bool curvilinear) const
Utility function to select appropriate transform.
A simple class to represent the positions of a step.
Definition: BDSStep.hh:33