BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSNavigatorPlacements.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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 BDSNAVIGATORPLACEMENTS_H
20#define BDSNAVIGATORPLACEMENTS_H
21
22#include "G4AffineTransform.hh"
23#include "G4Navigator.hh"
24#include "G4ThreeVector.hh"
25#include "G4Types.hh"
26
27#include <utility>
28
29class G4VPhysicalVolume;
30
50{
51public:
54
56 static void AttachWorldVolumeToNavigator(G4VPhysicalVolume* worldPVIn)
57 {navigator->SetWorldVolume(worldPVIn); worldPV = worldPVIn;}
58
59 static void ResetNavigatorStates();
60
64 G4ThreeVector ConvertToLocal(const G4ThreeVector& globalPosition,
65 G4bool& foundAPlacementVolume) const;
66
68 G4ThreeVector ConvertToLocalNoSetup(const G4ThreeVector& globalPosition) const;
69
74 G4ThreeVector ConvertAxisToGlobal(const G4ThreeVector& localAxis) const;
75
80 std::pair<G4ThreeVector, G4ThreeVector> ConvertAxisToGlobal(const std::pair<G4ThreeVector, G4ThreeVector>& localAxis) const;
81
82protected:
83 mutable G4AffineTransform globalToLocal;
84 mutable G4AffineTransform localToGlobal;
85
88 static G4Navigator* navigator;
89
90private:
92 inline const G4AffineTransform& GlobalToLocal() const {return globalToLocal;}
93 inline const G4AffineTransform& LocalToGlobal() const {return localToGlobal;}
95
102 G4bool InitialiseTransform(const G4ThreeVector& globalPosition) const;
103
106 static G4int numberOfInstances;
107
109 static G4VPhysicalVolume* worldPV;
110};
111
112
113#endif
Extra G4Navigator to get coordinate transforms for placement world.
const G4AffineTransform & LocalToGlobal() const
Utility function to select appropriate transform.
G4bool InitialiseTransform(const G4ThreeVector &globalPosition) const
static void AttachWorldVolumeToNavigator(G4VPhysicalVolume *worldPVIn)
Setup the navigator w.r.t. to a world volume - typically real world.
G4ThreeVector ConvertToLocalNoSetup(const G4ThreeVector &globalPosition) const
Similar to above function but does NOT initialise the transforms.
const G4AffineTransform & GlobalToLocal() const
Utility function to select appropriate transform.
G4ThreeVector ConvertToLocal(const G4ThreeVector &globalPosition, G4bool &foundAPlacementVolume) const
G4ThreeVector ConvertAxisToGlobal(const G4ThreeVector &localAxis) const
static G4VPhysicalVolume * worldPV
Cache of world PV to test if we're getting the wrong volume for the transform.
static G4Navigator * navigator