Jpp 19.3.0-rc.1
the software that should make you happy
Loading...
Searching...
No Matches
JPrintRange2D.cc File Reference

Auxiliary program to print the ranges of x, y and z values of 2D ROOT objects. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TClass.h"
#include "TKey.h"
#include "TRegexp.h"
#include "TH2.h"
#include "TGraph2D.h"
#include "TProfile2D.h"
#include "JTools/JRange.hh"
#include "JGizmo/JRootObjectID.hh"
#include "JGizmo/JGizmoToolkit.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Auxiliary program to print the ranges of x, y and z values of 2D ROOT objects.

The option -f corresponds to <file name>:<object name>.

Author
mdejong

Definition in file JPrintRange2D.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 29 of file JPrintRange2D.cc.

30{
31 using namespace std;
32 using namespace JPP;
33
34 typedef JRange<double> JRange_t;
35
36 vector<JRootObjectID> inputFile;
37 JRange_t x;
38 JRange_t y;
39 JRange_t z;
40 int debug;
41
42 try {
43
44 JParser<> zap("Auxiliary program to print the ranges of x, y and z values of 2D ROOT objects.");
45
46 zap['f'] = make_field(inputFile, "<input file>:<object name>");
47 zap['x'] = make_field(x, "x range") = JRange_t();
48 zap['y'] = make_field(y, "y range") = JRange_t();
49 zap['z'] = make_field(z, "z range") = JRange_t();
50 zap['d'] = make_field(debug) = 0;
51
52 zap(argc, argv);
53 }
54 catch(const exception &error) {
55 FATAL(error.what() << endl);
56 }
57
58
59 JRange_t X = JRange_t::DEFAULT_RANGE();
60 JRange_t Y = JRange_t::DEFAULT_RANGE();
61 JRange_t Z = JRange_t::DEFAULT_RANGE();
62
63 for (vector<JRootObjectID>::const_iterator input = inputFile.begin(); input != inputFile.end(); ++input) {
64
65 DEBUG("Input: " << *input << endl);
66
67 TDirectory* dir = getDirectory(*input);
68
69 if (dir == NULL) {
70 ERROR("File: " << input->getFullFilename() << " not opened." << endl);
71 continue;
72 }
73
74 const TRegexp regexp(input->getObjectName());
75
76 TIter iter(dir->GetListOfKeys());
77
78 for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
79
80 const TString tag(key->GetName());
81
82 DEBUG("Key: " << tag << " match = " << tag.Contains(regexp) << endl);
83
84 // option match
85
86 if (tag.Contains(regexp) && isTObject(key)) {
87
88 TObject* object = key->ReadObj();
89
90 try {
91
92 TH2& h2 = dynamic_cast<TH2&>(*object);
93
94 X.combine(x.join(JRange_t(h2.GetXaxis()->GetXmin(), h2.GetXaxis()->GetXmax())));
95 Y.combine(y.join(JRange_t(h2.GetYaxis()->GetXmin(), h2.GetYaxis()->GetXmax())));
96 Z.combine(z.join(JRange_t(h2.GetMinimum(), h2.GetMaximum())));
97 }
98 catch(exception&) {}
99
100 try {
101
102 TProfile2D& h2 = dynamic_cast<TProfile2D&>(*object);
103
104 X.combine(x.join(JRange_t(h2.GetXaxis()->GetXmin(), h2.GetXaxis()->GetXmax())));
105 Y.combine(y.join(JRange_t(h2.GetYaxis()->GetXmin(), h2.GetYaxis()->GetXmax())));
106 Z.combine(z.join(JRange_t(h2.GetMinimum(), h2.GetMaximum())));
107 }
108 catch(exception&) {}
109
110 try {
111
112 TGraph2D& g2 = dynamic_cast<TGraph2D&>(*object);
113
114 for (Int_t i = 0; i != g2.GetN(); ++i) {
115 if (x(g2.GetX()[i]) &&
116 y(g2.GetY()[i]) &&
117 z(g2.GetZ()[i])) {
118 X.include(g2.GetX()[i]);
119 Y.include(g2.GetY()[i]);
120 Z.include(g2.GetZ()[i]);
121 }
122 }
123 }
124 catch(exception&) {}
125 }
126 }
127 }
128
129
130 cout << SCIENTIFIC(15,5) << X.getLowerLimit() << ' '
131 << SCIENTIFIC(15,5) << Y.getLowerLimit() << ' '
132 << SCIENTIFIC(15,5) << Z.getLowerLimit() << ' '
133 << SCIENTIFIC(15,5) << X.getUpperLimit() << ' '
134 << SCIENTIFIC(15,5) << Y.getUpperLimit() << ' '
135 << SCIENTIFIC(15,5) << Z.getUpperLimit() << endl;
136}
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define ERROR(A)
Definition JMessage.hh:66
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Utility class to parse command line options.
Definition JParser.hh:1698
Range of values.
Definition JRange.hh:42
TDirectory * getDirectory(const JRootObjectID &id)
Get TDirectory pointer.
bool isTObject(const TKey *key)
Check if given key corresponds to a TObject.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Auxiliary data structure for floating point format specification.
Definition JManip.hh:488