BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldMagSkew.cc
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#include "BDSFieldMagSkew.hh"
20
21#include "globals.hh"
22#include "G4RotationMatrix.hh"
23
25 G4double angle):
26 field(fieldIn)
27{
28 rotation = new G4RotationMatrix();
29 antiRotation = new G4RotationMatrix();
30 rotation->rotateZ(angle);
31 antiRotation->rotateZ(-angle);
32
34}
35
36BDSFieldMagSkew::~BDSFieldMagSkew()
37{
38 delete rotation;
39 delete antiRotation;
40}
41
42G4ThreeVector BDSFieldMagSkew::GetField(const G4ThreeVector &position,
43 const G4double t) const
44{
45 G4ThreeVector rotatedPosition(position);
46 rotatedPosition = rotatedPosition.transform(*rotation);
47 G4ThreeVector normalField = field->GetField(rotatedPosition, t);
48 return (*antiRotation)*normalField;
49}
G4RotationMatrix * rotation
The rotation matrix used to rotate the coordinates.
G4RotationMatrix * antiRotation
The opposite rotation matrix used to transform the resultant field vector.
BDSFieldMagSkew()=delete
Private default constructor to force use of supplied ones.
BDSFieldMag *const field
virtual G4ThreeVector GetField(const G4ThreeVector &position, const G4double t=0) const
Get the field - local coordinates, and rotated.
Interface for static magnetic fields that may or may not be local.
Definition: BDSFieldMag.hh:39
G4bool FiniteStrength() const
Accessor.
Definition: BDSFieldMag.hh:80
virtual G4ThreeVector GetField(const G4ThreeVector &position, const G4double t=0) const =0
G4bool finiteStrength
Flag to cache whether finite nor not.
Definition: BDSFieldMag.hh:83