BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSOutput.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 "BDSBLMRegistry.hh"
23#include "BDSDebug.hh"
24#include "BDSEventInfo.hh"
25#include "BDSException.hh"
26#include "BDSGlobalConstants.hh"
27#include "BDSHistBinMapper.hh"
28#include "BDSHitApertureImpact.hh"
29#include "BDSHitCollimator.hh"
30#include "BDSHitEnergyDeposition.hh"
31#include "BDSHitEnergyDepositionGlobal.hh"
32#include "BDSHitSampler.hh"
33#include "BDSHitSamplerCylinder.hh"
34#include "BDSHitSamplerSphere.hh"
35#include "BDSHitSamplerLink.hh"
36#include "BDSOutput.hh"
37#include "BDSOutputROOTEventAperture.hh"
38#include "BDSOutputROOTEventBeam.hh"
39#include "BDSOutputROOTEventCollimator.hh"
40#include "BDSOutputROOTEventCollimatorInfo.hh"
41#include "BDSOutputROOTEventCoords.hh"
42#include "BDSOutputROOTEventLossWorld.hh"
43#include "BDSOutputROOTEventHeader.hh"
44#include "BDSOutputROOTEventHistograms.hh"
45#include "BDSOutputROOTEventInfo.hh"
46#include "BDSOutputROOTEventLoss.hh"
47#include "BDSOutputROOTEventModel.hh"
48#include "BDSOutputROOTEventOptions.hh"
49#include "BDSOutputROOTEventRunInfo.hh"
50#include "BDSOutputROOTEventSampler.hh"
51#include "BDSOutputROOTEventSamplerC.hh"
52#include "BDSOutputROOTEventSamplerS.hh"
53#include "BDSOutputROOTEventTrajectory.hh"
54#include "BDSOutputROOTParticleData.hh"
55#include "BDSParticleDefinition.hh"
56#include "BDSPrimaryVertexInformation.hh"
57#include "BDSPrimaryVertexInformationV.hh"
58#include "BDSScorerHistogramDef.hh"
59#include "BDSSDManager.hh"
60#include "BDSStackingAction.hh"
61#include "BDSTrajectoriesToStore.hh"
62#include "BDSTrajectoryPoint.hh"
63#include "BDSTrajectoryPointHit.hh"
64#include "BDSUtilities.hh"
65
66#include "globals.hh"
67#include "G4ParticleDefinition.hh"
68#include "G4PrimaryParticle.hh"
69#include "G4PrimaryVertex.hh"
70#include "G4THitsMap.hh"
71#include "G4Version.hh"
72
73#include "TH1D.h"
74#include "TH3D.h"
75
76#include "parser/beamBase.h"
77#include "parser/optionsBase.h"
78
79#include <algorithm>
80#include <cmath>
81#include <map>
82#include <ostream>
83#include <set>
84#include <utility>
85#include <vector>
86
87#include "CLHEP/Units/SystemOfUnits.h"
88
89const std::set<G4String> BDSOutput::protectedNames = {
90 "Event", "Histos", "Info", "Primary", "PrimaryGlobal",
91 "Eloss", "ElossVacuum", "ElossTunnel", "ElossWorld", "ElossWorldExit",
92 "ElossWorldContents",
93 "PrimaryFirstHit", "PrimaryLastHit", "Trajectory", "ApertureImpacts"
94};
95
96BDSOutput::BDSOutput(const G4String& baseFileNameIn,
97 const G4String& fileExtensionIn,
98 G4int fileNumberOffset):
100 baseFileName(baseFileNameIn),
101 fileExtension(fileExtensionIn),
102 outputFileNumber(fileNumberOffset),
103 sMaxHistograms(0),
104 nbins(0),
105 energyDeposited(0),
106 energyDepositedVacuum(0),
107 energyDepositedWorld(0),
108 energyDepositedWorldContents(0),
109 energyDepositedTunnel(0),
110 energyImpactingAperture(0),
111 energyImpactingApertureKinetic(0),
112 energyWorldExit(0),
113 energyWorldExitKinetic(0),
114 nCollimatorsInteracted(0)
115{
117 numberEventPerFile = g->NumberOfEventsPerNtuple();
118 useScoringMap = g->UseScoringMap();
119
120 storeApertureImpacts = g->StoreApertureImpacts();
121 storeApertureImpactsHistograms = g->StoreApertureImpactsHistograms();
122 storeCollimatorInfo = g->StoreCollimatorInfo();
123 storeCollimatorHitsLinks = g->StoreCollimatorHitsLinks();
124 storeCollimatorHitsIons = g->StoreCollimatorHitsIons();
125 // store primary hits if ion hits or links hits are turned on
126 storeCollimatorHits = g->StoreCollimatorHits() || storeCollimatorHitsLinks || storeCollimatorHitsIons || g->StoreCollimatorHitsAll();
127
129
130 storeELoss = g->StoreELoss();
131 // store histograms if storing general energy deposition as negligible in size
132 storeELossHistograms = g->StoreELossHistograms() || storeELoss;
133 storeELossTunnel = g->StoreELossTunnel();
134 storeELossTunnelHistograms = g->StoreELossTunnelHistograms() || storeELossTunnel;
135 storeELossVacuum = g->StoreELossVacuum();
136 storeELossVacuumHistograms = g->StoreELossVacuumHistograms() || storeELossVacuum;
137 storeELossWorld = g->StoreELossWorld();
138 storeELossWorldContents = g->StoreELossWorldContents() || g->UseImportanceSampling();
139 storeParticleData = g->StoreParticleData();
140 storeModel = g->StoreModel();
141 storePrimaries = g->StorePrimaries();
142 storePrimaryHistograms = g->StorePrimaryHistograms();
143 storeSamplerPolarCoords = g->StoreSamplerPolarCoords();
144 storeSamplerCharge = g->StoreSamplerCharge();
145 storeSamplerKineticEnergy = g->StoreSamplerKineticEnergy();
146 storeSamplerMass = g->StoreSamplerMass();
147 storeSamplerRigidity = g->StoreSamplerRigidity();
148 storeSamplerIon = g->StoreSamplerIon();
149 storeTrajectory = g->StoreTrajectory();
150 storeTrajectoryStepPoints = g->StoreTrajectoryStepPoints();
151 storeTrajectoryStepPointLast = g->StoreTrajectoryStepPointLast();
153
154 // easy option for everything - overwrite bools we've just set individually
155 if (g->StoreSamplerAll())
156 {
158 storeSamplerCharge = true;
160 storeSamplerMass = true;
162 storeSamplerIon = true;
163 }
164}
165
167{
169 {
170 PrepareCollimatorInformation(); // prepare names, offsets and indices
171 InitialiseCollimators(); // allocate local objects
172 }
176}
177
179{
180 headerOutput->Flush();
181 headerOutput->Fill(); // updates time stamp
182 WriteHeader();
184}
185
186void BDSOutput::FillParticleData(G4bool writeIons)
187{
188 // always prepare particle data and link to other classes, but optionally fill it
190 particleDataOutput->Fill(writeIons);
191
192#ifdef __ROOTDOUBLE__
194#else
196#endif
197 BDSOutputROOTEventCollimator::particleTable = particleDataOutput;
198 BDSOutputROOTEventAperture::particleTable = particleDataOutput;
199
202}
203
205{
207 WriteBeam();
209}
210
212{
214 WriteOptions();
216}
217
219{
220 if (storeModel)
221 {
222 const auto& smpm = BDSAcceleratorModel::Instance()->ScorerMeshPlacementsMap();
227 &smpm,
228 &materialIDToNameUnique,
230 WriteModel();
232 }
233}
234
235void BDSOutput::FillPrimary(const G4PrimaryVertex* vertex,
236 const G4int turnsTaken)
237{
238 const G4VUserPrimaryVertexInformation* vertexInfo = vertex->GetUserInformation();
239 if (const auto vertexInfoBDS = dynamic_cast<const BDSPrimaryVertexInformation*>(vertexInfo))
240 {
241 primary->Fill(vertexInfoBDS->primaryVertex.local,
242 vertexInfoBDS->momentum,
243 vertexInfoBDS->charge,
244 vertexInfoBDS->pdgID,
245 turnsTaken,
246 vertexInfoBDS->primaryVertex.beamlineIndex,
247 vertexInfoBDS->nElectrons,
248 vertexInfoBDS->mass,
249 vertexInfoBDS->rigidity);
250 primaryGlobal->Fill(vertexInfoBDS->primaryVertex.global);
251 }
252 else if (const auto vertexInfoBDSV = dynamic_cast<const BDSPrimaryVertexInformationV*>(vertexInfo))
253 {// vector version - multiple primaries at primary vertex
254 primary->Fill(vertexInfoBDSV,
255 turnsTaken);
256 primaryGlobal->Fill(vertexInfoBDSV);
257 }
258 auto nextLinkedVertex = vertex->GetNext();
259 if (nextLinkedVertex)
260 {FillPrimary(nextLinkedVertex, turnsTaken);}
261}
262
264 const BDSParticleDefinition* particle)
265{
266 G4bool isIon = particle->IsAnIon();
267 G4int ionA = 0;
268 G4int ionZ = 0;
269 if (isIon)
270 {// fill primary ion info correctly when we have no particle table available
271 ionA = particle->IonDefinition()->A();
272 ionZ = particle->IonDefinition()->Z();
273 }
274 primary->Fill(coords.local,
275 particle->Momentum(),
276 particle->Charge(),
277 particle->PDGID(),
278 0, 0,
279 particle->NElectrons(),
280 particle->Mass(),
281 particle->BRho(),
282 true,
283 &isIon, &ionA, &ionZ);
284 primaryGlobal->Fill(coords.global);
287}
288
290 const G4PrimaryVertex* vertex,
291 const std::vector<BDSHitsCollectionSampler*>& samplerHitsPlane,
292 const std::vector<BDSHitsCollectionSamplerCylinder*>& samplerHitsCylinder,
293 const std::vector<BDSHitsCollectionSamplerSphere*>& samplerHitsSphere,
294 const BDSHitsCollectionSamplerLink* samplerHitsLink,
295 const BDSHitsCollectionEnergyDeposition* energyLoss,
296 const BDSHitsCollectionEnergyDeposition* energyLossFull,
297 const BDSHitsCollectionEnergyDeposition* energyLossVacuum,
298 const BDSHitsCollectionEnergyDeposition* energyLossTunnel,
299 const BDSHitsCollectionEnergyDepositionGlobal* energyLossWorld,
300 const BDSHitsCollectionEnergyDepositionGlobal* energyLossWorldContents,
301 const BDSHitsCollectionEnergyDepositionGlobal* worldExitHits,
302 const std::vector<const BDSTrajectoryPointHit*>& primaryHits,
303 const std::vector<const BDSTrajectoryPointHit*>& primaryLosses,
304 const BDSTrajectoriesToStore* trajectories,
305 const BDSHitsCollectionCollimator* collimatorHits,
306 const BDSHitsCollectionApertureImpacts* apertureImpactHits,
307 const std::map<G4String, G4THitsMap<G4double>*>& scorerHits,
308 const G4int turnsTaken)
309{
310 // Clear integrals in this class -> here instead of BDSOutputStructures as
311 // looped over here -> do only once as expensive as lots of hits
312 energyDeposited = 0;
319 energyWorldExit = 0;
322
323 if (vertex && storePrimaries)
324 {FillPrimary(vertex, turnsTaken);}
325 FillSamplerHitsVector(samplerHitsPlane);
326 FillSamplerCylinderHitsVector(samplerHitsCylinder);
327 FillSamplerSphereHitsVector(samplerHitsSphere);
328 if (samplerHitsLink)
329 {FillSamplerHitsLink(samplerHitsLink);}
330 if (energyLoss)
331 {FillEnergyLoss(energyLoss, BDSOutput::LossType::energy);}
332 if (energyLossFull)
333 {FillEnergyLoss(energyLossFull, BDSOutput::LossType::energy);}
334 if (energyLossVacuum)
335 {FillEnergyLoss(energyLossVacuum, BDSOutput::LossType::vacuum);}
336 if (energyLossTunnel)
337 {FillEnergyLoss(energyLossTunnel, BDSOutput::LossType::tunnel);}
338 if (energyLossWorld)
339 {FillEnergyLoss(energyLossWorld, BDSOutput::LossType::world);}
340 if (worldExitHits)
341 {FillEnergyLoss(worldExitHits, BDSOutput::LossType::worldexit);}
342 if (energyLossWorldContents)
343 {FillEnergyLoss(energyLossWorldContents, BDSOutput::LossType::worldcontents);}
344 FillPrimaryHit(primaryHits);
345 FillPrimaryLoss(primaryLosses);
346 if (trajectories)
347 {FillTrajectories(trajectories);}
348 if (collimatorHits)
349 {FillCollimatorHits(collimatorHits, primaryLosses);}
350 if (apertureImpacts)
351 {FillApertureImpacts(apertureImpactHits);}
352 FillScorerHits(scorerHits); // map always exists
353
354 // we do this after energy loss and collimator hits as the energy loss
355 // is integrated for putting in event info and the number of collimators
356 // interacted with counted
357 if (info)
358 {FillEventInfo(info);}
359
362}
363
365{
366 CloseFile();
367 NewFile();
369}
370
372{
373 FillRunInfo(info);
376}
377
378G4bool BDSOutput::InvalidSamplerName(const G4String& samplerName)
379{
380 return protectedNames.find(samplerName) != protectedNames.end();
381}
382
383void BDSOutput::PrintProtectedNames(std::ostream& out)
384{
385 out << "Protected names for output " << G4endl;
386 for (const auto& key : protectedNames)
387 {out << "\"" << key << "\"" << G4endl;}
388}
389
391{
393 const BDSGlobalConstants* globalConstants = BDSGlobalConstants::Instance();
394
395 // Base root file name
396 G4String newFileName = baseFileName;
397
398 // if more than one file add number (starting at 0)
399 // of numberEventPerFile is specified and the number already generated exceeds that
400 if (numberEventPerFile > 0 && globalConstants->NGenerate() > numberEventPerFile)
401 {newFileName += "_" + std::to_string(outputFileNumber);} // note underscore
402
403 // policy: overwrite if output filename specifically set, otherwise increase number
404 // always check in interactive mode
405 if (!globalConstants->OutputFileNameSet() || !globalConstants->Batch())
406 {// check if file exists
407 G4String original = newFileName; // could have nper file number suffix too
408 G4int nTimeAppended = 1;
409 while (BDS::FileExists(newFileName + fileExtension)) // always test with extension
410 {// if exists increment suffix integer
411 newFileName = original + "-" + std::to_string(nTimeAppended);
412 nTimeAppended +=1;
413 }
414 }
415
416 // add extension now we've got the base part fixed
417 newFileName += fileExtension;
418
419 G4cout << __METHOD_NAME__ << "Setting up new file: " << newFileName << G4endl;
420
421 return newFileName;
422}
423
425{
426 // rounding up so last bin definitely covers smax
427 // (max - min) / bin width -> min = 0 here.
428 const G4double binWidth = BDSGlobalConstants::Instance()->ELossHistoBinWidth();
429 const BDSBeamline* flatBeamline = BDSAcceleratorModel::Instance()->BeamlineMain();
430 if (flatBeamline)
431 {// don't access a nullptr
432 if (!flatBeamline->empty())
433 {
434 G4double sMax = flatBeamline->GetLastItem()->GetSPositionEnd();
435 nbins = (int) std::ceil(sMax / binWidth); // round up to integer # of bins
436 }
437 }
438 else
439 {nbins = 1;} // can happen for generate primaries only
440
441 if (nbins == 0)
442 {nbins = 1;}
443
444 sMaxHistograms = nbins * binWidth;
445}
446
448{
449 // construct output histograms
450 // calculate histogram dimensions
452 const G4double smin = 0.0;
453 const G4double smax = sMaxHistograms / CLHEP::m;
454#ifdef BDSDEBUG
455 G4cout << __METHOD_NAME__ << "histogram parameters calculated to be: " << G4endl;
456 G4cout << "s minimum: " << smin << " m" << G4endl;
457 G4cout << "s maximum: " << smax << " m" << G4endl;
458 G4cout << "# of bins: " << nbins << G4endl;
459#endif
460 // create the histograms
462 {
463 histIndices1D["Phits"] = Create1DHistogram("PhitsHisto","Primary Hits",nbins,smin,smax);
464 histIndices1D["Ploss"] = Create1DHistogram("PlossHisto","Primary Loss",nbins,smin,smax);
465 }
467 {histIndices1D["Eloss"] = Create1DHistogram("ElossHisto","Energy Loss",nbins,smin,smax);}
468 // prepare bin edges for a by-element histogram
469 std::vector<G4double> binedges;
470 const BDSBeamline* flatBeamline = BDSAcceleratorModel::Instance()->BeamlineMain();
471 if (flatBeamline) // can be nullptr in case of generate primaries only
472 {binedges = flatBeamline->GetEdgeSPositions();}
473 else
474 {binedges = {0,1};}
475 // create per element ("pe") bin width histograms
477 {
478 histIndices1D["PhitsPE"] = Create1DHistogram("PhitsPEHisto",
479 "Primary Hits per Element",
480 binedges);
481 histIndices1D["PlossPE"] = Create1DHistogram("PlossPEHisto",
482 "Primary Loss per Element",
483 binedges);
484 }
486 {
487 histIndices1D["ElossPE"] = Create1DHistogram("ElossPEHisto",
488 "Energy Loss per Element" ,
489 binedges);
490 }
492 {
493 histIndices1D["ElossVacuum"] = Create1DHistogram("ElossVacuumHisto",
494 "Energy Loss in Vacuum",
495 nbins,smin,smax);
496 histIndices1D["ElossVacuumPE"] = Create1DHistogram("ElossVaccumPEHisto",
497 "Energy Loss in Vacuum per Element" ,
498 binedges);
499 }
500
502 {
503 histIndices1D["PFirstAI"] = Create1DHistogram("PFirstAIHisto",
504 "Primary aperture impacts",
505 nbins, smin, smax);
506 }
507
508 // only create tunnel histograms if we build the tunnel
509 const BDSBeamline* tunnelBeamline = BDSAcceleratorModel::Instance()->TunnelBeamline();
510 if (!tunnelBeamline)
511 {
512 storeELossTunnel = false;
514 }
515 if (storeELossTunnelHistograms && tunnelBeamline)
516 {
517 binedges = tunnelBeamline->GetEdgeSPositions();
518 histIndices1D["ElossTunnel"] = Create1DHistogram("ElossTunnelHisto",
519 "Energy Loss in Tunnel",
520 nbins, smin,smax);
521 histIndices1D["ElossTunnelPE"] = Create1DHistogram("ElossTunnelPEHisto",
522 "Energy Loss in Tunnel per Element",
523 binedges);
524 }
525
527 {
528 std::vector<G4String> collHistNames = {"CollPhitsPE",
529 "CollPlossPE",
530 "CollElossPE",
531 "CollPInteractedPE"};
532 std::vector<G4String> collHistDesciptions = {"Primary Hits per Coll",
533 "Primary Loss per Coll",
534 "Energy Loss per Collimator",
535 "Primary Interacted per Collimator"};
536 for (G4int i = 0; i < (G4int)collHistNames.size(); i++)
537 {
538 histIndices1D[collHistNames[i]] = Create1DHistogram(collHistNames[i],
539 collHistDesciptions[i],
541 }
542 }
543
544 // one unique 'scorer' - single 3d histogram 3d
546 {
548 if (!BDS::IsFinite(g->XMax() - g->XMin()))
549 {throw BDSException(__METHOD_NAME__, "0 width in general 3D scoring histogram in x dimension - check options, xmin and xmax");}
550 if (!BDS::IsFinite(g->YMax() - g->YMin()))
551 {throw BDSException(__METHOD_NAME__, "0 width in general 3D scoring histogram in y dimension - check options, ymin and ymax");}
552 if (!BDS::IsFinite(g->ZMax() - g->ZMin()))
553 {throw BDSException(__METHOD_NAME__, "0 width in general 3D scoring histogram in z dimension - check options, zmin and zmax");}
554 if (g->NBinsX() <= 0)
555 {throw BDSException(__METHOD_NAME__, "invalid number of bins in x dimension of 3D scoring histogram - check option, nbinsx");}
556 if (g->NBinsY() <= 0)
557 {throw BDSException(__METHOD_NAME__, "invalid number of bins in y dimension of 3D scoring histogram - check option, nbinsx");}
558 if (g->NBinsZ() <= 0)
559 {throw BDSException(__METHOD_NAME__, "invalid number of bins in z dimension of 3D scoring histogram - check option, nbinsx");}
560
561 G4int scInd = Create3DHistogram("ScoringMap", "Energy Deposition",
562 g->NBinsX(), g->XMin()/CLHEP::m, g->XMax()/CLHEP::m,
563 g->NBinsY(), g->YMin()/CLHEP::m, g->YMax()/CLHEP::m,
564 g->NBinsZ(), g->ZMin()/CLHEP::m, g->ZMax()/CLHEP::m);
565 histIndices3D["ScoringMap"] = scInd;
566 }
567
568 // scoring maps
569 const std::map<G4String, BDSScorerHistogramDef> scorerHistogramDefs = BDSAcceleratorModel::Instance()->ScorerHistogramDefinitionsMap();
570 if (!scorerHistogramDefs.empty())
571 {
572 for (const auto& nameDef : scorerHistogramDefs)
573 {
574 const auto def = nameDef.second;
575
576 // use safe output name without any slashes in the name
577 G4int histID = -1;
578
579 if (def.nBinsE <=1)
580 {
581
582 if (def.geometryType == "box"){
583 histID = Create3DHistogram(def.outputName, def.outputName,
584 def.nBinsX, def.xLow/CLHEP::m, def.xHigh/CLHEP::m,
585 def.nBinsY, def.yLow/CLHEP::m, def.yHigh/CLHEP::m,
586 def.nBinsZ, def.zLow/CLHEP::m, def.zHigh/CLHEP::m);}
587 else if (def.geometryType == "cylindrical"){
588 histID = Create3DHistogram(def.outputName, def.outputName,
589 def.nBinsZ, def.zLow/CLHEP::m, def.zHigh/CLHEP::m,
590 def.nBinsPhi, 0, 2*M_PI,
591 def.nBinsR, def.rLow/CLHEP::m, def.rHigh/CLHEP::m);}
592
593 histIndices3D[def.uniqueName] = histID;
594 histIndexToUnits3D[histID] = def.primitiveScorerUnitValue;
595 // avoid using [] operator for map as we have no default constructor for BDSHistBinMapper3D
596 }
597 else
598 {
599
600 if (def.geometryType == "box"){
601 histID = Create4DHistogram(def.outputName+"-"+def.eScale,def.outputName,def.eScale,def.eBinsEdges,
602 def.nBinsX, def.xLow/CLHEP::m, def.xHigh/CLHEP::m,
603 def.nBinsY, def.yLow/CLHEP::m, def.yHigh/CLHEP::m,
604 def.nBinsZ, def.zLow/CLHEP::m, def.zHigh/CLHEP::m,
605 def.nBinsE, def.eLow/CLHEP::GeV, def.eHigh/CLHEP::GeV);}
606 else if (def.geometryType == "cylindrical"){
607 histID = Create4DHistogram(def.outputName+"-"+def.eScale, def.outputName, def.eScale,def.eBinsEdges,
608 def.nBinsZ, def.zLow/CLHEP::m, def.zHigh/CLHEP::m,
609 def.nBinsPhi, 0, 2*M_PI,
610 def.nBinsR, def.rLow/CLHEP::m, def.rHigh/CLHEP::m,
611 def.nBinsE, def.eLow/CLHEP::GeV, def.eHigh/CLHEP::GeV);}
612
613 histIndices4D[def.uniqueName] = histID;
614 histIndexToUnits4D[histID] = def.primitiveScorerUnitValue;
615 }
616 scorerCoordinateMaps.insert(std::make_pair(def.uniqueName, def.coordinateMapper));
617 }
618 }
619
620 G4int nBLMs = BDSBLMRegistry::Instance()->NBLMs();
621 if (nBLMs > 0)
622 {
623 // the same primitive scorers for BLMs may be used in multiple SDs (each representing
624 // a unique combination of primitive scorers. However, we want 1 histogram per primitive
625 // scorer for all BLMs. We want to fill the same scoring quantity into the one histogram
626 // even from different collections ("SD/PS"). Determine 'set' of histogram names, which is
627 // set of primitive scorers use for BLMs.
628 // note there may be scorers from general 3d meshes and not just blms
629 std::vector<G4String> psnamesc = BDSSDManager::Instance()->PrimitiveScorerNamesComplete();
630 std::vector<G4String> psnames = BDSSDManager::Instance()->PrimitiveScorerNames();
631 std::map<G4String, G4double> scorerUnits = BDSSDManager::Instance()->PrimitiveScorerUnits();
632 std::set<G4String> blmHistoNames;
633 std::map<G4String, G4String> psFullNameToPS;
634 for (const auto& scorerNameComplete : psnamesc)
635 {
636 if (BDS::StrContains(scorerNameComplete, "blm_"))
637 {
638 for (const auto& scorerName : psnames)
639 {
640 if (BDS::StrContains(scorerNameComplete, "/"+scorerName)) // only match end of full name with '/'
641 {
642 blmHistoNames.insert(scorerName);
643 psFullNameToPS[scorerNameComplete] = scorerName;
644 }
645 }
646 }
647 }
648
649 // make BLM histograms and map the full collection name to that histogram ID for easy filling
650 // at the end of event. Note, multiple collections may feed into the same histogram.
651 for (const auto &hn : blmHistoNames)
652 {
653 G4String blmHistName = "BLM_" + hn;
654 G4int hind = Create1DHistogram(blmHistName, blmHistName, nBLMs, 0, nBLMs);
655 histIndices1D[blmHistName] = hind;
656 histIndexToUnits1D[hind] = scorerUnits[hn];
657 for (const auto& kv : psFullNameToPS)
658 {
659 if (hn == kv.second)
660 {blmCollectionNameToHistogramID[kv.first] = hind;}
661 }
662 }
663 }
664}
665
667{
668 if (info)
669 {*evtInfo = *(info->GetInfo());}
679 G4double ek = BDSStackingAction::energyKilled / CLHEP::GeV;
680 evtInfo->energyKilled = ek;
687 + ek;
688
690}
691
692void BDSOutput::FillSamplerHitsVector(const std::vector<BDSHitsCollectionSampler*>& hits)
693{
694 for (const auto& hc : hits)
695 {
696 if (!hc)
697 {continue;} // could be nullptr
698 if (!(hc->entries() > 0))
699 {continue;}
700 for (int i = 0; i < (int) hc->entries(); i++)
701 {
702 const BDSHitSampler* hit = (*hc)[i];
703 G4int samplerID = hit->samplerID;
704 G4int samplerVectorIndex = samplerIDToIndexPlane[samplerID];
705 samplerTrees[samplerVectorIndex]->Fill(hit, storeSamplerMass, storeSamplerCharge,
708 }
709 }
710 // extra information - do only once at the end
711 if (storeSamplerIon)
712 {
713 for (auto& sampler : samplerTrees)
714 {sampler->FillIon();}
715 }
716}
717
718void BDSOutput::FillSamplerCylinderHitsVector(const std::vector<BDSHitsCollectionSamplerCylinder*>& hits)
719{
720 for (const auto& hc : hits)
721 {
722 if (!hc)
723 {continue;} // could be nullptr
724 if (!(hc->entries() > 0))
725 {continue;}
726 for (int i = 0; i < (int) hc->entries(); i++)
727 {
728 const BDSHitSamplerCylinder* hit = (*hc)[i];
729 G4int samplerID = hit->samplerID;
730 G4int samplerVectorIndex = samplerIDToIndexCylinder[samplerID];
731 samplerCTrees[samplerVectorIndex]->Fill(hit, storeSamplerMass, storeSamplerCharge,
734 }
735 }
736 // extra information - do only once at the end
737 if (storeSamplerIon)
738 {
739 for (auto& sampler : samplerCTrees)
740 {sampler->FillIon();}
741 }
742}
743
744void BDSOutput::FillSamplerSphereHitsVector(const std::vector<BDSHitsCollectionSamplerSphere*>& hits)
745{
746 for (const auto& hc : hits)
747 {
748 if (!hc)
749 {continue;} // could be nullptr
750 if (!(hc->entries() > 0))
751 {continue;}
752 for (int i = 0; i < (int) hc->entries(); i++)
753 {
754 const BDSHitSamplerSphere* hit = (*hc)[i];
755 G4int samplerID = hit->samplerID;
756 G4int samplerVectorIndex = samplerIDToIndexSphere[samplerID];
757 samplerSTrees[samplerVectorIndex]->Fill(hit, storeSamplerMass, storeSamplerCharge,
760 }
761 }
762 // extra information - do only once at the end
763 if (storeSamplerIon)
764 {
765 for (auto& sampler : samplerSTrees)
766 {sampler->FillIon();}
767 }
768}
769
770
772{
773 if (!(hits->entries() > 0))
774 {return;}
775 for (int i = 0; i < (int)hits->entries(); i++)
776 {
777 const BDSHitSampler* hit = (*hits)[i];
778 G4int samplerID = hit->samplerID;
779 G4int samplerVectorIndex = samplerIDToIndexPlane[samplerID];
780 samplerTrees[samplerVectorIndex]->Fill(hit, storeSamplerMass, storeSamplerCharge,
783 }
784
785 // extra information
786 if (storeSamplerIon)
787 {
788 for (auto& sampler : samplerTrees)
789 {sampler->FillIon();}
790 }
791}
792
794{
795 G4int nHits = (G4int)hits->entries();
796 if (nHits == 0) // integer so ok to compare
797 {return;}
798 for (G4int i = 0; i < nHits; i++)
799 {
800 const BDSHitSamplerLink* hit = (*hits)[i];
801 G4int samplerID = hit->samplerID;
802 samplerID += 1; // offset index by one due to primary branch.
804 }
805 // extra information - do only once at the end
806 G4bool firstSampler = true;
807 for (auto& sampler : samplerTrees)
808 {
809 if (firstSampler) // skip primaries (1st sampler) as it always has extras filled in
810 {firstSampler = false; continue;}
811 if (storeSamplerIon)
812 {sampler->FillIon();}
813 }
814}
815
817 const LossType lossType)
818{
819 switch (lossType)
820 {
821 case BDSOutput::LossType::world:
822 case BDSOutput::LossType::worldexit:
823 case BDSOutput::LossType::worldcontents:
824 {break;}
825 default:
826 {return; break;} // don't fill for other types of hits
827 }
828
829 G4int nHits = (G4int)hits->entries();
830 if (nHits == 0) // integer so ok to compare
831 {return;}
832 switch (lossType)
833 {
834 case BDSOutput::LossType::world:
835 {
836 for (G4int i=0; i < nHits; i++)
837 {
838 BDSHitEnergyDepositionGlobal* hit = (*hits)[i];
839 energyDepositedWorld += hit->TotalEnergyWeighted()/CLHEP::GeV;
840 if (storeELossWorld)
841 {eLossWorld->Fill(hit);}
842 }
843 break;
844 }
845 case BDSOutput::LossType::worldexit:
846 {
847 for (G4int i = 0; i < nHits; i++)
848 {
849 BDSHitEnergyDepositionGlobal* hit = (*hits)[i];
850 energyWorldExit += hit->TotalEnergyWeighted()/CLHEP::GeV;
851 energyWorldExitKinetic += hit->KineticEnergyWeighted()/CLHEP::GeV;
852 if (storeELossWorld)
853 {eLossWorldExit->Fill(hit);}
854 }
855 break;
856 }
857 case BDSOutput::LossType::worldcontents:
858 {
859 for (G4int i = 0; i < nHits; i++)
860 {
861 BDSHitEnergyDepositionGlobal* hit = (*hits)[i];
862 energyDepositedWorldContents += hit->TotalEnergyWeighted()/CLHEP::GeV;
864 {eLossWorldContents->Fill(hit);}
865 }
866 break;
867 }
868 default:
869 {break;} // only to prevent compiler warning
870 }
871}
872
874 const LossType lossType)
875{
876 G4int nHits = (G4int)hits->entries();
877 if (nHits == 0)
878 {return;}
879 G4int indELoss = histIndices1D["Eloss"];
880 G4int indELossPE = histIndices1D["ElossPE"];
881 G4int indELossTunnel = -1;
882 G4int indELossTunnelPE = -1;
884 {
885 indELossTunnel = histIndices1D["ElossTunnel"];
886 indELossTunnelPE = histIndices1D["ElossTunnelPE"];
887 }
888 G4int indELossVacuum = -1;
889 G4int indELossVacuumPE = -1;
891 {
892 indELossVacuum = histIndices1D["ElossVacuum"];
893 indELossVacuumPE = histIndices1D["ElossVacuumPE"];
894 }
895 G4int indScoringMap = -1;
896 if (useScoringMap)
897 {indScoringMap = histIndices3D["ScoringMap"];}
898 for (G4int i=0; i < nHits; i++)
899 {
900 BDSHitEnergyDeposition* hit = (*hits)[i];
901 G4double sHit = hit->GetSHit()/CLHEP::m;
902 G4double eW = hit->GetEnergyWeighted()/CLHEP::GeV;
903 switch (lossType)
904 {
905 case BDSOutput::LossType::energy:
906 {
907 energyDeposited += eW;
908 if (storeELoss)
909 {eLoss->Fill(hit);}
911 {
912 runHistos->Fill1DHistogram(indELoss, sHit, eW);
913 evtHistos->Fill1DHistogram(indELoss, sHit, eW);
914 runHistos->Fill1DHistogram(indELossPE, sHit, eW);
915 evtHistos->Fill1DHistogram(indELossPE, sHit, eW);
916 }
917 break;
918 }
919 case BDSOutput::LossType::vacuum:
920 {
923 {eLossVacuum->Fill(hit);}
925 {
926 evtHistos->Fill1DHistogram(indELossVacuum, sHit, eW);
927 runHistos->Fill1DHistogram(indELossVacuumPE, sHit, eW);
928 }
929 break;
930 }
931 case BDSOutput::LossType::tunnel:
932 {
935 {eLossTunnel->Fill(hit);}
937 {
938 runHistos->Fill1DHistogram(indELossTunnel, sHit, eW);
939 evtHistos->Fill1DHistogram(indELossTunnel, sHit, eW);
940 runHistos->Fill1DHistogram(indELossTunnelPE, sHit, eW);
941 evtHistos->Fill1DHistogram(indELossTunnelPE, sHit, eW);
942 }
943 break;
944 }
945 default:
946 {break;} // only to prevent compiler warning
947 }
948
949 if (useScoringMap)
950 {
951 G4double x = hit->Getx()/CLHEP::m;
952 G4double y = hit->Gety()/CLHEP::m;
953 evtHistos->Fill3DHistogram(indScoringMap, x, y, sHit, eW);
954 runHistos->Fill3DHistogram(indScoringMap, x, y, sHit, eW);
955 }
956 }
957
959 nCollimators > 0 &&
960 (lossType == BDSOutput::LossType::energy) &&
962 {CopyFromHistToHist1D("ElossPE", "CollElossPE", collimatorIndices);}
963}
964
965void BDSOutput::FillPrimaryHit(const std::vector<const BDSTrajectoryPointHit*>& primaryHits)
966{
967 for (auto phit : primaryHits)
968 {
969 if (!phit)
970 {continue;}
971 pFirstHit->Fill(phit);
972 const G4double preStepSPosition = phit->point->GetPreS() / CLHEP::m;
974 {
975 runHistos->Fill1DHistogram(histIndices1D["Phits"], preStepSPosition);
976 evtHistos->Fill1DHistogram(histIndices1D["Phits"], preStepSPosition);
977 runHistos->Fill1DHistogram(histIndices1D["PhitsPE"], preStepSPosition);
978 evtHistos->Fill1DHistogram(histIndices1D["PhitsPE"], preStepSPosition);
979
981 {CopyFromHistToHist1D("PhitsPE", "CollPhitsPE", collimatorIndices);}
982 }
983 }
984}
985
986void BDSOutput::FillPrimaryLoss(const std::vector<const BDSTrajectoryPointHit*>& primaryLosses)
987{
988 for (auto ploss : primaryLosses)
989 {
990 if (!ploss)
991 {continue;}
992 pLastHit->Fill(ploss);
993 const G4double postStepSPosition = ploss->point->GetPostS() / CLHEP::m;
995 {
996 runHistos->Fill1DHistogram(histIndices1D["Ploss"], postStepSPosition);
997 evtHistos->Fill1DHistogram(histIndices1D["Ploss"], postStepSPosition);
998 runHistos->Fill1DHistogram(histIndices1D["PlossPE"], postStepSPosition);
999 evtHistos->Fill1DHistogram(histIndices1D["PlossPE"], postStepSPosition);
1000
1002 {CopyFromHistToHist1D("PlossPE", "CollPlossPE", collimatorIndices);}
1003 }
1004 }
1005}
1006
1008{
1009 if (storeTrajectory)
1011}
1012
1014 const std::vector<const BDSTrajectoryPointHit*>& primaryLossPoints)
1015{
1016 G4int nHits = (G4int)hits->entries();
1017 for (G4int i = 0; i < nHits; i++)
1018 {
1019 BDSHitCollimator* hit = (*hits)[i];
1020 G4int collimatorIndex = hit->collimatorIndex;
1021 collimators[collimatorIndex]->Fill(hit,
1022 collimatorInfo[collimatorIndex],
1023 collimatorDifferences[collimatorIndex],
1024 storeCollimatorHits); // this includes the || storeCollimatorHitsLinks || storeCollimatorHitsIons);
1025 }
1026
1027 // identify whether the primary loss point was in a collimator
1028 // only do if there's a beam line, ie finished in a beam line, and there are collimators
1029 if (!primaryLossPoints.empty())
1030 {
1031 for (auto primaryLossPoint : primaryLossPoints)
1032 {
1033 if (primaryLossPoint->point->GetBeamLine() && nCollimators > 0)
1034 {
1035 G4int lossPointBLInd = primaryLossPoint->point->GetBeamLineIndex(); // always the mass world index
1036 auto result = std::find(collimatorIndices.begin(), collimatorIndices.end(), lossPointBLInd);
1037 if (result != collimatorIndices.end())
1038 {
1039 G4int collIndex = (int) (result - collimatorIndices.begin());
1040 collimators[collIndex]->SetPrimaryStopped(true);
1041 collimators[collIndex]->primaryInteracted = true;
1042 // it must've interacted if it stopped - could be that we kill
1043 // secondaries and there's no energy deposition therefore not identified
1044 // as primaryInteracted=true in BDSOutputROOTEventCollimator::Fill()
1045 }
1046 }
1047 }
1048 }
1049
1050 // if required loop over collimators and get them to calculate and fill extra information
1052 {
1053 for (auto collimator : collimators)
1054 {collimator->FillExtras(storeCollimatorHitsIons, storeCollimatorHitsLinks);}
1055 }
1056
1057 // after all collimator hits have been filled, we summarise whether the primary
1058 // interacted in a histogram
1059 G4int histIndex = histIndices1D["CollPInteractedPE"];
1060 for (G4int i = 0; i < (G4int)collimators.size(); i++)
1061 {evtHistos->Fill1DHistogram(histIndex, i, (int)collimators[i]->primaryInteracted);}
1062
1063
1064 // loop over collimators and count the number that were interacted with in this event
1065 for (const auto collimator : collimators)
1066 {
1067 if (collimator->primaryInteracted)
1069 }
1070}
1071
1073{
1075 {return;}
1076
1077 G4int nPrimaryImpacts = 0;
1078 G4int nHits = (G4int)hits->entries();
1079 G4int histIndex = histIndices1D["PFirstAI"];
1080 for (G4int i = 0; i < nHits; i++)
1081 {
1082 const BDSHitApertureImpact* hit = (*hits)[i];
1083 G4double eW = (hit->totalEnergy / CLHEP::GeV) * hit->weight;
1085 energyImpactingApertureKinetic += (hit->preStepKineticEnergy * hit->weight ) / CLHEP::GeV;
1086 if (hit->parentID == 0)
1087 {
1088 nPrimaryImpacts += 1;
1089 // only store one primary aperture hit in this histogram even if they were multiple
1090 if (storeApertureImpactsHistograms && nPrimaryImpacts == 1)
1091 {evtHistos->Fill1DHistogram(histIndex, hit->S / CLHEP::m);}
1092 }
1093 // hits are generated in order as the particle progresses
1094 // through the model, so the first one in the collection
1095 // for the primary is the first one in S.
1097 {apertureImpacts->Fill(hit, nPrimaryImpacts==1);}
1098 }
1099}
1100
1101void BDSOutput::FillScorerHits(const std::map<G4String, G4THitsMap<G4double>*>& scorerHitsMap)
1102{
1103 for (const auto& nameHitsMap : scorerHitsMap)
1104 {
1105#if G4VERSION < 1039
1106 if (nameHitsMap.second->GetSize() == 0)
1107#else
1108 if (nameHitsMap.second->size() == 0)
1109#endif
1110#ifdef BDSDEBUG
1111 {G4cout << nameHitsMap.first << " empty" << G4endl; continue;}
1112#else
1113 {continue;}
1114#endif
1115 FillScorerHitsIndividual(nameHitsMap.first, nameHitsMap.second);
1116 }
1117}
1118
1119void BDSOutput::FillScorerHitsIndividual(const G4String& histogramDefName,
1120 const G4THitsMap<G4double>* hitMap)
1121{
1122 if (BDS::StrContains(histogramDefName, "blm_"))
1123 {return FillScorerHitsIndividualBLM(histogramDefName, hitMap);}
1124
1125 if (!(histIndices3D.find(histogramDefName) == histIndices3D.end()))
1126 {
1127 G4int histIndex = histIndices3D[histogramDefName];
1128 G4double unit = BDS::MapGetWithDefault(histIndexToUnits3D, histIndex, 1.0);
1129 // avoid using [] operator for map as we have no default constructor for BDSHistBinMapper3D
1130 const BDSHistBinMapper& mapper = scorerCoordinateMaps.at(histogramDefName);
1131 TH3D* hist = evtHistos->Get3DHistogram(histIndex);
1132 G4int x,y,z,e;
1133#if G4VERSION < 1039
1134 for (const auto& hit : *hitMap->GetMap())
1135#else
1136 for (const auto& hit : *hitMap)
1137#endif
1138 {
1139 // convert from scorer global index to 3d i,j,k index of 3d scorer
1140 mapper.IJKLFromGlobal(hit.first, x,y,z,e);
1141 G4int rootGlobalIndex = (hist->GetBin(x + 1, y + 1, z + 1)); // convert to root system (add 1 to avoid underflow bin)
1142 evtHistos->Set3DHistogramBinContent(histIndex, rootGlobalIndex, *hit.second / unit);
1143 }
1145 }
1146
1147 if (!(histIndices4D.find(histogramDefName) == histIndices4D.end()))
1148 {
1149 G4int histIndex = histIndices4D[histogramDefName];
1150 G4double unit = BDS::MapGetWithDefault(histIndexToUnits4D, histIndex, 1.0);
1151 // avoid using [] operator for map as we have no default constructor for BDSHistBinMapper3D
1152 const BDSHistBinMapper& mapper = scorerCoordinateMaps.at(histogramDefName);
1153 G4int x,y,z,e;
1154#if G4VERSION < 1039
1155 for (const auto& hit : *hitMap->GetMap())
1156#else
1157 for (const auto& hit : *hitMap)
1158#endif
1159 {
1160 // convert from scorer global index to 4d i,j,k,e index of 4d scorer
1161 mapper.IJKLFromGlobal(hit.first, x,y,z,e);
1162 evtHistos->Set4DHistogramBinContent(histIndex, x, y, z, e - 1, *hit.second / unit); // - 1 to go back to the Boost Histogram indexing (-1 for the underflow bin)
1163 }
1164 runHistos->AccumulateHistogram4D(histIndex, evtHistos->Get4DHistogram(histIndex));
1165 }
1166}
1167
1168void BDSOutput::FillScorerHitsIndividualBLM(const G4String& histogramDefName,
1169 const G4THitsMap<G4double>* hitMap)
1170{
1171 G4int histIndex = blmCollectionNameToHistogramID[histogramDefName];
1172#if G4VERSION < 1039
1173 for (const auto& hit : *hitMap->GetMap())
1174#else
1175 for (const auto& hit : *hitMap)
1176#endif
1177 {
1178#ifdef BDSDEBUG
1179 G4cout << "Filling hist " << histIndex << ", bin: " << hit.first+1 << " value: " << *hit.second << G4endl;
1180#endif
1181 G4double unit = BDS::MapGetWithDefault(histIndexToUnits1D, histIndex, 1.0);
1182 evtHistos->Fill1DHistogram(histIndex,hit.first, *hit.second / unit);
1183 runHistos->Fill1DHistogram(histIndex,hit.first, *hit.second / unit);
1184 }
1185}
1186
1188{
1189 if (info)
1191}
1192
1193void BDSOutput::CopyFromHistToHist1D(const G4String& sourceName,
1194 const G4String& destinationName,
1195 const std::vector<G4int>& indices)
1196{
1197 TH1D* sourceEvt = evtHistos->Get1DHistogram(histIndices1D[sourceName]);
1198 TH1D* destinationEvt = evtHistos->Get1DHistogram(histIndices1D[destinationName]);
1199 // for the run ones we are overwriting but this is ok
1200 TH1D* sourceRun = runHistos->Get1DHistogram(histIndices1D[sourceName]);
1201 TH1D* destinationRun = runHistos->Get1DHistogram(histIndices1D[destinationName]);
1202 G4int binIndex = 1; // starts at 1 for TH1; 0 is underflow
1203 for (const auto index : indices)
1204 {
1205 destinationEvt->SetBinContent(binIndex, sourceEvt->GetBinContent(index + 1));
1206 destinationEvt->SetBinError(binIndex, sourceEvt->GetBinError(index + 1));
1207 destinationRun->SetBinContent(binIndex, sourceRun->GetBinContent(index + 1));
1208 destinationRun->SetBinError(binIndex, sourceRun->GetBinError(index + 1));
1209 binIndex++;
1210 }
1211}
const std::map< G4String, BDSScorerHistogramDef > & ScorerHistogramDefinitionsMap() const
Access all scorer histogram definitions.
const std::map< G4String, G4Transform3D > & ScorerMeshPlacementsMap() const
Access all scorer histogram definitions.
BDSBeamline * TunnelBeamline() const
Access the beam line containing all the tunnel segments.
const BDSBeamline * BeamlineMain() const
Accessor.
static BDSBLMRegistry * Instance()
Accessor for registry.
G4double GetSPositionEnd() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
std::vector< G4double > GetEdgeSPositions() const
Definition: BDSBeamline.cc:906
const BDSBeamlineElement * GetLastItem() const
Return a reference to the last element.
Definition: BDSBeamline.hh:126
G4bool empty() const
Iterator mechanics.
Definition: BDSBeamline.hh:175
Interface to store event information use G4 hooks.
Definition: BDSEventInfo.hh:42
const BDSOutputROOTEventInfo * GetInfo() const
Accessor.
Definition: BDSEventInfo.hh:64
General exception with possible name of object and message.
Definition: BDSException.hh:35
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
BDS::TrajectoryOptions StoreTrajectoryOptions() const
options that require some implementation.
G4bool UseImportanceSampling() const
Is importance sampling being used.
Mapping from axis indices to 1D index.
Snapshot of information for particle passing through a collimator.
Snapshot of information for particle passing through a collimator.
G4int collimatorIndex
Index of collimator the hit is in.
Information recorded for a step leaving a volume.
Information recorded for a single piece of energy deposition.
G4double Getx() const
Accessor for extra piece of information.
G4double Gety() const
Accessor for extra piece of information.
The information recorded from a particle impacting a sampler.
The information recorded from a particle impacting a sampler.
The information recorded from a particle impacting a sampler.
G4int A() const
Accessor.
G4int Z() const
Accessor.
Class to store all beam options for a BDSIM run.
void Fill(const BDSOutputROOTEventCoords *other)
Fill from another instance.
void Fill(const std::vector< std::string > &analysedFilesIn=std::vector< std::string >(), const std::vector< std::string > &combinedFilesIn=std::vector< std::string >())
void AccumulateHistogram3D(G4int histoId, TH3D *otherHistogram)
Add the values from one supplied 3D histogram to another. Uses TH3-Add().
BDSBH4DBase * Get4DHistogram(int iHisto) const
Accessors.
TH1D * Get1DHistogram(int iHisto) const
Accessors.
void Set3DHistogramBinContent(G4int histoId, G4int globalBinID, G4double value)
TH3D * Get3DHistogram(int iHisto) const
Accessors.
double energyKilled
Total energy of killed particles that weren't in a sensitive volume.
double energyWorldExitKinetic
Total kinetic energy leaving the world.
double energyWorldExit
Total energy leaving the world.
double energyDepositedWorldContents
Total energy deposited in the world contents for this event.
double energyImpactingApertureKinetic
Total kinetic energy impacting the aperture.
double energyDepositedWorld
Total energy deposited in the world for this event.
double energyDepositedTunnel
Total energy deposited in the tunnel for this event.
double energyDepositedVacuum
Total energy deposited in vacuum volumes.
int nCollimatorsInteracted
Number of collimators primary interacted with.
double energyTotal
Sum of above 5 variables that totals all energy.
double energyDeposited
Total energy deposited in machine (not world or tunnel).
double energyImpactingAperture
Total energy impacting the aperture.
void Fill(const BDSOutputROOTEventLoss *other)
Fill from another instance.
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.
Class to store all options for a BDSIM run.
Information pertaining to a run.
Information stored per sampler per event.
void Fill(G4bool fillIons)
Fill maps of particle information from Geant4.
virtual void Flush()
Clear maps.
Holder for output information.
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.
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.
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 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
BDSOutputROOTEventAperture * apertureImpacts
Impacts on the aperture.
BDSOutputROOTParticleData * particleDataOutput
Geant4 information / particle tables.
BDSOutputROOTEventOptions * optionsOutput
Options output.
virtual void WriteModel()=0
Write a representation of the samplers and beamline.
void FillEnergyLoss(const BDSHitsCollectionEnergyDeposition *loss, const LossType type)
Fill a collection of energy hits into the appropriate output structure.
Definition: BDSOutput.cc:873
G4double energyDepositedVacuum
Integral when filling hit.
Definition: BDSOutput.hh:327
virtual void WriteParticleData()=0
Write the geant4 information.
G4bool createCollimatorOutputStructures
Definition: BDSOutput.hh:323
virtual void WriteBeam()=0
Write the beam.
virtual void NewFile()=0
Open a new file. This should call WriteHeader() in it.
G4double sMaxHistograms
Definition: BDSOutput.hh:294
G4bool storeCollimatorHitsIons
Storage option.
Definition: BDSOutput.hh:303
G4bool storeSamplerCharge
Storage option.
Definition: BDSOutput.hh:311
G4bool storeSamplerRigidity
Storage option.
Definition: BDSOutput.hh:314
virtual void WriteFileEventLevel()=0
BDSOutput()=delete
No default constructor.
G4bool storeCollimatorHitsLinks
Storage option.
Definition: BDSOutput.hh:302
G4bool storeELossHistograms
Storage option.
Definition: BDSOutput.hh:304
G4bool storePrimaryHistograms
Storage option.
Definition: BDSOutput.hh:308
G4double energyDepositedTunnel
Integral when filling hit.
Definition: BDSOutput.hh:330
G4bool storeSamplerKineticEnergy
Storage option.
Definition: BDSOutput.hh:312
void CloseAndOpenNewFile()
Close a file and open a new one.
Definition: BDSOutput.cc:364
void FillEventInfo(const BDSEventInfo *info)
Fill event summary information.
Definition: BDSOutput.cc:666
virtual void WriteOptions()=0
Write the options.
std::map< G4String, G4int > blmCollectionNameToHistogramID
Definition: BDSOutput.hh:178
std::map< G4int, G4double > histIndexToUnits1D
Definition: BDSOutput.hh:347
G4bool storePrimaries
Options for dynamic bits of output.
Definition: BDSOutput.hh:170
G4bool storeSamplerPolarCoords
Storage option.
Definition: BDSOutput.hh:310
void FillScorerHitsIndividual(const G4String &hsitogramDefName, const G4THitsMap< G4double > *hitMap)
Fill an individual scorer hits map into a particular output histogram.
Definition: BDSOutput.cc:1119
const G4String baseFileName
Base file name.
Definition: BDSOutput.hh:281
void FillSamplerHitsVector(const std::vector< BDSHitsCollectionSampler * > &hits)
Fill sampler hits from a vector<sampler hits collection>.
Definition: BDSOutput.cc:692
G4bool storeTrajectory
Options for dynamic bits of output.
Definition: BDSOutput.hh:171
void FillSamplerHitsLink(const BDSHitsCollectionSamplerLink *hits)
Fill sampler link hits into output structures.
Definition: BDSOutput.cc:793
G4int nbins
Number of bins for each histogram required.
Definition: BDSOutput.hh:297
G4bool useScoringMap
Whether the single 3D histogram will be built.
Definition: BDSOutput.hh:290
void CopyFromHistToHist1D(const G4String &sourceName, const G4String &destinationName, const std::vector< G4int > &indices)
Definition: BDSOutput.cc:1193
G4String GetNextFileName()
Get the next file name based on the base file name and the accrued number of files.
Definition: BDSOutput.cc:390
virtual void InitialiseGeometryDependent()
Definition: BDSOutput.cc:166
void FillPrimary(const G4PrimaryVertex *vertex, const G4int turnsTaken)
Definition: BDSOutput.cc:235
G4double energyDepositedWorldContents
Integral when filling hit.
Definition: BDSOutput.hh:329
const G4String fileExtension
File extension to add to each file.
Definition: BDSOutput.hh:282
G4bool storeSamplerIon
Storage option.
Definition: BDSOutput.hh:315
G4bool storeModel
Storage option.
Definition: BDSOutput.hh:309
void FillRunInfo(const BDSEventInfo *info)
Fill run level summary information.
Definition: BDSOutput.cc:1187
G4double energyWorldExitKinetic
Integral when filling hit.
Definition: BDSOutput.hh:334
std::map< G4String, G4int > histIndices1D
Map of histogram name (short) to index of histogram in output.
Definition: BDSOutput.hh:339
void FillTrajectories(const BDSTrajectoriesToStore *trajectories)
Copy a set of trajectories to the output structure.
Definition: BDSOutput.cc:1007
G4bool storeELoss
Options for dynamic bits of output.
Definition: BDSOutput.hh:163
void FillParticleData(G4bool writeIons)
Fill the local structure particle data with information. Also calls WriteParticleData().
Definition: BDSOutput.cc:186
G4bool storeParticleData
Storage option.
Definition: BDSOutput.hh:307
G4bool storeApertureImpactsHistograms
Options for dynamic bits of output.
Definition: BDSOutput.hh:169
G4double energyImpactingApertureKinetic
Integral when filling hit.
Definition: BDSOutput.hh:332
G4double energyImpactingAperture
Integral when filling hit.
Definition: BDSOutput.hh:331
void FillEventPrimaryOnly(const BDSParticleCoordsFullGlobal &coords, const BDSParticleDefinition *particle)
Definition: BDSOutput.cc:263
G4bool storeELossWorldContents
Options for dynamic bits of output.
Definition: BDSOutput.hh:167
G4bool storeCollimatorInfo
Storage option.
Definition: BDSOutput.hh:300
virtual void WriteFileRunLevel()=0
void FillRun(const BDSEventInfo *info)
Copy run information to output structure.
Definition: BDSOutput.cc:371
void FillBeam(const GMAD::BeamBase *beam)
Definition: BDSOutput.cc:204
void FillHeader()
Fill the local structure header with information - updates time stamp.
Definition: BDSOutput.cc:178
void FillPrimaryLoss(const std::vector< const BDSTrajectoryPointHit * > &primaryLosses)
Fill a collection volume exit hits into the appropriate output structure.
Definition: BDSOutput.cc:986
static void PrintProtectedNames(std::ostream &out)
Feedback for protected names.
Definition: BDSOutput.cc:383
std::map< G4String, G4int > histIndices4D
Map of histogram name (short) to index of histogram in output.
Definition: BDSOutput.hh:341
void FillPrimaryHit(const std::vector< const BDSTrajectoryPointHit * > &primaryHits)
Fill the hit where the primary particle impact.
Definition: BDSOutput.cc:965
G4int outputFileNumber
Number of output file.
Definition: BDSOutput.hh:284
void FillApertureImpacts(const BDSHitsCollectionApertureImpacts *hits)
Fill aperture impact hits.
Definition: BDSOutput.cc:1072
void FillSamplerHits(const BDSHitsCollectionSampler *hits)
Fill sampler hits into output structures.
Definition: BDSOutput.cc:771
G4bool storeELossWorld
Options for dynamic bits of output.
Definition: BDSOutput.hh:166
void FillEvent(const BDSEventInfo *info, const G4PrimaryVertex *vertex, const std::vector< BDSHitsCollectionSampler * > &samplerHitsPlane, const std::vector< BDSHitsCollectionSamplerCylinder * > &samplerHitsCylinder, const std::vector< BDSHitsCollectionSamplerSphere * > &samplerHitsSphere, const BDSHitsCollectionSamplerLink *samplerHitsLink, const BDSHitsCollectionEnergyDeposition *energyLoss, const BDSHitsCollectionEnergyDeposition *energyLossFull, const BDSHitsCollectionEnergyDeposition *energyLossVacuum, const BDSHitsCollectionEnergyDeposition *energyLossTunnel, const BDSHitsCollectionEnergyDepositionGlobal *energyLossWorld, const BDSHitsCollectionEnergyDepositionGlobal *energyLossWorldContents, const BDSHitsCollectionEnergyDepositionGlobal *worldExitHits, const std::vector< const BDSTrajectoryPointHit * > &primaryHits, const std::vector< const BDSTrajectoryPointHit * > &primaryLosses, const BDSTrajectoriesToStore *trajectories, const BDSHitsCollectionCollimator *collimatorHits, const BDSHitsCollectionApertureImpacts *apertureImpactHits, const std::map< G4String, G4THitsMap< G4double > * > &scorerHitsMap, const G4int turnsTaken)
Copy event information from Geant4 simulation structures to output structures.
Definition: BDSOutput.cc:289
void CreateHistograms()
Create histograms.
Definition: BDSOutput.cc:447
G4int storeTrajectoryStepPoints
Storage option.
Definition: BDSOutput.hh:316
G4bool storeTrajectoryStepPointLast
Storage option.
Definition: BDSOutput.hh:317
G4bool storeELossVacuum
Options for dynamic bits of output.
Definition: BDSOutput.hh:165
std::map< G4String, BDSHistBinMapper > scorerCoordinateMaps
Map of histogram name (short) to index of histogram in output.
Definition: BDSOutput.hh:342
G4bool storeSamplerMass
Storage option.
Definition: BDSOutput.hh:313
void FillModel()
Definition: BDSOutput.cc:218
G4double energyWorldExit
Integral when filling hit.
Definition: BDSOutput.hh:333
G4int nCollimatorsInteracted
Integral when filling hit.
Definition: BDSOutput.hh:335
void FillOptions(const GMAD::OptionsBase *options)
Definition: BDSOutput.cc:211
G4bool storeELossVacuumHistograms
Storage option.
Definition: BDSOutput.hh:306
LossType
Enum for different types of energy loss that can be written out.
Definition: BDSOutput.hh:182
void CalculateHistogramParameters()
Calculate the number of bins and required maximum s.
Definition: BDSOutput.cc:424
G4bool storeELossTunnel
Options for dynamic bits of output.
Definition: BDSOutput.hh:164
G4double energyDeposited
Integral when filling hit.
Definition: BDSOutput.hh:326
static const std::set< G4String > protectedNames
Invalid names for samplers - kept here as this is where the output structures are created.
Definition: BDSOutput.hh:287
G4int numberEventPerFile
Number of events stored per file.
Definition: BDSOutput.hh:283
G4double energyDepositedWorld
Integral when filling hit.
Definition: BDSOutput.hh:328
G4bool storeCollimatorHits
Storage option.
Definition: BDSOutput.hh:301
void FillScorerHits(const std::map< G4String, G4THitsMap< G4double > * > &scorerHitsMap)
Fill a map of scorer hits into the output.
Definition: BDSOutput.cc:1101
static G4bool InvalidSamplerName(const G4String &samplerName)
Test whether a sampler name is invalid or not.
Definition: BDSOutput.cc:378
void FillCollimatorHits(const BDSHitsCollectionCollimator *hits, const std::vector< const BDSTrajectoryPointHit * > &primaryLossPoints)
Fill collimator hits.
Definition: BDSOutput.cc:1013
BDS::TrajectoryOptions storeTrajectoryOptions
Storage option.
Definition: BDSOutput.hh:318
std::map< G4String, G4int > histIndices3D
Map of histogram name (short) to index of histogram in output.
Definition: BDSOutput.hh:340
G4bool storeApertureImpacts
Options for dynamic bits of output.
Definition: BDSOutput.hh:168
virtual void CloseFile()=0
G4bool storeELossTunnelHistograms
Storage option.
Definition: BDSOutput.hh:305
virtual void WriteHeader()=0
Write the header.
A set of particle coordinates in both local and global.
Wrapper for particle definition.
G4double Mass() const
Accessor.
G4bool NElectrons() const
Accessor.
G4double Momentum() const
Accessor.
G4double Charge() const
Accessor.
G4double BRho() const
Accessor.
BDSIonDefinition * IonDefinition() const
Accessor.
G4bool IsAnIon() const
Accessor.
Full set of coordinates for association with primary vertex. Vector version.
Full set of coordinates for association with primary vertex.
const std::map< G4String, G4double > & PrimitiveScorerUnits() const
Access the map of units for primitive scorers.
const std::vector< G4String > & PrimitiveScorerNames() const
Access a vector of the just primitive scorer part of the names.
const std::vector< G4String > & PrimitiveScorerNamesComplete() const
Access a vector the full primitive scorer names as registered.
Double map of trajectories to bitset of which filters matched whether to store them.
Options for a beam distribution.
Definition: beamBase.h:35
Basic options class independent of Geant4.
Definition: optionsBase.h:36
G4bool StrContains(const G4String &str, const G4String &test)
Utility function to simplify lots of syntax changes for pedantic g4 changes.
Definition: BDSUtilities.cc:66
V MapGetWithDefault(const std::map< K, V > &m, const K &key, const V &defaultValue)
G4bool FileExists(const G4String &filename)
Checks if filename exists.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())