1 |
// these includes are FWLite-safe
|
2 |
#include "DataFormats/FWLite/interface/Handle.h"
|
3 |
#include "DataFormats/FWLite/interface/Event.h"
|
4 |
// these are from ROOT, so they're safe too
|
5 |
#include <TString.h>
|
6 |
#include <TObjString.h>
|
7 |
#include <TObjArray.h>
|
8 |
#include <TDirectory.h>
|
9 |
#include <TEnv.h>
|
10 |
|
11 |
#if !defined(__CINT__) && !defined(__MAKECINT__)
|
12 |
#include "UserCode/GPetrucc/interface/fwliteHelpers.h"
|
13 |
#else
|
14 |
#ifndef UserCode_GPetrucc_CINT_load_library
|
15 |
#define UserCode_GPetrucc_CINT_load_library
|
16 |
// load the library that contains the dictionaries
|
17 |
int _load_UserCodeGPetrucc = gSystem->Load("libUserCodeGPetrucc.so");
|
18 |
#endif
|
19 |
#endif
|
20 |
|
21 |
#include "UserCode/GPetrucc/interface/fwlite/EventSelectors.h"
|
22 |
|
23 |
namespace fwlite {
|
24 |
template<typename T>
|
25 |
class Scanner {
|
26 |
public:
|
27 |
typedef fwlite::Handle<std::vector<T> > HandleT;
|
28 |
Scanner(fwlite::EventBase *ev, const char *label, const char *instance = "", const char *process="") :
|
29 |
event_(ev), label_(label), instance_(instance),
|
30 |
printFullEventId_(ev->isRealData()),
|
31 |
ignoreExceptions_(false),
|
32 |
exprSep_(":"),
|
33 |
maxLinesToPrint_(50)
|
34 |
{
|
35 |
objType = helper::Parser::elementType(Reflex::Type::ByTypeInfo(HandleT::TempWrapT::typeInfo()));
|
36 |
eventSelectors_.SetOwner();
|
37 |
}
|
38 |
|
39 |
|
40 |
void scan(const char *exprs, const char *cut="", int nmax=-1) {
|
41 |
helper::ScannerBase scanner(objType);
|
42 |
scanner.setIgnoreExceptions(ignoreExceptions_);
|
43 |
|
44 |
TObjArray *exprArray = TString(exprs).Tokenize(exprSep_);
|
45 |
int rowline = 0;
|
46 |
if (printFullEventId_) {
|
47 |
printf(" : %9s : %4s : %9s : %3s", "RUN", "LUMI", "EVENT", "#IT");
|
48 |
rowline += 3*4+9+4+9+3-1; // -1 as first char remain blank
|
49 |
} else {
|
50 |
printf(" : %5s : %3s", "EVENT", "#IT");
|
51 |
rowline += 3+6+3+3-1; // -1 as first char remain blank
|
52 |
}
|
53 |
for (int i = 0; i < exprArray->GetEntries(); ++i) {
|
54 |
const char *ex = ((TObjString *)(*exprArray)[i])->GetString();
|
55 |
scanner.addExpression(ex);
|
56 |
printf(" : %8s", (strlen(ex)>8 ? ex+strlen(ex)-8 : ex)); // the rightmost part is usually the more interesting one
|
57 |
rowline += 3+8;
|
58 |
}
|
59 |
std::cout << " :" << std::endl;
|
60 |
rowline += 2;
|
61 |
delete exprArray;
|
62 |
|
63 |
TString rule('-', rowline);
|
64 |
std::cout << " " << rule << " " << std::endl;
|
65 |
|
66 |
if (strlen(cut)) scanner.setCut(cut);
|
67 |
|
68 |
int iev = 0, line = 0;
|
69 |
for (event_->toBegin(); (iev != nmax) && !event_->atEnd(); ++iev, ++(*event_)) {
|
70 |
if (!selectEvent(*event_)) continue;
|
71 |
handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
|
72 |
const std::vector<T> & vals = *handle_;
|
73 |
for (size_t j = 0, n = vals.size(); j < n; ++j) {
|
74 |
if (!scanner.test(&vals[j])) continue;
|
75 |
if (printFullEventId_) {
|
76 |
const edm::EventAuxiliary &id = event_->eventAuxiliary();
|
77 |
printf(" : %9d : %4d : %9d : %3d", id.run(), id.luminosityBlock(), id.event(), j);
|
78 |
} else {
|
79 |
printf(" : %5d : %3d", iev, j);
|
80 |
}
|
81 |
scanner.print(&vals[j]);
|
82 |
std::cout << " :" << std::endl;
|
83 |
if (++line == maxLinesToPrint_) {
|
84 |
line = 0;
|
85 |
if (!wantMore()) {
|
86 |
iev = nmax-1; // this is to exit the outer loop
|
87 |
break; // and this to exit the inner one
|
88 |
}
|
89 |
}
|
90 |
}
|
91 |
}
|
92 |
std::cout << std::endl;
|
93 |
}
|
94 |
|
95 |
TH1 * draw(const char *expr, const char *cut = "", TString drawopt = "", TH1 *hist = 0) {
|
96 |
// prep the machinery
|
97 |
helper::ScannerBase scanner(objType);
|
98 |
scanner.setIgnoreExceptions(ignoreExceptions_);
|
99 |
if (!scanner.addExpression(expr)) return 0;
|
100 |
if (strlen(cut)) scanner.setCut(cut);
|
101 |
|
102 |
// make histo, if needed
|
103 |
if (hist == 0) {
|
104 |
htempDelete();
|
105 |
hist = new TH1F("htemp",
|
106 |
(strlen(cut) ? TString(expr)+"{"+cut+"}" : TString(expr)),
|
107 |
gEnv->GetValue("Hist.Binning.1D.x",100), 0, 0);
|
108 |
hist->SetBit(TH1::kCanRebin);
|
109 |
}
|
110 |
|
111 |
// fill histogram
|
112 |
for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
|
113 |
if (!selectEvent(*event_)) continue;
|
114 |
handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
|
115 |
const std::vector<T> & vals = *handle_;
|
116 |
for (size_t j = 0, n = vals.size(); j < n; ++j) {
|
117 |
scanner.fill1D(&vals[j], hist);
|
118 |
}
|
119 |
}
|
120 |
|
121 |
if (!TString(drawopt).Contains("goff",TString::kIgnoreCase)) hist->Draw(drawopt);
|
122 |
return hist;
|
123 |
}
|
124 |
TH1 * draw(const char *expr, int bins, double xlow, double xhigh, const char *cut = "", const char *drawopt = "") {
|
125 |
htempDelete();
|
126 |
TH1 * htemp = new TH1F("htemp", (strlen(cut) ? TString(expr)+"{"+cut+"}" : TString(expr)), bins, xlow, xhigh);
|
127 |
return draw(expr,cut,drawopt,htemp);
|
128 |
}
|
129 |
|
130 |
TProfile * drawProf(TString xexpr, TString yexpr, const char *cut = "", TString drawopt = "", TProfile *hist = 0) {
|
131 |
// prep the machinery
|
132 |
helper::ScannerBase scanner(objType);
|
133 |
scanner.setIgnoreExceptions(ignoreExceptions_);
|
134 |
if (!scanner.addExpression(xexpr)) return 0;
|
135 |
if (!scanner.addExpression(yexpr)) return 0;
|
136 |
if (strlen(cut)) scanner.setCut(cut);
|
137 |
|
138 |
// make histo, if needed
|
139 |
if (hist == 0) {
|
140 |
htempDelete();
|
141 |
hist = new TProfile("htemp",
|
142 |
(strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr),
|
143 |
gEnv->GetValue("Hist.Binning.1D.x",100), 0., 0.);
|
144 |
hist->SetBit(TH1::kCanRebin);
|
145 |
}
|
146 |
|
147 |
// fill histogram
|
148 |
for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
|
149 |
if (!selectEvent(*event_)) continue;
|
150 |
handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
|
151 |
const std::vector<T> & vals = *handle_;
|
152 |
for (size_t j = 0, n = vals.size(); j < n; ++j) {
|
153 |
scanner.fillProf(&vals[j], hist);
|
154 |
}
|
155 |
}
|
156 |
|
157 |
if (!TString(drawopt).Contains("goff",TString::kIgnoreCase)) hist->Draw(drawopt);
|
158 |
return hist;
|
159 |
}
|
160 |
TProfile * drawProf(TString xexpr, int bins, double xlow, double xhigh, TString yexpr, const char *cut = "", const char *drawopt = "") {
|
161 |
htempDelete();
|
162 |
TProfile * htemp = new TProfile("htemp", (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr), bins, xlow, xhigh);
|
163 |
return drawProf(xexpr,yexpr,cut,drawopt,htemp);
|
164 |
}
|
165 |
|
166 |
TH2 * draw2D(TString xexpr, TString yexpr, const char *cut = "", TString drawopt = "", TH2 *hist = 0) {
|
167 |
// prep the machinery
|
168 |
helper::ScannerBase scanner(objType);
|
169 |
scanner.setIgnoreExceptions(ignoreExceptions_);
|
170 |
if (!scanner.addExpression((const char *)xexpr)) return 0;
|
171 |
if (!scanner.addExpression((const char *)yexpr)) return 0;
|
172 |
if (strlen(cut)) scanner.setCut(cut);
|
173 |
|
174 |
// make histo, if needed
|
175 |
if (hist == 0) {
|
176 |
// ok this is much more a hack than for the 1D case
|
177 |
double xmin = 0, xmax = -1, ymin = 0, ymax = -1;
|
178 |
for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
|
179 |
if (!selectEvent(*event_)) continue;
|
180 |
handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
|
181 |
const std::vector<T> & vals = *handle_;
|
182 |
for (size_t j = 0, n = vals.size(); j < n; ++j) {
|
183 |
if (!scanner.test(&vals[j])) continue;
|
184 |
double x = scanner.eval(&vals[j],0);
|
185 |
double y = scanner.eval(&vals[j],1);
|
186 |
if ((xmax == -1) || (x >= xmax)) xmax = x;
|
187 |
if ((xmin == 0) || (x <= xmin)) xmin = x;
|
188 |
if ((ymax == -1) || (y >= ymax)) ymax = y;
|
189 |
if ((ymin == 0) || (y <= ymin)) ymin = y;
|
190 |
}
|
191 |
}
|
192 |
htempDelete();
|
193 |
hist = new TH2F("htemp",
|
194 |
(strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr),
|
195 |
gEnv->GetValue("Hist.Binning.2D.x",20), xmin, xmax,
|
196 |
gEnv->GetValue("Hist.Binning.2D.y",20), ymin, ymax);
|
197 |
}
|
198 |
|
199 |
// fill histogram
|
200 |
for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
|
201 |
if (!selectEvent(*event_)) continue;
|
202 |
handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
|
203 |
const std::vector<T> & vals = *handle_;
|
204 |
for (size_t j = 0, n = vals.size(); j < n; ++j) {
|
205 |
scanner.fill2D(&vals[j], hist);
|
206 |
}
|
207 |
}
|
208 |
|
209 |
if (!strlen(hist->GetXaxis()->GetTitle())) hist->GetXaxis()->SetTitle(xexpr);
|
210 |
if (!strlen(hist->GetYaxis()->GetTitle())) hist->GetYaxis()->SetTitle(yexpr);
|
211 |
if (!TString(drawopt).Contains("goff",TString::kIgnoreCase)) hist->Draw(drawopt);
|
212 |
return hist;
|
213 |
}
|
214 |
TH2 * draw2D(TString xexpr, int xbins, double xlow, double xhigh,
|
215 |
TString yexpr, int ybins, double ylow, double yhigh,
|
216 |
const char *cut = "", const char *drawopt = "") {
|
217 |
htempDelete();
|
218 |
TH2 * htemp = new TH2F("htemp", (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr),
|
219 |
xbins, xlow, xhigh, ybins,ylow,yhigh);
|
220 |
return draw2D(xexpr,yexpr,cut,drawopt,htemp);
|
221 |
}
|
222 |
|
223 |
|
224 |
TGraph * drawGraph(TString xexpr, TString yexpr, const char *cut = "", TString drawopt = "AP", TGraph *graph = 0) {
|
225 |
// prep the machinery
|
226 |
helper::ScannerBase scanner(objType);
|
227 |
scanner.setIgnoreExceptions(ignoreExceptions_);
|
228 |
if (!scanner.addExpression((const char *)xexpr)) return 0;
|
229 |
if (!scanner.addExpression((const char *)yexpr)) return 0;
|
230 |
if (strlen(cut)) scanner.setCut(cut);
|
231 |
|
232 |
// make graph, if needed
|
233 |
if (graph == 0) {
|
234 |
htempDelete();
|
235 |
graph = new TGraph();
|
236 |
graph->SetNameTitle("htemp", (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
|
237 |
}
|
238 |
|
239 |
// fill graph
|
240 |
for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
|
241 |
if (!selectEvent(*event_)) continue;
|
242 |
handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
|
243 |
const std::vector<T> & vals = *handle_;
|
244 |
for (size_t j = 0, n = vals.size(); j < n; ++j) {
|
245 |
scanner.fillGraph(&vals[j], graph);
|
246 |
}
|
247 |
}
|
248 |
|
249 |
if (!strlen(graph->GetXaxis()->GetTitle())) graph->GetXaxis()->SetTitle(xexpr);
|
250 |
if (!strlen(graph->GetYaxis()->GetTitle())) graph->GetYaxis()->SetTitle(yexpr);
|
251 |
if (!TString(drawopt).Contains("goff",TString::kIgnoreCase)) graph->Draw(drawopt);
|
252 |
return graph;
|
253 |
}
|
254 |
|
255 |
void setPrintFullEventId(bool printIt=true) { printFullEventId_ = printIt; }
|
256 |
void setExpressionSeparator(TString separator) { exprSep_ = separator; }
|
257 |
void setIgnoreExceptions(bool ignoreThem) { ignoreExceptions_ = ignoreThem; }
|
258 |
void setMaxLinesToPrint(int lines) { maxLinesToPrint_ = (lines > 0 ? lines : 2147483647); }
|
259 |
|
260 |
void addEventSelector(fwlite::EventSelector *selector) { eventSelectors_.Add(selector); }
|
261 |
void clearEventSelector() { eventSelectors_.Clear(); }
|
262 |
TObjArray & eventSelectors() { return eventSelectors_; }
|
263 |
bool selectEvent(const fwlite::EventBase &ev) const {
|
264 |
for (int i = 0, n = eventSelectors_.GetEntries(); i < n; ++i) {
|
265 |
if (!((fwlite::EventSelector *)(eventSelectors_[i]))->accept(ev)) return false;
|
266 |
}
|
267 |
return true;
|
268 |
}
|
269 |
private:
|
270 |
fwlite::EventBase *event_;
|
271 |
std::string label_, instance_, process_;
|
272 |
bool printFullEventId_;
|
273 |
bool ignoreExceptions_;
|
274 |
TString exprSep_;
|
275 |
HandleT handle_;
|
276 |
Reflex::Type objType;
|
277 |
|
278 |
TObjArray eventSelectors_;
|
279 |
|
280 |
int maxLinesToPrint_;
|
281 |
bool wantMore() const {
|
282 |
// ask if user wants more
|
283 |
fprintf(stderr,"Type <CR> to continue or q to quit ==> ");
|
284 |
// read first char
|
285 |
int readch = getchar(), answer = readch;
|
286 |
// poll out remaining chars from buffer
|
287 |
while (readch != '\n' && readch != EOF) readch = getchar();
|
288 |
// check first char
|
289 |
return !(answer == 'q' || answer == 'Q');
|
290 |
}
|
291 |
|
292 |
void htempDelete() {
|
293 |
if (gDirectory) {
|
294 |
TObject *obj = gDirectory->Get("htemp");
|
295 |
if (obj) obj->Delete();
|
296 |
}
|
297 |
}
|
298 |
|
299 |
};
|
300 |
}
|