BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
Compare.hh
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#ifndef COMPCOMPARE_H
20#define COMPCOMPARE_H
21
22#include <vector>
23
24#include "ResultEvent.hh"
25
26#include "BDSOutputROOTEventSampler.hh"
27
28class Result;
29
30class TDirectory;
31class TFile;
32class TH1;
33class TTree;
34
40namespace Compare
41{
42 const static double Chi2Tolerance = 40;
43 const static double TreeTolerance = 0.05;
44 const static double OpticsSimgaTolerance = 10;
45 const static double EventTreeTolerance = 1e-10;
46
48 std::vector<Result*> Files(TFile* f1, TFile* f2);
49
53 void Directories(TDirectory* d1,
54 TDirectory* d2,
55 std::vector<Result*>& results);
56
58 void Histograms(TH1* h1, TH1* h2, std::vector<Result*>& results);
59
61 void Trees(TTree* t1, TTree* t2, std::vector<Result*>& results);
62
65 void Optics(TTree* t1, TTree* t2, std::vector<Result*>& results);
66
67 void EventTree(TTree* t1, TTree* t2, std::vector<Result*>& results,
68 const std::vector<std::string>& samplerNames,
69 const std::vector<std::string>& samplerCNames,
70 const std::vector<std::string>& samplerSNames);
71
72#ifdef __ROOTDOUBLE__
73 void Sampler(BDSOutputROOTEventSampler<double>* e1,
75 ResultEvent* results);
76#else
77 void Sampler(BDSOutputROOTEventSampler<float>* e1,
79 ResultEvent* results);
80#endif
81
83 template <typename T>
84 bool Diff(const std::vector<T>& v1, std::vector<T>& v2, int i) {return std::abs(v1[i] - v2[i]) > (T)EventTreeTolerance;}
85 template <>
86 bool Diff(const std::vector<bool>& v1, std::vector<bool>& v2, int i) {return v1[i] != v2[i];}
87
88 template <typename T>
89 bool Diff(T v1, T v2) {return std::abs(v1 - v2) > (T)Compare::EventTreeTolerance;}
90 template <>
91 bool Diff(bool v1, bool v2) {return v1 != v2;}
92
95 void PrintNoMatching(std::string className, std::string objectName);
96
98 bool Summarise(const std::vector<Result*>& results);
99
101 bool StringStartsWith(std::string aString, std::string prefix);
102
103 bool IsInVector(std::string key, const std::vector<std::string>& vec);
104
105 inline bool NanOrInf(const double& val) {return std::isnan(val) || std::isinf(val);}
106
107 inline bool GTEZero(const double& val) {return val >= 0;}
108 inline bool LTZero(const double& val) {return val < 0;}
109
112}
113
114#endif
Information stored per sampler per event.
Result of comparison of single entry of an Event tree in BDSIM output.
Definition: ResultEvent.hh:35
Base class of a result.
Definition: Result.hh:37
Comparison utility functions.
Definition: Compare.hh:41
std::vector< Result * > Files(TFile *f1, TFile *f2)
Compare two files.
Definition: Compare.cc:53
void Trees(TTree *t1, TTree *t2, std::vector< Result * > &results)
Compare two TTrees.
Definition: Compare.cc:183
bool Summarise(const std::vector< Result * > &results)
Loop over results and print any failed ones. Returns true if all passed.
bool Diff(const std::vector< T > &v1, std::vector< T > &v2, int i)
Return true if they're different.
Definition: Compare.hh:84
bool hasPrimaries
Bool for checking if primaries were written to file.
Definition: Compare.hh:111
void PrintNoMatching(std::string className, std::string objectName)
Definition: Compare.cc:576
void Histograms(TH1 *h1, TH1 *h2, std::vector< Result * > &results)
Compare two histogams.
Definition: Compare.cc:134
void Directories(TDirectory *d1, TDirectory *d2, std::vector< Result * > &results)
Definition: Compare.cc:61
void Optics(TTree *t1, TTree *t2, std::vector< Result * > &results)
Definition: Compare.cc:265
bool StringStartsWith(std::string aString, std::string prefix)
Check whether a string is prefixed with another string.
Definition: Compare.cc:581