BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSSurvey.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 "BDSSurvey.hh"
20
21#include "BDSAcceleratorComponent.hh"
22#include "BDSBeamline.hh"
23#include "BDSBeamlineElement.hh"
24#include "BDSBeamPipeInfo.hh"
25#include "BDSDebug.hh"
26#include "BDSMagnet.hh"
27#include "BDSMagnetStrength.hh"
28
29#include <fstream>
30#include <iomanip>
31#include <ctime>
32#include <string>
33
34using std::setw;
35
36BDSSurvey::BDSSurvey(G4String filename):
37 nullStrength(new BDSMagnetStrength()),
38 gp(14)
39{
40 magnetKeys = nullStrength->AllKeys();
41
42 G4cout << __METHOD_NAME__ << "Generating Survey: " << filename << " ..." << G4endl;
43 survey.open(filename);
44 WriteHeader();
45}
46
47BDSSurvey::~BDSSurvey()
48{
49 if (survey.is_open())
50 {survey.close();}
51 delete nullStrength;
52}
53
55{
56 time_t currenttime;
57 time(&currenttime);
58 std::string timestring = asctime(localtime(&currenttime));
59 timestring = timestring.substr(0,timestring.size()-1);
60
61 survey << "### BDSIM output - created "<< timestring << G4endl;
62 survey << std::left
63 << setw(15) << "Type " << " "
64 << setw(40) << "Name " << " "
65 << setw(gp) << "SStart[m] " << " "
66 << setw(gp) << "SMid[m] " << " "
67 << setw(gp) << "SEnd[m] " << " "
68 << setw(gp) << "ChordLength[m]" << " "
69 << setw(gp) << "ArcLength[m] " << " "
70 << setw(gp) << "X[m] " << " "
71 << setw(gp) << "Y[m] " << " "
72 << setw(gp) << "Z[m] " << " "
73 << setw(gp) << "Phi[rad] " << " "
74 << setw(gp) << "Theta[rad] " << " "
75 << setw(gp) << "Psi[rad] " << " "
76 << setw(gp) << "Aper1[m] " << " "
77 << setw(gp) << "Aper2[m] " << " "
78 << setw(gp) << "Aper3[m] " << " "
79 << setw(gp) << "Aper4[m] " << " "
80 << setw(15) << "Aper_Type " << " "
81 << setw(gp) << "Angle[rad] ";
82
83 for (auto const& key : magnetKeys)
84 {
85 const G4String unit = nullStrength->UnitName(key);
86 if (!unit.empty())
87 {survey << " " << setw(18) << key + "[" + unit + "]";}
88 else
89 {survey << " " << setw(18) << key;}
90 }
91 survey << G4endl;
92
93 survey.setf(std::ios::fixed, std::ios::floatfield);
94 survey.setf(std::ios::showpoint);
95 // print up to 7 significant numbers
96 survey.precision(7);
97}
98
100{
101 for (auto element : *beamline)
102 {Write(element);}
103
104 survey << "### Total length = " << beamline->GetTotalChordLength()/CLHEP::m << "m" << G4endl;
105 survey << "### Total arc length = " << beamline->GetTotalArcLength()/CLHEP::m << "m" << G4endl;
106 survey.close();
107}
108
110{
111 BDSAcceleratorComponent* acceleratorComponent = beamlineElement->GetAcceleratorComponent();
112
113 G4RotationMatrix* rm = beamlineElement->GetRotationMiddle();
114 G4double phi = rm->getPhi();
115 G4double theta = rm->getTheta();
116 G4double psi = rm->getPsi();
117
118 G4double sStart = beamlineElement->GetSPositionStart() /CLHEP::m;
119 G4double sMiddle = beamlineElement->GetSPositionMiddle()/CLHEP::m;
120 G4double sEnd = beamlineElement->GetSPositionEnd() /CLHEP::m;
121 G4ThreeVector pos = beamlineElement->GetPositionMiddle();
122
123 BDSBeamPipeInfo* beamPipeInfo = acceleratorComponent->GetBeamPipeInfo();
124
125 survey << std::left << std::setprecision(6) << std::fixed
126 << setw(15) << acceleratorComponent->GetType() << " "
127 << setw(40) << acceleratorComponent->GetName() << " "
128 << setw(gp) << sStart << " "
129 << setw(gp) << sMiddle << " "
130 << setw(gp) << sEnd << " "
131 << setw(gp) << acceleratorComponent->GetChordLength()/CLHEP::m << " "
132 << setw(gp) << acceleratorComponent->GetArcLength()/CLHEP::m << " "
133 << setw(gp) << pos.x()/CLHEP::m << " "
134 << setw(gp) << pos.y()/CLHEP::m << " "
135 << setw(gp) << pos.z()/CLHEP::m << " "
136 << setw(gp) << phi/CLHEP::radian << " "
137 << setw(gp) << theta/CLHEP::radian << " "
138 << setw(gp) << psi/CLHEP::radian << " "
139 << setw(gp) << (beamPipeInfo ? beamPipeInfo->aper1/CLHEP::m : 0) << " "
140 << setw(gp) << (beamPipeInfo ? beamPipeInfo->aper2/CLHEP::m : 0) << " "
141 << setw(gp) << (beamPipeInfo ? beamPipeInfo->aper3/CLHEP::m : 0) << " "
142 << setw(gp) << (beamPipeInfo ? beamPipeInfo->aper4/CLHEP::m : 0) << " "
143 << setw(15) << (beamPipeInfo ? beamPipeInfo->beamPipeType : 0) << " "
144 << setw(gp) << acceleratorComponent->GetAngle();
145
146 BDSMagnetStrength const* st;
147 if (BDSMagnet* magCast = dynamic_cast<BDSMagnet*>(acceleratorComponent))
148 {
149 st = magCast->MagnetStrength();
150 if (!st)
151 {st = nullStrength;}
152 }
153 else
154 {st = nullStrength;}
155
157 survey << G4endl;
158}
Abstract class that represents a component of an accelerator.
virtual G4String GetName() const
The name of the component without modification.
virtual G4double GetChordLength() const
virtual G4double GetArcLength() const
BDSBeamPipeInfo * GetBeamPipeInfo() const
G4String GetType() const
Get a string describing the type of the component.
Holder class for all information required to describe a beam pipe model.
G4double aper3
Public member for direct access.
G4double aper1
Public member for direct access.
G4double aper4
Public member for direct access.
G4double aper2
Public member for direct access.
BDSBeamPipeType beamPipeType
Public member for direct access.
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4ThreeVector GetPositionMiddle() const
Accessor.
G4double GetSPositionEnd() const
Accessor.
BDSAcceleratorComponent * GetAcceleratorComponent() const
Accessor.
G4double GetSPositionMiddle() const
Accessor.
G4RotationMatrix * GetRotationMiddle() const
Accessor.
G4double GetSPositionStart() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
G4double GetTotalChordLength() const
Get the total length of the beamline - the sum of the chord length of each element.
Definition: BDSBeamline.hh:136
G4double GetTotalArcLength() const
Get the total ARC length for the beamline - ie the maximum s position.
Definition: BDSBeamline.hh:139
Efficient storage of magnet strengths.
std::ostream & WriteValuesInSIUnitsForSuvey(std::ostream &out, const G4int precision=14) const
static G4String UnitName(const G4String &key)
Access a unit name for a given key.
Abstract base class that implements features common to all magnets.
Definition: BDSMagnet.hh:45
const int gp
General precision - number of characters / item.
Definition: BDSSurvey.hh:69
BDSMagnetStrength const * nullStrength
Definition: BDSSurvey.hh:64
std::ofstream survey
Output file stream.
Definition: BDSSurvey.hh:59
std::vector< G4String > magnetKeys
Cache of all the possible magnet strength parameters.
Definition: BDSSurvey.hh:67
void WriteHeader()
Write header.
Definition: BDSSurvey.cc:54
void Write(BDSBeamlineElement *beamlineElement)
write line
Definition: BDSSurvey.cc:109