Tools for Reactor-Antineutrino Studies in RAT-PAC  1.2.0
Code for running and analyzing reactor-IBD simulations
TRefMatch.h
1 // TRefMatch
3 
16 // ~ Mark J. Duvall ~ mjduvall@hawaii.edu ~ 10/2022 ~ //
17 
18 //Copyright (C) 2022 Mark J. Duvall / T. Rocks Science
19 //
20 // This program is free software: you can redistribute it and/or modify
21 // it under the terms of the GNU General Public License as published by
22 // the Free Software Foundation, either version 3 of the License, or
23 // (at your option) any later version.
24 //
25 // This program is distributed in the hope that it will be useful,
26 // but WITHOUT ANY WARRANTY; without even the implied warranty of
27 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
28 // GNU General Public License for more details.
29 //
30 // You should have received a copy of the GNU General Public License
31 // along with this program. If not, see <http://www.gnu.org/licenses/>.
32 
33 #ifndef TRefMatch
34 #define TRefMatch
35 
36 class TRefMatch : public TClass {
37 
38 private:
39  // members
40  // NOTE: There is no individual datamember for the best-match angle; it is stored in fResults[0]
41  const char* fTestFileName;
42  const char* fTestTreeName;
43  const char* fTestVarName;
44  TFile* fTestSampleFile;
45  TTree* fTestSampleTree;
46  TBranch* fTestSampleBranch;
47  const char* fReferenceTreeName;
49  TFile* fOutFile;
51  TSystemDirectory* fReferenceFileDir;
52  Long64_t fnTestSampleEvents;
53  Long64_t fnReferenceEvents;
54  Bool_t fkAnderson;
55  Double_t fTestSamplePhiTrue;
56  Double_t fProb;
57  Double_t fSig;
58  TMatrixD fResultsMatrix;
59  TVectorD fResults;
60  TFile* fBestMatchFile;
61  Bool_t fkHasInit;
62  Bool_t fkHasRun;
63  TCanvas* fCanvas;
65  TGraph* fResultsGraph;
66 
67 private:
68  // methods
69  void Init();
70  void SetTestSamplePhiTrue( Double_t phi ) { fTestSamplePhiTrue = phi; }
71  void SetProb( Double_t prob ) { fProb = prob; }
72  void SetSig( Double_t sig ) { fSig = sig; }
73  void SetResultsMatrix( TMatrixD resmat ) { fResultsMatrix = resmat; }
74  void SetResults( TVectorD res ) { fResults = res; }
75  void SetBestMatchFile( TFile* bestMatchFile ) { fBestMatchFile = bestMatchFile; }
76  void SetCanvas( TCanvas *c ) { fCanvas = c; }
77  void SetTestSampleHist( TH1D* h ) { fTestSampleHist = h; }
78  void SetResultsGraph( TGraph *g ) { fResultsGraph = g; }
79 
80 public:
81  // methods
82  // ctors and inits
83  TRefMatch();
84  TRefMatch( const char* fileName, const char* treeName = "T", const char* branchVarName = "phi" );
85  void Init( const char* fileName, const char* treeName, const char* branchVarName );
86  // getters:
87  const char* GetTestFileName() { return fTestFileName; }
88  const char* GetTestTreeName() { return fTestTreeName; }
89  const char* GetTestVarName() { return fTestVarName; }
90  TFile* GetFile() { return fTestSampleFile; }
91  TTree* GetTree() { return fTestSampleTree; }
92  TBranch* GetBranch() { return fTestSampleBranch; }
93  const char* GetReferenceTreeName() { return fReferenceTreeName; }
94  TList* GetReferenceFileList() { return fReferenceFileList; }
95  TFile* GetOutFile() { return fOutFile; }
96  TRegexp GetReferenceFilePattern() { return fReferenceFilePattern; }
97  TSystemDirectory* GetReferenceFileDir() { return fReferenceFileDir; }
98  Long64_t GetnTestSampleEvents() { return fnTestSampleEvents; }
99  Long64_t GetnReferenceEvents() { return fnReferenceEvents; }
100  Bool_t GetAnderson() { return fkAnderson; }
101  Double_t GetTestSamplePhiTrue() { return fTestSamplePhiTrue; }
102  Double_t GetProb() { return fProb; }
103  Double_t GetSig() { return fSig; }
104  TMatrixD GetResultsMatrix() { return fResultsMatrix; }
105  TVectorD GetResults() { return fResults; }
106  TFile* GetBestMatchFile() { return fBestMatchFile; }
107  TCanvas* GetCanvas() { return fCanvas; }
108  TH1D* GetTestSampleHist() { return fTestSampleHist; }
109  TGraph* GetResultsGraph() { return fResultsGraph; }
110  // setters:
111  void SetReferenceFileDir(TSystemDirectory* refFileDir);
112  void SetReferenceFileDir(const char* refFileDirName);
113  void SetReferenceFilePattern(TRegexp patternRE);
114  void SetReferenceFilePattern(const char* pattern);
115  void SetReferenceFileList(TList *fileList) { fReferenceFileList = fileList; }
116  void SetAnderson(Bool_t bAnderson) { fkAnderson = bAnderson; }
117  void SetReferenceTreeName(const char* treeName) { fReferenceTreeName = treeName; }
118  void FillReferenceFileList();
119  void SetnEvents(Long64_t nTestSampleEvents=0, Long64_t nReferenceEvents=0) { fnTestSampleEvents=nTestSampleEvents; fnReferenceEvents=nReferenceEvents; }
120  // utility:
121  Double_t Prob2Sig( Double_t prob );
122  Double_t Sig2Prob( Double_t sig );
123  void ExtractRef( const char* runName, const char* treeName = "T", const char* branchName = "phi" );
124  void ExtractTest( const char* runName, const char* treeName = "T_ncap", const char* branchName = "phi" );
125  // MAIN:
126  Double_t UnbinnedKSTest(TTree *T1, TTree *T2, const char* branchName1, const char* branchName2="", Long64_t nEvents1=0, Long64_t nEvents2=0 );
127  Double_t AndersonDarlingTest(TTree *T1, TTree *T2, const char* branchName1, const char* branchName2="", Long64_t nEvents1=0, Long64_t nEvents2=0 );
128  void RefCompare();
129  // plots:
130  void DrawResults(Bool_t kDrawFit=kFALSE);
131  // print info:
132  void PrintVerbose();
133  void PrintResults();
134  // save and / or close
135  void SaveResults();
136  void Close();
137  // do everything
138  void Run(Long64_t nTestSampleEvents=0);
139 
140 
141 //Integrating the TRefMatch class to ROOT.
142 ClassDef(TRefMatch,4)
143 
144 }; //end class
145 
146 // all pau! )
147 #endif
Class implementing KS/AD-test reference-matching algorithm. //DOC//.
Definition: TRefMatch.h:36
const char * fReferenceTreeName
Name of tree in reference files.
Definition: TRefMatch.h:47
Double_t Sig2Prob(Double_t sig)
Convert significance to probability.
Definition: TRefMatch.cxx:252
TFile * fBestMatchFile
Address of best-match reference file.
Definition: TRefMatch.h:60
void SetReferenceFilePattern(TRegexp patternRE)
Set search pattern (by TRegexp object) for reference filenames.
Definition: TRefMatch.cxx:171
TTree * fTestSampleTree
Address of test-sample tree.
Definition: TRefMatch.h:45
TH1D * fTestSampleHist
TH1 for test-sample data.
Definition: TRefMatch.h:64
TCanvas * fCanvas
Canvas for drawing results.
Definition: TRefMatch.h:63
void SaveResults()
Print canvas and save object.
Definition: TRefMatch.cxx:856
TVectorD fResults
3x1 vector containing best match value and corresponding probability and significance
Definition: TRefMatch.h:59
const char * fTestFileName
Name of file containing data sample to be tested.
Definition: TRefMatch.h:41
Long64_t fnTestSampleEvents
Number of events to use from test-sample tree (default value 0 will use all events)
Definition: TRefMatch.h:52
const char * fTestVarName
Name of branch/variable containing data sample to be tested.
Definition: TRefMatch.h:43
TFile * fOutFile
Output file.
Definition: TRefMatch.h:49
Bool_t fkHasRun
Status indicator for whether RefComp has been called yet.
Definition: TRefMatch.h:62
TRefMatch()
Default ctor.
Definition: TRefMatch.cxx:33
Long64_t fnReferenceEvents
Number of events to use from reference trees (default value 0 will use all events)
Definition: TRefMatch.h:53
Double_t Prob2Sig(Double_t prob)
Convert probability to significance.
Definition: TRefMatch.cxx:229
Bool_t fkAnderson
If true, use Anderson-Darling rather than Kolmogorov-Smirnov.
Definition: TRefMatch.h:54
void SetReferenceFileDir(TSystemDirectory *refFileDir)
Set directory (by pointer) to search for reference files.
Definition: TRefMatch.cxx:143
TBranch * fTestSampleBranch
Address of test-sample branch.
Definition: TRefMatch.h:46
void Run(Long64_t nTestSampleEvents=0)
Execute a typical analysis.
Definition: TRefMatch.cxx:825
void Init()
Initialize: validate and fill members.
Definition: TRefMatch.cxx:77
TMatrixD fResultsMatrix
Matrix of comparison results.
Definition: TRefMatch.h:58
void FillReferenceFileList()
Populate the reference-file list based on the current reference-file directory and pattern.
Definition: TRefMatch.cxx:192
void ExtractRef(const char *runName, const char *treeName="T", const char *branchName="phi")
Extract distribution from reference run.
Definition: TRefMatch.cxx:718
Double_t fProb
Match probability (current entry)
Definition: TRefMatch.h:56
const char * fTestTreeName
Name of tree containing data sample to be tested.
Definition: TRefMatch.h:42
void ExtractTest(const char *runName, const char *treeName="T_ncap", const char *branchName="phi")
Extract distribution from experimental run.
Definition: TRefMatch.cxx:769
TFile * fTestSampleFile
Address of test-sample file.
Definition: TRefMatch.h:44
void PrintVerbose()
Mostly settings.
Definition: TRefMatch.cxx:932
Double_t fSig
Match significance (current entry)
Definition: TRefMatch.h:57
void DrawResults(Bool_t kDrawFit=kFALSE)
Plot sample distribution and algorithm results.
Definition: TRefMatch.cxx:594
Double_t fTestSamplePhiTrue
The true value of phi for the test sample.
Definition: TRefMatch.h:55
Bool_t fkHasInit
Status indicator for whether Init has been called yet.
Definition: TRefMatch.h:61
Double_t AndersonDarlingTest(TTree *T1, TTree *T2, const char *branchName1, const char *branchName2="", Long64_t nEvents1=0, Long64_t nEvents2=0)
Apply (unbinned) Anderson-Darling test.
Definition: TRefMatch.cxx:279
void Close()
Close pads, files, etc.
Definition: TRefMatch.cxx:888
TSystemDirectory * fReferenceFileDir
Directory containing reference files.
Definition: TRefMatch.h:51
void PrintResults()
Results summary.
Definition: TRefMatch.cxx:946
void RefCompare()
Perform reference-comparison algorithm.
Definition: TRefMatch.cxx:477
TGraph * fResultsGraph
TGraph for algorithm results.
Definition: TRefMatch.h:65
Double_t UnbinnedKSTest(TTree *T1, TTree *T2, const char *branchName1, const char *branchName2="", Long64_t nEvents1=0, Long64_t nEvents2=0)
Apply unbinned Kolmogorov-Smirnov test.
Definition: TRefMatch.cxx:363
TRegexp fReferenceFilePattern
Regex describing reference files.
Definition: TRefMatch.h:50
TList * fReferenceFileList
List of files containing reference distributions.
Definition: TRefMatch.h:48