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