ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/GPetrucc/interface/fwlite/Scanner.h
Revision: 1.9
Committed: Thu Feb 11 11:15:18 2010 UTC (15 years, 3 months ago) by gpetrucc
Content type: text/plain
Branch: MAIN
Changes since 1.8: +3 -3 lines
Log Message:
- Don't own event selectors
  (negligent users will get leaks, but nobody will get crashes)
- fix drawProfile

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