BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
HeaderAnalysis.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 "BDSOutputROOTEventHeader.hh"
20#include "Header.hh"
21#include "HeaderAnalysis.hh"
22
23#include "TChain.h"
24#include "TFile.h"
25#include "TTree.h"
26
27#include <set>
28#include <stdexcept>
29#include <string>
30#include <vector>
31
32HeaderAnalysis::HeaderAnalysis(const std::vector<std::string>& filenamesIn,
33 Header* headerIn,
34 TChain* chainIn):
35 filenames(filenamesIn),
36 header(headerIn),
37 chain(chainIn)
38{;}
39
40HeaderAnalysis::~HeaderAnalysis() noexcept
41{;}
42
43unsigned long long int HeaderAnalysis::CountNOriginalEvents(unsigned long long int& nEventsInFileIn,
44 unsigned long long int& nEventsInFileSkippedIn,
45 unsigned long long int& nEventsRequestedIn,
46 unsigned int& distrFileLoopNTimesIn)
47{
48 unsigned long long int nOriginalEvents = 0;
49 nEventsInFileIn = 0;
50 nEventsInFileSkippedIn = 0;
51 nEventsRequestedIn = 0;
52 distrFileLoopNTimesIn = 0;
53
54 for (const auto& fn : filenames)
55 {
56 TFile* ft = nullptr;
57 try
58 {ft = new TFile(fn.c_str(), "READ");}
59 catch (const std::exception& e)
60 {continue;}
61 if (ft->IsZombie())
62 {
63 delete ft;
64 continue;
65 }
66 Header* ha = new Header();
67 TTree* ht = dynamic_cast<TTree*>(ft->Get("Header"));
68 if (!ht)
69 {
70 delete ft;
71 continue;
72 }
73 ha->SetBranchAddress(ht);
74 ht->GetEntry(ht->GetEntries()-1); // get the last entry
75
76 nOriginalEvents += ha->header->nOriginalEvents;
77 nEventsInFileIn += ha->header->nEventsInFile;
78 nEventsInFileSkippedIn += ha->header->nEventsInFileSkipped;
79 nEventsRequestedIn += ha->header->nEventsRequested;
80 distrFileLoopNTimesIn += ha->header->distrFileLoopNTimes;
81 ft->Close();
82 delete ft;
83 delete ha;
84 }
85 return nOriginalEvents;
86}
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.
unsigned long long int nEventsInFile
Number of events in the input distribution file irrespective of filters.
unsigned long long int CountNOriginalEvents(unsigned long long int &nEventsInFileIn, unsigned long long int &nEventsInFileSkippedIn, unsigned long long int &nEventsRequestedIn, unsigned int &distrFileLoopNTimesIn)
HeaderAnalysis(const std::vector< std::string > &filenamesIn, Header *headerIn, TChain *chainIn)
Constructor intended for use to construct a model analysis object.
Header loader.
Definition: Header.hh:34
BDSOutputROOTEventHeader * header
Member that ROOT can map file data to locally.
Definition: Header.hh:44
void SetBranchAddress(TTree *t)
Set the branch addresses to address the contents of the file.
Definition: Header.cc:41