Jpp test-rotations-old
the software that should make you happy
Loading...
Searching...
No Matches
JScale1D.cc File Reference

Auxiliary program to scale ROOT histograms. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <sstream>
#include <cmath>
#include "TROOT.h"
#include "TFile.h"
#include "TObject.h"
#include "TKey.h"
#include "TH1.h"
#include "TGraph.h"
#include "TGraphErrors.h"
#include "TString.h"
#include "TRegexp.h"
#include "JTools/JRange.hh"
#include "JGizmo/JRootObjectID.hh"
#include "JGizmo/JGizmoToolkit.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 scale ROOT histograms.

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

The scale factor (option -F <factor>) refers to a ROOT TFormula. The expression may contain member methods of the corresponding object.

Author
mdejong

Definition in file JScale1D.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 35 of file JScale1D.cc.

36{
37 using namespace std;
38 using namespace JPP;
39
40 typedef JRange<double> JRange_t;
41
42 vector<JRootObjectID> inputFile;
43 string outputFile;
44 TString formula;
45 JRange_t X;
46 JRange_t Y;
47 string option;
48 int debug;
49
50 try {
51
52 JParser<> zap("Auxiliary program to scale ROOT histograms."\
53 "\nNote that the formula may contain method names of the specified object.");
54
55 zap['f'] = make_field(inputFile, "<input file>:<object name>");
56 zap['o'] = make_field(outputFile, "ROOT file with histogram(s)") = "scale.root";
57 zap['F'] = make_field(formula, "ROOT TFormula (may contain method names of object)");
58 zap['x'] = make_field(X, "abscissa range (only for Integral)") = JRange_t();
59 zap['y'] = make_field(Y, "abscissa range (only for Integral)") = JRange_t();
60 zap['O'] = make_field(option) = "", "nosw2", "width", "nosw2 width";
61 zap['d'] = make_field(debug) = 0;
62
63 zap(argc, argv);
64 }
65 catch(const exception &error) {
66 FATAL(error.what() << endl);
67 }
68
69
70 TFile out(outputFile.c_str(), "recreate");
71
72 for (vector<JRootObjectID>::const_iterator input = inputFile.begin(); input != inputFile.end(); ++input) {
73
74 DEBUG("Input: " << *input << endl);
75
76 TDirectory* dir = getDirectory(*input);
77
78 if (dir == NULL) {
79 ERROR("File: " << input->getFullFilename() << " not opened." << endl);
80 continue;
81 }
82
83 const TRegexp regexp(input->getObjectName());
84
85 TIter iter(dir->GetListOfKeys());
86
87 for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
88
89 const TString tag(key->GetName());
90
91 DEBUG("Key: " << tag << " match = " << tag.Contains(regexp) << endl);
92
93 // option match
94
95 if (tag.Contains(regexp) && isTObject(key)) {
96
97 TObject* object = key->ReadObj();
98
99 TH1* h1 = dynamic_cast<TH1*> (object);
100 TGraphErrors* ge = dynamic_cast<TGraphErrors*>(object);
101 TGraph* g1 = dynamic_cast<TGraph*> (object);
102
103 if (h1 != NULL ||
104 ge != NULL ||
105 g1 != NULL) {
106
107 const double factor = getResult(formula, object);
108
109 if (h1 != NULL) {
110
111 if (X != JRange_t()) { h1->GetXaxis()->SetRangeUser(X.getLowerLimit(), X.getUpperLimit()); }
112 if (Y != JRange_t()) { h1->GetYaxis()->SetRangeUser(Y.getLowerLimit(), Y.getUpperLimit()); }
113
114 h1->Scale(factor, option.c_str());
115
116 } else if (ge != NULL) {
117
118 for (Int_t i = 0; i != ge->GetN(); ++i) {
119 ge->GetY() [i] *= factor;
120 ge->GetEY()[i] *= factor;
121 }
122
123 } else if (g1 != NULL) {
124
125 for (Int_t i = 0; i != g1->GetN(); ++i) {
126 g1->GetY() [i] *= factor;
127 }
128 }
129
130 out.WriteTObject(object);
131 }
132 }
133 }
134 }
135
136 out.Write();
137 out.Close();
138}
string outputFile
#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
Double_t g1(const Double_t x)
Function.
Definition JQuantiles.cc:25
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.
Double_t getResult(const TString &text, TObject *object=NULL)
Get result of given textual formula.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).