BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSOutputROOTEventModel.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 "BDSOutputROOTEventModel.hh"
20
21#ifndef __ROOTBUILD__
22#include "BDSAcceleratorModel.hh"
23#include "BDSBeamline.hh"
24#include "BDSBeamlineElement.hh"
25#include "BDSBeamPipeInfo.hh"
26#include "BDSMagnet.hh"
27#include "BDSMagnetStrength.hh"
28#include "BDSSamplerRegistry.hh"
29
30#include "G4RotationMatrix.hh"
31#include "G4String.hh"
32#include "G4Types.hh"
33
34#include <algorithm>
35#include <map>
36#include <string>
37#include <vector>
38#include <utility>
39#endif
40
42
44 n(0),
45 storeCollimatorInfo(false),
46 nCollimators(0)
47{
48 Flush();
49}
50
52{;}
53
54int BDSOutputROOTEventModel::findNearestElement(const TVector3& point) const
55{
56 // we make a vector of pairs of the distance between the mid of each element
57 // and the specified point along with the index. We then sort this vector by
58 // the distance (first part of pair only). Could use a map as it's sorted but
59 // can't guarantee access with a double or float as a key due to binary representation.
60 std::vector<std::pair<float, int> > distanceAndIndex;
61 distanceAndIndex.reserve(midRefPos.size());
62 for (int i = 0; i < (int)midRefPos.size(); i++)
63 {distanceAndIndex.emplace_back(std::make_pair((midRefPos[i]-point).Mag(), i));}
64
65 struct customLess
66 {bool operator()(std::pair<float, int>& a, std::pair<float, int>& b) const {return a.first < b.first;};};
67 std::sort(distanceAndIndex.begin(), distanceAndIndex.end(), customLess());
68 return distanceAndIndex[0].second;
69}
70
72{
73 n = 0;
74 componentName.clear();
75 placementName.clear();
76 componentType.clear();
77 length.clear();
78 angle.clear();
79 staPos.clear();
80 midPos.clear();
81 endPos.clear();
82 staRot.clear();
83 midRot.clear();
84 endRot.clear();
85 staRefPos.clear();
86 midRefPos.clear();
87 endRefPos.clear();
88 staRefRot.clear();
89 midRefRot.clear();
90 endRefRot.clear();
91 tilt.clear();
92 offsetX.clear();
93 offsetY.clear();
94 staS.clear();
95 midS.clear();
96 endS.clear();
97 beamPipeType.clear();
98 beamPipeAper1.clear();
99 beamPipeAper2.clear();
100 beamPipeAper3.clear();
101 beamPipeAper4.clear();
102 material.clear();
103 k1.clear();
104 k2.clear();
105 k3.clear();
106 k4.clear();
107 k5.clear();
108 k6.clear();
109 k7.clear();
110 k8.clear();
111 k9.clear();
112 k10.clear();
113 k11.clear();
114 k12.clear();
115 k1s.clear();
116 k2s.clear();
117 k3s.clear();
118 k4s.clear();
119 k5s.clear();
120 k6s.clear();
121 k7s.clear();
122 k8s.clear();
123 k9s.clear();
124 k10s.clear();
125 k11s.clear();
126 k12s.clear();
127 ks.clear();
128 hkick.clear();
129 vkick.clear();
130 bField.clear();
131 eField.clear();
132 e1.clear();
133 e2.clear();
134 fint.clear();
135 fintx.clear();
136 fintk2.clear();
137 fintxk2.clear();
138
139 storeCollimatorInfo = false;
140 collimatorIndices.clear();
142 nCollimators = 0;
143 collimatorInfo.clear();
145 scoringMeshTranslation.clear();
146 scoringMeshRotation.clear();
147 scoringMeshName.clear();
148
149 materialIDToName.clear();
150 materialNameToID.clear();
151
152 samplerNamesUnique.clear();
153 samplerSPosition.clear();
154 samplerCNamesUnique.clear();
155 samplerSNamesUnique.clear();
156 samplerCRadius.clear();
157 samplerSRadius.clear();
158}
159
160#ifndef __ROOTBUILD__
162 n(0),
163 storeCollimatorInfo(storeCollimatorInfoIn),
164 nCollimators(0)
165{;}
166
167TRotation BDSOutputROOTEventModel::ConvertToROOT(const G4RotationMatrix* rm) const
168{
169 TRotation rr = TRotation();
170 rr.SetToIdentity();
171 if (!rm)
172 {return rr;}
173
174 double rotAngle;
175 CLHEP::Hep3Vector axis;
176
177 rm->getAngleAxis(rotAngle, axis);
178 rr.Rotate(rotAngle, TVector3(axis.x(), axis.y(), axis.z()));
179 return rr;
180}
181
182TRotation BDSOutputROOTEventModel::ConvertToROOT(const G4RotationMatrix& rm) const
183{
184 return ConvertToROOT(&rm);
185}
186
187TVector3 BDSOutputROOTEventModel::ConvertToROOT(const G4ThreeVector& v) const
188{
189 return TVector3(v.x() / CLHEP::m, v.y() / CLHEP::m, v.z() / CLHEP::m);
190}
191
192void BDSOutputROOTEventModel::Fill(const std::vector<G4int>& collimatorIndicesIn,
193 const std::map<G4String, G4int>& collimatorIndicesByNameIn,
194 const std::vector<BDSOutputROOTEventCollimatorInfo>& collimatorInfoIn,
195 const std::vector<G4String>& collimatorBranchNamesIn,
196 const std::map<G4String, G4Transform3D>* scorerMeshPlacements,
197 const std::map<short int, G4String>* materialIDToNameUnique,
198 G4bool storeTrajectory)
199{
201 for (const auto& nameSPos : sr->GetUniquePlaneNamesAndSPosition())
202 {
203 samplerNamesUnique.push_back(std::string(nameSPos.first) + ".");
204 samplerSPosition.push_back((double) nameSPos.second / CLHEP::m);
205 }
206 for (const auto& name : sr->GetUniqueNamesCylinder())
207 {samplerCNamesUnique.push_back(std::string(name) + ".");}
208 for (const auto& name : sr->GetUniqueNamesSphere())
209 {samplerSNamesUnique.push_back(std::string(name) + ".");}
210 samplerCRadius = sr->GetUniqueNameToRadiusCylinder();
211 samplerSRadius = sr->GetUniqueNameToRadiusSphere();
212
213 for (const auto& name : collimatorBranchNamesIn)
214 {collimatorBranchNamesUnique.push_back(std::string(name) + ".");}
215
216 if (scorerMeshPlacements)
217 {
218 for (const auto& kv : *scorerMeshPlacements)
219 {
220 const G4String& name = kv.first;
221 scoringMeshTranslation[name] = ConvertToROOT(kv.second.getTranslation());
222 scoringMeshRotation[name] = ConvertToROOT(kv.second.getRotation());
223 scoringMeshName.emplace_back(name);
224 }
225 }
226
227 const BDSBeamline* beamline = BDSAcceleratorModel::Instance()->BeamlineMain();
228 if (!beamline)
229 {return;} // in case of generatePrimariesOnly there is no model - return
230
232 {
233 for (const auto value : collimatorIndicesIn)
234 {collimatorIndices.push_back((int)value);}
235
236 nCollimators = (int)collimatorIndices.size();
237
238 for (const auto& kv : collimatorIndicesByNameIn)
239 {collimatorIndicesByName[(std::string)kv.first] = (int)kv.second;}
240
241 collimatorInfo = collimatorInfoIn;
242 }
243
244 if (materialIDToNameUnique && storeTrajectory)
245 {
246 for (const auto& kv : *materialIDToNameUnique)
247 {
248 materialIDToName[kv.first] = (std::string)kv.second;
249 materialNameToID[(std::string)kv.second] = kv.first;
250 }
251 }
252
253 n = (int)beamline->size();
254
255 for (auto i = beamline->begin(); i != beamline->end(); ++i)
256 {
257 componentName.push_back((*i)->GetName());
258 placementName.push_back((*i)->GetPlacementName());
259 componentType.push_back((*i)->GetType());
260 length.push_back((float)((*i)->GetAcceleratorComponent()->GetArcLength() / CLHEP::m));
261 angle.push_back((float)((*i)->GetAcceleratorComponent()->GetAngle() / CLHEP::radian));
262 staPos.push_back(ConvertToROOT((*i)->GetPositionStart()));
263 midPos.push_back(ConvertToROOT((*i)->GetPositionMiddle()));
264 endPos.push_back(ConvertToROOT((*i)->GetPositionEnd()));
265 staRot.push_back(ConvertToROOT((*i)->GetRotationStart()));
266 midRot.push_back(ConvertToROOT((*i)->GetRotationMiddle()));
267 endRot.push_back(ConvertToROOT((*i)->GetRotationEnd()));
268 staRefPos.emplace_back(ConvertToROOT((*i)->GetReferencePositionStart()));
269 midRefPos.emplace_back(ConvertToROOT((*i)->GetReferencePositionMiddle()));
270 endRefPos.emplace_back(ConvertToROOT((*i)->GetReferencePositionEnd()));
271 staRefRot.push_back(ConvertToROOT((*i)->GetReferenceRotationStart()));
272 midRefRot.push_back(ConvertToROOT((*i)->GetReferenceRotationMiddle()));
273 endRefRot.push_back(ConvertToROOT((*i)->GetReferenceRotationEnd()));
274
275 // tilt and offset
276 BDSTiltOffset* to = (*i)->GetTiltOffset();
277 if (to)
278 {
279 tilt.push_back((float)(to->GetTilt() / CLHEP::rad));
280 offsetX.push_back((float)(to->GetXOffset() / CLHEP::m));
281 offsetY.push_back((float)(to->GetYOffset() / CLHEP::m));
282 }
283 else
284 {
285 tilt.push_back(0);
286 offsetX.push_back(0);
287 offsetY.push_back(0);
288 }
289
290 // S positions
291 staS.push_back((float)((*i)->GetSPositionStart() / CLHEP::m));
292 midS.push_back((float)((*i)->GetSPositionMiddle() / CLHEP::m));
293 endS.push_back((float)((*i)->GetSPositionEnd() / CLHEP::m));
294
295 // beam pipe
296 BDSBeamPipeInfo* beampipeinfo = (*i)->GetBeamPipeInfo();
297 beamPipeType.push_back(beampipeinfo ? beampipeinfo->beamPipeType.ToString() : "");
298 beamPipeAper1.push_back(beampipeinfo ? beampipeinfo->aper1 / CLHEP::m : 0);
299 beamPipeAper2.push_back(beampipeinfo ? beampipeinfo->aper2 / CLHEP::m : 0);
300 beamPipeAper3.push_back(beampipeinfo ? beampipeinfo->aper3 / CLHEP::m : 0);
301 beamPipeAper4.push_back(beampipeinfo ? beampipeinfo->aper4 / CLHEP::m : 0);
302
303 // associated material if any
304 const auto accComp = (*i)->GetAcceleratorComponent();
305 material.push_back(accComp->Material());
306
307 // helper shortcuts to all the member vectors
308 std::vector<std::vector<float>*> localNorm = {&k1,&k2,&k3,&k4,&k5,&k6,&k7,&k8,&k9,&k10,&k11,&k12};
309 std::vector<std::vector<float>*> localSkew = {&k1s,&k2s,&k3s,&k4s,&k5s,&k6s,&k7s,&k8s,&k9s,&k10s,&k11s,&k12s};
310
311 // helper lambda to avoid duplication
312 auto fillzero = [&]
313 {
314 for (int j = 0; j < (int)localNorm.size(); j++)
315 {localNorm[j]->push_back(0);}
316 for (int j = 0; j < (int)localSkew.size(); j++)
317 {localSkew[j]->push_back(0);}
318 ks.push_back(0);
319 hkick.push_back(0);
320 vkick.push_back(0);
321 bField.push_back(0);
322 eField.push_back(0);
323 e1.push_back(0);
324 e2.push_back(0);
325 hgap.push_back(0);
326 fint.push_back(0);
327 fintx.push_back(0);
328 fintk2.push_back(0);
329 fintxk2.push_back(0);
330 };
331
332 // fill magnet strength data
333 if (BDSMagnet* mag = dynamic_cast<BDSMagnet*>(accComp))
334 {
335 const BDSMagnetStrength* ms = mag->MagnetStrength();
336 if (!ms)
337 {
338 fillzero();
339 continue;
340 }
341 // assume localNorm and normComponents are same size
342 std::vector<G4double> normComponents = ms->NormalComponents();
343 for (int j = 0; j < (int)localNorm.size(); j++)
344 {localNorm[j]->push_back((float)normComponents[j]);}
345 std::vector<G4double> skewComponents = ms->SkewComponents();
346 for (int j = 0; j < (int)localSkew.size(); j++)
347 {localSkew[j]->push_back((float)skewComponents[j]);}
348
349 ks.push_back((float)((*ms)["ks"]/BDSMagnetStrength::Unit("ks")));
350 hkick.push_back( (float)((*ms)["hkick"]/BDSMagnetStrength::Unit("hkick")));
351 vkick.push_back( (float)((*ms)["vkick"]/BDSMagnetStrength::Unit("vkick")));
352 bField.push_back((float)((*ms)["field"]/BDSMagnetStrength::Unit("field")));
353 eField.push_back((float)((*ms)["efield"]/BDSMagnetStrength::Unit("efield")));
354 e1.push_back((float)((*ms)["e1"]/BDSMagnetStrength::Unit("e1")));
355 e2.push_back((float)((*ms)["e2"]/BDSMagnetStrength::Unit("e2")));
356 hgap.push_back((float)((*ms)["hgap"]/BDSMagnetStrength::Unit("hgap")));
357 fint.push_back((float)(*ms)["fint"]);
358 fintx.push_back((float)(*ms)["fintx"]);
359 fintk2.push_back((float)(*ms)["fintk2"]);
360 fintxk2.push_back((float)(*ms)["fintxk2"]);
361 }
362 else
363 {// not a magnet
364 fillzero();
365 }
366 }
367}
368#endif
const BDSBeamline * BeamlineMain() const
Accessor.
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 vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
BeamlineVector::size_type size() const
Get the number of elements.
Definition: BDSBeamline.hh:148
iterator begin()
Iterator mechanics.
Definition: BDSBeamline.hh:167
iterator end()
Iterator mechanics.
Definition: BDSBeamline.hh:168
Efficient storage of magnet strengths.
std::vector< G4double > NormalComponents() const
Accessor for all normal components - k1 - k12.
std::vector< G4double > SkewComponents() const
Accessor for all skew components - k1 - k12.
static G4double Unit(const G4String &key)
Access a unit factor for a given key.
Abstract base class that implements features common to all magnets.
Definition: BDSMagnet.hh:45
Information stored per model representing accelerator.
std::vector< BDSOutputROOTEventCollimatorInfo > collimatorInfo
Collimator information explicitly.
bool storeCollimatorInfo
Whether optional collimator information was stored.
std::vector< float > vkick
Vertical fractional momentum kick.
std::vector< float > eField
E field in V/m.
std::vector< std::string > collimatorBranchNamesUnique
Vector of all collimator branch names in event tree used to load data.
int nCollimators
Number of collimators in beam line.
virtual ~BDSOutputROOTEventModel()
Destructor.
int findNearestElement(const TVector3 &point) const
Find element index closest to point.
std::vector< float > bField
B field in T.
std::vector< float > hkick
Horizontal fractional momentum kick.
std::vector< int > collimatorIndices
std::vector< std::string > material
Material associated with element if any.
virtual void Fill(const std::vector< G4int > &collimatorIndicesIn={}, const std::map< G4String, G4int > &collimatorIndicesByNameIn={}, const std::vector< BDSOutputROOTEventCollimatorInfo > &collimatorInfoIn={}, const std::vector< G4String > &collimatorBranchNamesIn={}, const std::map< G4String, G4Transform3D > *scorerMeshPlacements=nullptr, const std::map< short int, G4String > *materialIDToNameUnique=nullptr, G4bool storeTrajectory=false)
Fill root output.
void Flush()
Initialise all members.
BDSOutputROOTEventModel()
Default constructor.
std::map< std::string, int > collimatorIndicesByName
Similar cache but by name of collimator as built by BDSIM.
std::vector< float > ks
Solenoid strength.
TRotation ConvertToROOT(const G4RotationMatrix *rm) const
Utility function.
static BDSSamplerRegistry * Instance()
Accessor for registry.
A holder for any placement offsets and rotations for a BDSAcceleratorComponent.
G4double GetYOffset() const
Accessor.
G4double GetXOffset() const
Accessor.
G4double GetTilt() const
Accessor.