BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSParticleCoords.cc
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#include "BDSParticleCoords.hh"
20
21#include "G4ThreeVector.hh"
22#include "G4Transform3D.hh"
23
24#include "CLHEP/Geometry/Point3D.h"
25
26#include <ostream>
27
28BDSParticleCoords::BDSParticleCoords():
29 x(0),
30 y(0),
31 z(0),
32 xp(0),
33 yp(0),
34 zp(0),
35 T(0)
36{;}
37
38BDSParticleCoords::BDSParticleCoords(G4double xIn,
39 G4double yIn,
40 G4double zIn,
41 G4double xpIn,
42 G4double ypIn,
43 G4double zpIn,
44 G4double tIn):
45 x(xIn),
46 y(yIn),
47 z(zIn),
48 xp(xpIn),
49 yp(ypIn),
50 zp(zpIn),
51 T(tIn)
52{;}
53
54BDSParticleCoords::BDSParticleCoords(G4ThreeVector pos,
55 G4ThreeVector mom,
56 G4double tIn):
57 x(pos.x()),
58 y(pos.y()),
59 z(pos.z()),
60 xp(mom.x()),
61 yp(mom.y()),
62 zp(mom.z()),
63 T(tIn)
64{;}
65
66BDSParticleCoords BDSParticleCoords::ApplyTransform(const G4Transform3D& transform) const
67{
68 // if no transform, return copy of self
69 if (transform == G4Transform3D::Identity)
70 {return BDSParticleCoords(*this);}
71
72 G4ThreeVector originalPos = G4ThreeVector(x,y,z);
73 G4ThreeVector newPos = transform * (HepGeom::Point3D<G4double>)originalPos;
74 G4ThreeVector originalMom = G4ThreeVector(xp,yp,zp);
75 G4ThreeVector newMom = transform * (HepGeom::Vector3D<G4double>)originalMom;
76 return BDSParticleCoords(newPos.x(), newPos.y(), newPos.z(),
77 newMom.x(), newMom.y(), newMom.z(),
78 T);
79}
80
81BDSParticleCoords BDSParticleCoords::ApplyOffset(const G4ThreeVector& offset) const
82{
83 return BDSParticleCoords(x + offset.x(),
84 y + offset.y(),
85 z + offset.z(),
86 xp, yp, zp, T);
87}
88
89void BDSParticleCoords::AddOffset(const G4ThreeVector& offset)
90{
91 x += offset.x();
92 y += offset.y();
93 z += offset.z();
94}
95
96std::ostream& operator<< (std::ostream& out, BDSParticleCoords const& p)
97{
98 p.Print(out);
99 return out;
100}
101
102void BDSParticleCoords::Print(std::ostream& out) const
103{
104 out << "Position: (" << x << ", " << y << ", " << z << ")";
105 out << ", Momentum: (" << xp << ", " << yp << ", " << zp << ")";
106 out << ", t: " << T << G4endl;
107}
A set of particle coordinates.
virtual void Print(std::ostream &out) const
Actual print out method so it can be called from a derived class.
void AddOffset(const G4ThreeVector &offset)
Apply an offset to the spatial coordinates only - assignment.
BDSParticleCoords ApplyTransform(const G4Transform3D &transform) const
Apply a transform to the coordinates and return a copy of them transformed.
BDSParticleCoords ApplyOffset(const G4ThreeVector &offset) const
Apply an offset to the spatial coordinates only and return a copy.