BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSOutputStructures.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 "BDSAcceleratorModel.hh"
20#include "BDSBeamline.hh"
21#include "BDSBeamlineElement.hh"
22#include "BDSDebug.hh"
23#include "BDSHitEnergyDeposition.hh"
24#include "BDSEventInfo.hh"
25#include "BDSGlobalConstants.hh"
26#include "BDSOutputStructures.hh"
27#include "BDSOutputROOTEventAperture.hh"
28#include "BDSOutputROOTEventBeam.hh"
29#include "BDSOutputROOTEventCollimator.hh"
30#include "BDSOutputROOTEventCoords.hh"
31#include "BDSOutputROOTEventHeader.hh"
32#include "BDSOutputROOTEventHistograms.hh"
33#include "BDSOutputROOTEventInfo.hh"
34#include "BDSOutputROOTEventLoss.hh"
35#include "BDSOutputROOTEventLossWorld.hh"
36#include "BDSOutputROOTEventModel.hh"
37#include "BDSOutputROOTEventOptions.hh"
38#include "BDSOutputROOTEventRunInfo.hh"
39#include "BDSOutputROOTEventSampler.hh"
40#include "BDSOutputROOTEventSamplerC.hh"
41#include "BDSOutputROOTEventSamplerS.hh"
42#include "BDSOutputROOTEventTrajectory.hh"
43#include "BDSOutputROOTParticleData.hh"
44#include "BDSHitSampler.hh"
45#include "BDSSamplerRegistry.hh"
46#include "BDSTrajectoryPoint.hh"
47
48#include "globals.hh"
49#include "G4Material.hh"
50#include "G4MaterialTable.hh"
51
52#include "CLHEP/Units/SystemOfUnits.h"
53
54#include <algorithm>
55#include <map>
56#include <string>
57#include <vector>
58#include <utility>
59
61 nCollimators(0),
62 localSamplersInitialised(false),
63 localCollimatorsInitialised(false)
64{
65 G4bool storeCollimatorInfo = globals->StoreCollimatorInfo();
66 G4bool storeTurn = globals->StoreELossTurn();
67 G4bool storeLinks = globals->StoreELossLinks();
68 G4bool storeLocal = globals->StoreELossLocal();
69 G4bool storeGlobal = globals->StoreELossGlobal();
70 G4bool storeTime = globals->StoreELossTime();
71 G4bool storeStepLength = globals->StoreELossStepLength();
72 G4bool storePreStepKineticEnergy = globals->StoreELossPreStepKineticEnergy();
73 G4bool storeModelID = globals->StoreELossModelID();
74 G4bool storeELPhysics = globals->StoreELossPhysicsProcesses();
75 // store the model id if either modelID requested or store links
76 storeModelID = storeModelID || storeLinks;
77
78 particleDataOutput = new BDSOutputROOTParticleData();
79 headerOutput = new BDSOutputROOTEventHeader();
80 beamOutput = new BDSOutputROOTEventBeam();
81 optionsOutput = new BDSOutputROOTEventOptions();
82 modelOutput = new BDSOutputROOTEventModel(storeCollimatorInfo);
83
84 eLoss = new BDSOutputROOTEventLoss(storeTurn, storeLinks, storeModelID, storeLocal,
85 storeGlobal, storeTime, storeStepLength,
86 storePreStepKineticEnergy, storeELPhysics);
87 eLossVacuum = new BDSOutputROOTEventLoss(storeTurn, storeLinks, storeModelID, storeLocal,
88 storeGlobal, storeTime, storeStepLength,
89 storePreStepKineticEnergy, storeELPhysics);
90 eLossTunnel = new BDSOutputROOTEventLoss(storeTurn, storeLinks, storeModelID, storeLocal,
91 storeGlobal, storeTime, storeStepLength,
92 storePreStepKineticEnergy, storeELPhysics);
93 eLossWorld = new BDSOutputROOTEventLossWorld();
94 eLossWorldExit = new BDSOutputROOTEventLossWorld();
95 eLossWorldContents = new BDSOutputROOTEventLossWorld();
96
97 pFirstHit = new BDSOutputROOTEventLoss(true, true, true, true, true, true, false, true, true);
98 pLastHit = new BDSOutputROOTEventLoss(true, true, true, true, true, true, false, true, true);
99
100 apertureImpacts = new BDSOutputROOTEventAperture();
101
102 traj = new BDSOutputROOTEventTrajectory();
103 evtHistos = new BDSOutputROOTEventHistograms();
104 evtInfo = new BDSOutputROOTEventInfo();
105 runHistos = new BDSOutputROOTEventHistograms();
106 runInfo = new BDSOutputROOTEventRunInfo();
107
108#ifndef __ROOTDOUBLE__
109 primary = new BDSOutputROOTEventSampler<float>("Primary");
110#else
111 primary = new BDSOutputROOTEventSampler<double>("Primary");
112#endif
113 primaryGlobal = new BDSOutputROOTEventCoords();
114}
115
116BDSOutputStructures::~BDSOutputStructures()
117{
118 delete particleDataOutput;
119 delete headerOutput;
120 delete beamOutput;
121 delete optionsOutput;
122 delete modelOutput;
123 delete primaryGlobal;
124 delete eLoss;
125 delete eLossVacuum;
126 delete eLossTunnel;
127 delete eLossWorld;
128 delete eLossWorldExit;
129 delete eLossWorldContents;
130 delete pFirstHit;
131 delete pLastHit;
132 delete apertureImpacts;
133 delete traj;
134 delete evtHistos;
135 delete evtInfo;
136 delete runHistos;
137 delete runInfo;
138 for (auto sampler : samplerTrees)
139 {delete sampler;}
140 for (auto sampler : samplerCTrees)
141 {delete sampler;}
142 for (auto sampler : samplerSTrees)
143 {delete sampler;}
144 for (auto collimator : collimators)
145 {delete collimator;}
146 delete primary;
147}
148
149G4int BDSOutputStructures::Create1DHistogram(G4String name, G4String title,
150 G4int nbins, G4double xmin, G4double xmax)
151{
152 G4int result = evtHistos->Create1DHistogram(name, title, nbins, xmin, xmax);
153 // index from runHistos will be the same as used only through interfaces in this class
154 runHistos->Create1DHistogram(name, title, nbins, xmin, xmax);
155 return result;
156}
157
158G4int BDSOutputStructures::Create1DHistogram(G4String name, G4String title,
159 std::vector<double>& edges)
160{
161 G4int result = evtHistos->Create1DHistogram(name,title,edges);
162 runHistos->Create1DHistogram(name,title,edges);
163 return result;
164}
165
166G4int BDSOutputStructures::Create3DHistogram(G4String name, G4String title,
167 G4int nBinsX, G4double xMin, G4double xMax,
168 G4int nBinsY, G4double yMin, G4double yMax,
169 G4int nBinsZ, G4double zMin, G4double zMax)
170{
171 G4int result = evtHistos->Create3DHistogram(name, title,
172 nBinsX, xMin, xMax,
173 nBinsY, yMin, yMax,
174 nBinsZ, zMin, zMax);
175 // index from runHistos will be the same as used only through interfaces in this class
176 runHistos->Create3DHistogram(name, title,
177 nBinsX, xMin, xMax,
178 nBinsY, yMin, yMax,
179 nBinsZ, zMin, zMax);
180 return result;
181}
182
183G4int BDSOutputStructures::Create4DHistogram(const G4String& name,
184 const G4String& title,
185 const G4String& eScale,
186 const std::vector<double>& eBinsEdges,
187 G4int nBinsX, G4double xMin, G4double xMax,
188 G4int nBinsY, G4double yMin, G4double yMax,
189 G4int nBinsZ, G4double zMin, G4double zMax,
190 G4int nBinsE, G4double eMin, G4double eMax)
191{
192 G4int result = evtHistos->Create4DHistogram(name, title, eScale, eBinsEdges,
193 nBinsX, xMin, xMax,
194 nBinsY, yMin, yMax,
195 nBinsZ, zMin, zMax,
196 nBinsE, eMin, eMax);
197
198 runHistos->Create4DHistogram(name, title, eScale, eBinsEdges,
199 nBinsX, xMin, xMax,
200 nBinsY, yMin, yMax,
201 nBinsZ, zMin, zMax,
202 nBinsE, eMin, eMax);
203 return result;
204}
205
207{
209 {
210 auto samplerRegistry = BDSSamplerRegistry::Instance();
211 const auto sNames = samplerRegistry->GetUniqueNamesPlane();
212#ifdef USE_SIXTRACKLINK
213 // TODO hardcoded because of sixtrack dynamic buildup
214 // Sixtrack does lazy initialisation for collimators in link to Geant4 so we don't know
215 // a priori how many link elements there'll be. If we allow it to be dynamically built up
216 // we risk the vector expanding and moving in memory, therefore breaking all the && (address
217 // of pointer) links of TTree::SetBranchAddress in the output. Therefore, we reserve a size
218 // of 300 in the hope that this is less than the LHC 120 collimators for both beams. Ideally,
219 // the sixtrack interface should be rewritten so we know at construction time how many will
220 // be built.
221 samplerTrees.reserve(300);
222#else
223 samplerTrees.reserve(sNames.size());
224#endif
226 for (const auto& samplerName : sNames)
227 {
228#ifndef __ROOTDOUBLE__
230#else
232#endif
233 samplerTrees.push_back(res);
234 samplerNames.push_back(samplerName);
235 }
236 const auto planeIDs = samplerRegistry->GetSamplerIDsPlane();
237 G4int i = 0;
238 for (const auto& ID : planeIDs)
239 {samplerIDToIndexPlane[ID] = i; i++;}
240
241 // cylindrical samplers
242 const auto scNames = samplerRegistry->GetUniqueNamesCylinder();
243 samplerCTrees.reserve(scNames.size());
244 for (const auto& samplerName : scNames)
245 {
246 samplerCTrees.emplace_back(new BDSOutputROOTEventSamplerC(samplerName));
247 samplerCNames.emplace_back(samplerName);
248 }
249 const auto cylinderIDs = samplerRegistry->GetSamplerIDsCylinder();
250 i = 0;
251 for (const auto& ID : cylinderIDs)
252 {samplerIDToIndexCylinder[ID] = i; i++;}
253
254 // spherical samplers
255 const auto ssNames = samplerRegistry->GetUniqueNamesSphere();
256 samplerSTrees.reserve(ssNames.size());
257 for (const auto& samplerName : ssNames)
258 {
259 samplerSTrees.emplace_back(new BDSOutputROOTEventSamplerS(samplerName));
260 samplerSNames.emplace_back(samplerName);
261 }
262 const auto sphereIDs = samplerRegistry->GetSamplerIDsSphere();
263 i = 0;
264 for (const auto& ID : sphereIDs)
265 {samplerIDToIndexSphere[ID] = i; i++;}
266 }
267}
268
270{
271 materialToID.clear();
272 materialIDToNameUnique.clear();
273
274 const auto materialTable = G4Material::GetMaterialTable(); // should be an std::vector<G4Material*>*
275
276 // It's totally permitted to use degenerate material names as the geometry is done by pointer
277 // We need a way to sort the materials for a given input irrespective of pointer or memory
278 // location so the result is the same for multiple runs of bdsim.
279 // Use a pair of <name, density>. A c++ map will be internally sorted by keys and the various
280 // comparison operators are defined by pairs in <utility>.
281 // Once sorted, by a map, we then loop over that map and generate integer IDs for each
282 // material
283 // This is a little overkill really as we ensure in BDSMaterials we don't make materials
284 // with degenerate names and ultimately, we can't define degenerate materials in GMAD so
285 // this shouldn't happen. Perhaps it could from GDML.
286 std::map<std::pair<G4String, G4double>, G4Material*> sortingMap;
287 std::map<G4String, int> nameCount;
288 std::map<G4Material*, G4String> matToUniqueName;
289 for (const auto& mat : *materialTable)
290 {
291 G4String matName = mat->GetName();
292 sortingMap[std::make_pair(matName, mat->GetDensity())] = mat;
293
294 auto search = nameCount.find(matName);
295 if (search != nameCount.end())
296 {
297 search->second += 1;
298 matToUniqueName[mat] = matName + std::to_string(search->second);
299 }
300 else
301 {
302 nameCount[matName] = 0;
303 matToUniqueName[mat] = matName;
304 }
305 }
306 short int i = 0;
307 for (const auto& kv : sortingMap)
308 {
309 materialToID[kv.second] = i;
310 materialIDToNameUnique[i] = matToUniqueName[kv.second];
311 i++;
312 }
313}
314
316{
317 G4int result = 0;
318 for (auto const& samplerName : BDSSamplerRegistry::Instance()->GetUniqueNames())
319 {// only put it in if it doesn't exist already
320 if (std::find(samplerNames.begin(), samplerNames.end(), samplerName) == samplerNames.end())
321 {
322 result++;
323#ifndef __ROOTDOUBLE__
325#else
327#endif
328 samplerTrees.push_back(res);
329 samplerNames.push_back(samplerName);
330 }
331 }
333 return result;
334}
335
337{
338 const G4String collimatorPrefix = "COLL_";
339 const BDSBeamline* flatBeamline = BDSAcceleratorModel::Instance()->BeamlineMain();
341 nCollimators = (G4int)collimatorIndices.size();
342
343 for (auto index : collimatorIndices)
344 {
345 // prepare output structure name
346 const BDSBeamlineElement* el = flatBeamline->at(index);
347 // use the 'placement' name for a unique name (with copy number included)
348 G4String collimatorName = collimatorPrefix + el->GetPlacementName();
349 collimatorNames.push_back(collimatorName);
350 collimatorIndicesByName[el->GetName()] = index;
352
354 info.Fill(el);
355 collimatorInfo.push_back(info);
356
357 // cache difference in apertures for efficient interpolation and avoid
358 // repeated calculation. not required in info for output though.
359 G4double xDiff = info.xSizeOut - info.xSizeIn;
360 G4double yDiff = info.ySizeOut - info.ySizeIn;
361 collimatorDifferences.emplace_back(xDiff, yDiff); // construct in place
362 }
363}
364
366{
368 {
370 for (int i = 0; i < (int)collimatorIndices.size(); i++)
371 {collimators.push_back(new BDSOutputROOTEventCollimator());}
372 }
373}
374
376{
378}
379
381{
382 headerOutput->Flush();
383}
384
386{
388}
389
391{
392 *beamOutput = BDSOutputROOTEventBeam(); // default
393}
394
396{
398}
399
401{
402 primary->Flush();
403 for (auto sampler : samplerTrees)
404 {sampler->Flush();}
405 for (auto sampler : samplerCTrees)
406 {sampler->Flush();}
407 for (auto sampler : samplerSTrees)
408 {sampler->Flush();}
409 for (auto collimator : collimators)
410 {collimator->Flush();}
411 primaryGlobal->Flush();
412 eLoss->Flush();
413 eLossVacuum->Flush();
414 eLossTunnel->Flush();
415 eLossWorld->Flush();
416 eLossWorldExit->Flush();
417 eLossWorldContents->Flush();
418 pFirstHit->Flush();
419 pLastHit->Flush();
420 apertureImpacts->Flush();
421 traj->Flush();
422 evtHistos->Flush();
423 evtInfo->Flush();
424}
425
427{
428 runInfo->Flush();
429}
const BDSBeamline * BeamlineMain() const
Accessor.
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4String GetPlacementName() const
Accessor.
G4String GetName() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
const BDSBeamlineElement * at(int iElement) const
Return a reference to the element at i.
Definition: BDSBeamline.hh:120
std::vector< G4int > GetIndicesOfCollimators() const
Return indices in order of ecol, rcol, jcol and crystalcol elements.
Definition: BDSBeamline.cc:954
A class that holds global options and constants.
Data stored for energy deposition hits per event.
Class to store all beam options for a BDSIM run.
Data stored for each collimator in the model.
Data stored for each collimator per event.
Information about the software and the file.
Holder for a set of histograms to be stored.
virtual void Flush()
Flush the contents.
Information pertaining to an individual event.
Data stored for world hits per event.
Data stored for energy deposition hits per event.
Information stored per model representing accelerator.
void Flush()
Initialise all members.
Class to store all options for a BDSIM run.
Information pertaining to a run.
Information stored per cylindrical sampler per event.
Information stored per spherical sampler per event.
virtual void Flush()
Clean Sampler.
Structure to record a trajectory.
void Flush()
add comment to avoid warning (no need to make persistent, see issue #191)
Geant4 particle data for particles used in simulation.
virtual void Flush()
Clear maps.
void InitialiseCollimators()
Construct collimtors.
BDSOutputROOTEventModel * modelOutput
Model output.
G4int Create1DHistogram(G4String name, G4String title, G4int nbins, G4double xmin, G4double xmax)
Create histograms for both evtHistos and runHistos. Return index from evtHistos.
std::vector< G4String > collimatorNames
Names of collimators in output structures.
std::vector< BDSOutputROOTEventCollimator * > collimators
Collimator output structures.
std::vector< std::string > samplerNames
Sampler names to use.
BDSOutputROOTEventSampler< double > * primary
Primary sampler structure.
BDSOutputROOTEventRunInfo * runInfo
Run information.
void InitialiseMaterialMap()
Construct a map of material pointer to integer ID and name.
void ClearStructuresOptions()
Clear the local options structure.
BDSOutputROOTEventHistograms * runHistos
Run level histograms.
BDSOutputROOTEventLossWorld * eLossWorldExit
World exit hits.
void ClearStructuresRunLevel()
Clear the local structures in this class in preparation for a new run.
BDSOutputStructures()=delete
Unused default constructors.
void InitialiseSamplers()
Construct samplers.
void ClearStructuresBeam()
Clear the local beam structure.
std::map< G4String, G4int > collimatorIndicesByName
Indices mapped to their name.
BDSOutputROOTEventLossWorld * eLossWorldContents
Externally supplied world contents hits.
G4int nCollimators
Number of collimators in beam line.
std::vector< BDSOutputROOTEventSampler< double > * > samplerTrees
Sampler structures.
std::vector< std::pair< G4double, G4double > > collimatorDifferences
void ClearStructuresModel()
Clear the local model structure.
void ClearStructuresHeader()
Clear the local header structure.
G4int Create4DHistogram(const G4String &name, const G4String &title, const G4String &eScale, const std::vector< double > &eBinsEdges, G4int nBinsX, G4double xMin, G4double xMax, G4int nBinsY, G4double yMin, G4double yMax, G4int nBinsZ, G4double zMin, G4double zMax, G4int nBinsE, G4double eMin, G4double eMax)
Create histograms for both evtHistos and runHistos. Return index from evtHistos.
BDSOutputROOTEventHeader * headerOutput
Information about the file.
BDSOutputROOTEventLossWorld * eLossWorld
World energy deposition.
BDSOutputROOTEventLoss * eLossVacuum
General energy deposition.
std::vector< G4int > collimatorIndices
Indices in beam line that are collimators.
BDSOutputROOTEventInfo * evtInfo
Event information.
BDSOutputROOTEventLoss * pLastHit
Primary loss point.
BDSOutputROOTEventHistograms * evtHistos
Event level histograms.
BDSOutputROOTEventLoss * eLossTunnel
Tunnel energy deposition.
G4int UpdateSamplerStructures()
Interface to allow setting up samplers later for dynamic geometry construction a la SixTrack....
G4int Create3DHistogram(G4String name, G4String title, G4int nBinsX, G4double xMin, G4double xMax, G4int nBinsY, G4double yMin, G4double yMax, G4int nBinsZ, G4double zMin, G4double zMax)
Create histograms for both evtHistos and runHistos. Return index from evtHistos.
void ClearStructuresEventLevel()
Clear the local structures in this class in preparation for a new event.
BDSOutputROOTEventLoss * pFirstHit
Primary hit point.
BDSOutputROOTEventTrajectory * traj
Trajectories.
BDSOutputROOTEventBeam * beamOutput
Beam output.
BDSOutputROOTEventLoss * eLoss
General energy deposition.
std::vector< BDSOutputROOTEventCollimatorInfo > collimatorInfo
void ClearStructuresParticleData()
Clear the local particle data structure.
G4bool localCollimatorsInitialised
Whether we've setup the member vector of collimators. Similarly to localSamplersInitialised.
BDSOutputROOTEventAperture * apertureImpacts
Impacts on the aperture.
BDSOutputROOTParticleData * particleDataOutput
Geant4 information / particle tables.
BDSOutputROOTEventOptions * optionsOutput
Options output.
static BDSSamplerRegistry * Instance()
Accessor for registry.