ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/GPetrucc/interface/fwlite/Scanner.h
Revision: 1.7
Committed: Thu Feb 4 10:53:26 2010 UTC (15 years, 3 months ago) by gpetrucc
Content type: text/plain
Branch: MAIN
Changes since 1.6: +95 -14 lines
Log Message:
- "press <CR> to continue" in scan()
- no more warning for existing histo "htemp" in draw
- draw2d without range now works
- added drawGraph

File Contents

# User Rev Content
1 gpetrucc 1.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 gpetrucc 1.7 #include <TDirectory.h>
9 gpetrucc 1.4 #include <TEnv.h>
10 gpetrucc 1.1
11     #if !defined(__CINT__) && !defined(__MAKECINT__)
12     #include "UserCode/GPetrucc/interface/fwliteHelpers.h"
13     #else
14 gpetrucc 1.6 #ifndef UserCode_GPetrucc_CINT_load_library
15     #define UserCode_GPetrucc_CINT_load_library
16 gpetrucc 1.1 // load the library that contains the dictionaries
17     int _load_UserCodeGPetrucc = gSystem->Load("libUserCodeGPetrucc.so");
18     #endif
19 gpetrucc 1.6 #endif
20    
21     #include "UserCode/GPetrucc/interface/fwlite/EventSelectors.h"
22 gpetrucc 1.1
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 gpetrucc 1.4 event_(ev), label_(label), instance_(instance),
30     printFullEventId_(ev->isRealData()),
31 gpetrucc 1.5 ignoreExceptions_(false),
32 gpetrucc 1.7 exprSep_(":"),
33     maxLinesToPrint_(50)
34 gpetrucc 1.1 {
35 gpetrucc 1.6 objType = helper::Parser::elementType(Reflex::Type::ByTypeInfo(HandleT::TempWrapT::typeInfo()));
36     eventSelectors_.SetOwner();
37 gpetrucc 1.1 }
38    
39    
40 gpetrucc 1.7 void scan(const char *exprs, const char *cut="", int nmax=-1) {
41 gpetrucc 1.1 helper::ScannerBase scanner(objType);
42 gpetrucc 1.5 scanner.setIgnoreExceptions(ignoreExceptions_);
43 gpetrucc 1.1
44 gpetrucc 1.4 TObjArray *exprArray = TString(exprs).Tokenize(exprSep_);
45 gpetrucc 1.1 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 gpetrucc 1.2 printf(" : %5s : %3s", "EVENT", "#IT");
51 gpetrucc 1.1 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 gpetrucc 1.3 printf(" : %8s", (strlen(ex)>8 ? ex+strlen(ex)-8 : ex)); // the rightmost part is usually the more interesting one
57 gpetrucc 1.1 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 gpetrucc 1.7 int iev = 0, line = 0;
69 gpetrucc 1.6 for (event_->toBegin(); (iev != nmax) && !event_->atEnd(); ++iev, ++(*event_)) {
70     if (!selectEvent(*event_)) continue;
71 gpetrucc 1.1 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 gpetrucc 1.7 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 gpetrucc 1.1 }
91     }
92     std::cout << std::endl;
93     }
94 gpetrucc 1.7
95 gpetrucc 1.4 TH1 * draw(const char *expr, const char *cut = "", TString drawopt = "", TH1 *hist = 0) {
96     // prep the machinery
97     helper::ScannerBase scanner(objType);
98 gpetrucc 1.5 scanner.setIgnoreExceptions(ignoreExceptions_);
99 gpetrucc 1.7 if (!scanner.addExpression(expr)) return 0;
100 gpetrucc 1.4 if (strlen(cut)) scanner.setCut(cut);
101    
102     // make histo, if needed
103     if (hist == 0) {
104 gpetrucc 1.7 htempDelete();
105 gpetrucc 1.4 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 gpetrucc 1.6 if (!selectEvent(*event_)) continue;
114 gpetrucc 1.4 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 gpetrucc 1.7 htempDelete();
126 gpetrucc 1.4 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 gpetrucc 1.5 scanner.setIgnoreExceptions(ignoreExceptions_);
134 gpetrucc 1.7 if (!scanner.addExpression(xexpr)) return 0;
135     if (!scanner.addExpression(yexpr)) return 0;
136 gpetrucc 1.4 if (strlen(cut)) scanner.setCut(cut);
137    
138     // make histo, if needed
139     if (hist == 0) {
140 gpetrucc 1.7 htempDelete();
141 gpetrucc 1.4 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 gpetrucc 1.6 if (!selectEvent(*event_)) continue;
150 gpetrucc 1.4 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 gpetrucc 1.7 htempDelete();
162 gpetrucc 1.4 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 gpetrucc 1.5 scanner.setIgnoreExceptions(ignoreExceptions_);
170 gpetrucc 1.7 if (!scanner.addExpression((const char *)xexpr)) return 0;
171     if (!scanner.addExpression((const char *)yexpr)) return 0;
172 gpetrucc 1.4 if (strlen(cut)) scanner.setCut(cut);
173    
174     // make histo, if needed
175     if (hist == 0) {
176 gpetrucc 1.7 // 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 gpetrucc 1.4 hist = new TH2F("htemp",
194     (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr),
195 gpetrucc 1.7 gEnv->GetValue("Hist.Binning.2D.x",20), xmin, xmax,
196     gEnv->GetValue("Hist.Binning.2D.y",20), ymin, ymax);
197 gpetrucc 1.4 }
198    
199     // fill histogram
200     for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
201 gpetrucc 1.6 if (!selectEvent(*event_)) continue;
202 gpetrucc 1.4 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 gpetrucc 1.7 htempDelete();
218 gpetrucc 1.4 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 gpetrucc 1.7
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 gpetrucc 1.1 void setPrintFullEventId(bool printIt=true) { printFullEventId_ = printIt; }
256 gpetrucc 1.4 void setExpressionSeparator(TString separator) { exprSep_ = separator; }
257 gpetrucc 1.5 void setIgnoreExceptions(bool ignoreThem) { ignoreExceptions_ = ignoreThem; }
258 gpetrucc 1.7 void setMaxLinesToPrint(int lines) { maxLinesToPrint_ = (lines > 0 ? lines : 2147483647); }
259    
260 gpetrucc 1.6 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 gpetrucc 1.1 private:
270 gpetrucc 1.2 fwlite::EventBase *event_;
271     std::string label_, instance_, process_;
272 gpetrucc 1.1 bool printFullEventId_;
273 gpetrucc 1.5 bool ignoreExceptions_;
274 gpetrucc 1.4 TString exprSep_;
275 gpetrucc 1.1 HandleT handle_;
276 gpetrucc 1.4 Reflex::Type objType;
277 gpetrucc 1.1
278 gpetrucc 1.6 TObjArray eventSelectors_;
279    
280 gpetrucc 1.7 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 gpetrucc 1.1 };
300     }