BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
HeaderAnalysis.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 "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 <stdexcept>
28#include <string>
29#include <vector>
30
31HeaderAnalysis::HeaderAnalysis(const std::vector<std::string>& filenamesIn,
32 Header* headerIn,
33 TChain* chainIn):
34 filenames(filenamesIn),
35 header(headerIn),
36 chain(chainIn)
37{;}
38
39HeaderAnalysis::~HeaderAnalysis() noexcept
40{;}
41
42unsigned long long int HeaderAnalysis::CountNOriginalEvents()
43{
44 unsigned long long int nOriginalEvents = 0;
45 for (int i = 0; i < chain->GetEntries(); i++) // assumes 1 header entry per file - fine
46 {
47 chain->GetEntry(i);
48 if (header->header->nOriginalEvents > 0)
49 {nOriginalEvents += header->header->nOriginalEvents;}
50 else
51 {// not a skimmed file, so nOriginalEvents in header is 0 -> get info from that individual file
52 TFile* ftemp;
53 try
54 {ftemp = new TFile(filenames[i].c_str());}
55 catch (const std::exception& e)
56 {continue;}
57 TTree* eventTree = (TTree*)ftemp->Get("Event");
58 if (!eventTree)
59 {
60 ftemp->Close();
61 delete ftemp;
62 continue;
63 }
64 Long64_t nentries = eventTree->GetEntries();
65 if (nentries > 0)
66 {nOriginalEvents += (unsigned long long int)nentries;}
67 ftemp->Close();
68 delete ftemp;
69 }
70 }
71 return nOriginalEvents;
72}
unsigned long long int nOriginalEvents
Number of original events if skimmed.
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