Jpp 20.0.0-rc.2
the software that should make you happy
Loading...
Searching...
No Matches
JZebraMantis.cc File Reference

Program to compare histograms in root files that have same directory structure,
and where the histograms have the same names. More...

#include <iostream>
#include <fstream>
#include "TString.h"
#include "TRegexp.h"
#include "TObjArray.h"
#include "TObjString.h"
#include "TFile.h"
#include "TKey.h"
#include "Jeep/JParser.hh"
#include "JSupport/JMeta.hh"
#include "JLang/JPredicate.hh"
#include "JLang/JVectorize.hh"
#include "JCompareHistograms/JTest_t.hh"
#include "JCompareHistograms/JTestDictionary.hh"
#include "JGizmo/JRootObjectID.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Program to compare histograms in root files that have same directory structure,
and where the histograms have the same names.

The input histograms and the test to be applied for each histogram are specified
in an ASCII formatted steering file which is passed by the command line.

Each row of the steering file should have multiple columns.

Column 1 is the name of the histogram to be compared (including the full path inside the root file)
Column 2 is an integer value that indicates the test to be performed. See JTestDictionary()
Columns 3..n are reserved for the different parameters of the test.

Author
rgruiz, bjung

Definition in file JZebraMantis.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 83 of file JZebraMantis.cc.

83 {
84
85 using namespace JPP;
86 using namespace std;
87
88 string steeringFile;
89
90 string file1;
91 string file2;
92
93 string output;
94 string ascii;
95
96 vector<string> keys;
97
98 bool onlyFailures;
99
100 int debug;
101
102 const array_type<string> listOfKeys = get_keys(JTestSummary().getProperties());
103
104 try {
105
106 const string keysExplainer = MAKE_STRING("Terminal output:" << endl << listOfKeys);
107
108 JParser<> zap("\nProgram to compare histograms in root files that have the same directory structure. See the link below this usage for further details.\n");
109
110 zap['s'] = make_field(steeringFile , "ASCII steering file with list of histograms and tests");
111 zap['a'] = make_field(file1 , "input file 1");
112 zap['b'] = make_field(file2 , "input file 2");
113 zap['o'] = make_field(output , "output file root") = "zebramantis.root";
114 zap['t'] = make_field(ascii , "output file txt" ) = "zebramantis.txt";
115 zap['k'] = make_field(keys , keysExplainer ) = JPARSER::initialised();
116 zap['w'] = make_field(onlyFailures , "write only failed tests" );
117 zap['d'] = make_field(debug) = 2;
118
119 zap(argc,argv);
120 }
121 catch(const exception &error) {
122 ERROR(error.what() << endl);
123 }
124
125 if (keys.empty()) {
126 keys = listOfKeys;
127 }
128
129 TFile* f1 = TFile::Open(file1.c_str());
130 TFile* f2 = TFile::Open(file2.c_str());
131
132 TFile out(output.c_str(),"recreate");
133
134 ofstream results;
135 results.open (ascii);
136 results << "# " << listOfKeys << endl;
137
138 const vector<JRootObjectID> objectIDs = readDir(f1);
139
140 std::ifstream infile(steeringFile);
141
143
144 size_t npassed = 0;
145 size_t nfailed = 0;
146
147 for (string line; getline(infile, line); ) {
148
149 istringstream iss(line);
150
151 TString name;
152 int testID;
153
154 if (!(iss >> name >> testID)) {
155 continue;
156 }
157
158 DEBUG("Input: " << name << ' ' << testID << endl);
159
160 const TRegexp regexp(name);
161
162 for (vector<JRootObjectID>::const_iterator objectID = objectIDs.cbegin() ; objectID != objectIDs.cend() ; ++objectID) {
163
164 const TString& dirName = objectID->getDirectory();
165 const TString& fullName = objectID->getFullObjectName();
166
167 DEBUG("Key: " << fullName << " match = " << fullName.Contains(regexp) << endl);
168
169 if ((fullName.Index(regexp) != -1)) {
170
171 TObject* obj1 = (TObject*)f1->Get(fullName);
172 TObject* obj2 = (TObject*)f2->Get(fullName);
173
174 if (!obj1 || !obj2) {
175 DEBUG("Could not retrieve " << fullName << endl);
176 continue;
177 }
178
179 d[testID]->read(iss);
180 d[testID]->test(obj1,obj2);
181
182 if (dirName.Length() > 0 && !out.GetDirectory(dirName)) {
183
184 if (dirName[0] == JRootObjectID::PATHNAME_SEPARATOR) { // Remove leading forward slash
185 out.mkdir(TString(dirName(1, dirName.Length() - 1)));
186 } else {
187 out.mkdir(dirName);
188 }
189 }
190
191 out.cd(dirName);
192
193 for (vector<JTestResult>::iterator r = d[testID]->begin() ; r != d[testID]->end() ; ++r) {
194
195 if (onlyFailures && r->passed) {
196 continue;
197 }
198
199 print(cout, *r, keys.cbegin(), keys.cend(), false);
200 print(results, *r, listOfKeys.cbegin(), listOfKeys.cend(), true);
201
202 r->obj->Write();
203 }
204
205 const size_t Npass = count_if(d[testID]->cbegin(), d[testID]->cend(),
206 make_predicate(&JTestResult::passed, true));
207
208 npassed += Npass;
209 nfailed += (d[testID]->size() - Npass);
210
211 d[testID]->clear();
212 }
213 }
214 }
215
216 infile.close();
217
218 results << WHITE << "# PASSED: " << npassed << " " << " FAILED: " << nfailed << " FAILURE FRACTION: " << float (nfailed)/(nfailed+npassed) << endl;
219
220 putObject(&out, JMeta(argc, argv));
221 JMeta::copy(file1.c_str(), out);
222 JMeta::copy(file2.c_str(), out);
223
224 results.close();
225 out .Close();
226
227 return 0;
228}
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define ERROR(A)
Definition JMessage.hh:66
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
void readDir(TDirectory *dir, std::vector< TString > &v)
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
void print(const TH1 &h1, std::ostream &out)
Print histogram parameters.
Dictionary to map different tests to unique integer indices.
Class dedicated to standardize the title of the graphical objects produced by the JTest_t() derived c...
Utility class to parse command line options.
Definition JParser.hh:1698
std::istream & getline(std::istream &in, JString &object)
Read string from input stream until end of line.
Definition JString.hh:478
const array_type< JKey_t > & get_keys(const std::map< JKey_t, JValue_t, JComparator_t, JAllocator_t > &data)
Method to create array of keys of map.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
Auxiliary data structure for return type of make methods.
Definition JVectorize.hh:28
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
Auxiliary class for ROOT I/O of application specific meta data.
Definition JMeta.hh:72