Jpp master_rocky-44-g75b7c4f75
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 "TApplication.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 30 of file JPrintRange2D.cc.

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