BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSOutputLoaderSampler.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 "BDSDebug.hh"
20#include "BDSException.hh"
21#include "BDSOutputLoaderSampler.hh"
22#include "BDSOutputROOTEventOptions.hh"
23#include "BDSOutputROOTEventSampler.hh"
24#include "BDSUtilities.hh"
25
26#include "globals.hh"
27
28#include "RtypesCore.h"
29#include "TFile.h"
30#include "TTree.h"
31
32#include <string>
33
34BDSOutputLoaderSampler::BDSOutputLoaderSampler(const G4String& filePath,
35 const G4String& samplerName):
36 BDSOutputLoader(filePath),
37 doublePrecision(false),
38 nEvents(0),
39 localSamplerFloat(nullptr),
40 localSamplerDouble(nullptr)
41{
42 G4String samplerNameLocal = samplerName;
43 if (!BDS::EndsWith(samplerNameLocal, "."))
44 {samplerNameLocal += ".";}
45
46 TBranch* search = eventTree->GetBranch(samplerNameLocal);
47 if (!search)
48 {throw BDSException(__METHOD_NAME__, "no such sampler name \"" + samplerName + "\"");}
49
50 localSamplerDouble = new BDSOutputROOTEventSampler<double>();
51 localSamplerFloat = new BDSOutputROOTEventSampler<float>();
52 doublePrecision = localOptions->outputDoublePrecision;
53 if (localOptions->outputDoublePrecision)
54 {eventTree->SetBranchAddress(samplerNameLocal, &localSamplerDouble);}
55 else
56 {eventTree->SetBranchAddress(samplerNameLocal, &localSamplerFloat);}
57 nEvents = (G4long)eventTree->GetEntries();
58}
59
60BDSOutputLoaderSampler::~BDSOutputLoaderSampler()
61{
62 delete localSamplerFloat;
63 delete localSamplerDouble;
64}
65
66const BDSOutputROOTEventSampler<float>* BDSOutputLoaderSampler::SamplerDataFloat(G4long eventNumber)
67{
68 localSamplerFloat->Flush();
69 Common(eventNumber);
70 return localSamplerFloat;
71}
72
73const BDSOutputROOTEventSampler<double>* BDSOutputLoaderSampler::SamplerDataDouble(G4long eventNumber)
74{
75 localSamplerDouble->Flush();
76 Common(eventNumber);
77 return localSamplerDouble;
78}
79
80void BDSOutputLoaderSampler::Common(G4long eventNumber)
81{
82 if (eventNumber > nEvents)
83 {
84 G4String msg = "requested event number " + std::to_string(eventNumber);
85 msg += " is beyond number of entries in the file (" + std::to_string(nEvents) + ")";
86 throw BDSException(__METHOD_NAME__, msg);
87 }
88 // always change back to this file - assuming other root files could be open
89 file->cd();
90
91 // cannot retrieve a seed state beyond that in the file - protection here to
92 // make life simpler elsewhere
93 if (eventNumber > eventTree->GetEntries())
94 {
95 G4cout << __METHOD_NAME__ << "event index beyond number stored in file - no seed state loaded" << G4endl;
96 return;
97 }
98 eventTree->GetEntry((Long64_t)eventNumber);
99}
General exception with possible name of object and message.
Definition: BDSException.hh:35
Loader of ROOT Event output for recreating events.
Information stored per sampler per event.
virtual void Flush()
Clean Sampler.
G4bool EndsWith(const std::string &expression, const std::string &suffix)
Return true if a string "expression" ends with "suffix".