ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/GPetrucc/interface/fwlite/Scanner.h
Revision: 1.5
Committed: Sun Jan 31 19:44:41 2010 UTC (15 years, 3 months ago) by gpetrucc
Content type: text/plain
Branch: MAIN
Changes since 1.4: +7 -0 lines
Log Message:
allow ignoring exceptions when scanning

File Contents

# Content
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 <TEnv.h>
9
10 #if !defined(__CINT__) && !defined(__MAKECINT__)
11 #include "UserCode/GPetrucc/interface/fwliteHelpers.h"
12 #else
13 // load the library that contains the dictionaries
14 int _load_UserCodeGPetrucc = gSystem->Load("libUserCodeGPetrucc.so");
15 #endif
16
17 namespace fwlite {
18 template<typename T>
19 class Scanner {
20 public:
21 typedef fwlite::Handle<std::vector<T> > HandleT;
22 Scanner(fwlite::EventBase *ev, const char *label, const char *instance = "", const char *process="") :
23 //helper::ScannerBase(Reflex::Type::ByTypeInfo(typeid(T))),
24 event_(ev), label_(label), instance_(instance),
25 printFullEventId_(ev->isRealData()),
26 ignoreExceptions_(false),
27 exprSep_(":")
28 {
29 Reflex::Type wrapperType = Reflex::Type::ByTypeInfo(HandleT::TempWrapT::typeInfo());
30 objType = helpme.elementType(Reflex::Type::ByTypeInfo(HandleT::TempWrapT::typeInfo()));
31 }
32
33
34 void scan(const char *exprs, const char *cut="", int nmax=20) {
35 helper::ScannerBase scanner(objType);
36 scanner.setIgnoreExceptions(ignoreExceptions_);
37
38 TObjArray *exprArray = TString(exprs).Tokenize(exprSep_);
39 int rowline = 0;
40 if (printFullEventId_) {
41 printf(" : %9s : %4s : %9s : %3s", "RUN", "LUMI", "EVENT", "#IT");
42 rowline += 3*4+9+4+9+3-1; // -1 as first char remain blank
43 } else {
44 printf(" : %5s : %3s", "EVENT", "#IT");
45 rowline += 3+6+3+3-1; // -1 as first char remain blank
46 }
47 for (int i = 0; i < exprArray->GetEntries(); ++i) {
48 const char *ex = ((TObjString *)(*exprArray)[i])->GetString();
49 scanner.addExpression(ex);
50 printf(" : %8s", (strlen(ex)>8 ? ex+strlen(ex)-8 : ex)); // the rightmost part is usually the more interesting one
51 rowline += 3+8;
52 }
53 std::cout << " :" << std::endl;
54 rowline += 2;
55 delete exprArray;
56
57 TString rule('-', rowline);
58 std::cout << " " << rule << " " << std::endl;
59
60 if (strlen(cut)) scanner.setCut(cut);
61
62 int iev = 0;
63 for (event_->toBegin(); (iev < nmax) && !event_->atEnd(); ++iev, ++(*event_)) {
64 handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
65 const std::vector<T> & vals = *handle_;
66 for (size_t j = 0, n = vals.size(); j < n; ++j) {
67 if (!scanner.test(&vals[j])) continue;
68 if (printFullEventId_) {
69 const edm::EventAuxiliary &id = event_->eventAuxiliary();
70 printf(" : %9d : %4d : %9d : %3d", id.run(), id.luminosityBlock(), id.event(), j);
71 } else {
72 printf(" : %5d : %3d", iev, j);
73 }
74 scanner.print(&vals[j]);
75 std::cout << " :" << std::endl;
76 }
77 }
78 std::cout << std::endl;
79 }
80
81
82 TH1 * draw(const char *expr, const char *cut = "", TString drawopt = "", TH1 *hist = 0) {
83 // prep the machinery
84 helper::ScannerBase scanner(objType);
85 scanner.setIgnoreExceptions(ignoreExceptions_);
86 scanner.addExpression(expr);
87 if (strlen(cut)) scanner.setCut(cut);
88
89 // make histo, if needed
90 if (hist == 0) {
91 hist = new TH1F("htemp",
92 (strlen(cut) ? TString(expr)+"{"+cut+"}" : TString(expr)),
93 gEnv->GetValue("Hist.Binning.1D.x",100), 0, 0);
94 hist->SetBit(TH1::kCanRebin);
95 }
96
97 // fill histogram
98 for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
99 handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
100 const std::vector<T> & vals = *handle_;
101 for (size_t j = 0, n = vals.size(); j < n; ++j) {
102 scanner.fill1D(&vals[j], hist);
103 }
104 }
105
106 if (!TString(drawopt).Contains("goff",TString::kIgnoreCase)) hist->Draw(drawopt);
107 return hist;
108 }
109 TH1 * draw(const char *expr, int bins, double xlow, double xhigh, const char *cut = "", const char *drawopt = "") {
110 TH1 * htemp = new TH1F("htemp", (strlen(cut) ? TString(expr)+"{"+cut+"}" : TString(expr)), bins, xlow, xhigh);
111 return draw(expr,cut,drawopt,htemp);
112 }
113
114 TProfile * drawProf(TString xexpr, TString yexpr, const char *cut = "", TString drawopt = "", TProfile *hist = 0) {
115 // prep the machinery
116 helper::ScannerBase scanner(objType);
117 scanner.setIgnoreExceptions(ignoreExceptions_);
118 scanner.addExpression(xexpr);
119 scanner.addExpression(yexpr);
120 if (strlen(cut)) scanner.setCut(cut);
121
122 // make histo, if needed
123 if (hist == 0) {
124 hist = new TProfile("htemp",
125 (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr),
126 gEnv->GetValue("Hist.Binning.1D.x",100), 0., 0.);
127 hist->SetBit(TH1::kCanRebin);
128 }
129
130 // fill histogram
131 for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
132 handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
133 const std::vector<T> & vals = *handle_;
134 for (size_t j = 0, n = vals.size(); j < n; ++j) {
135 scanner.fillProf(&vals[j], hist);
136 }
137 }
138
139 if (!TString(drawopt).Contains("goff",TString::kIgnoreCase)) hist->Draw(drawopt);
140 return hist;
141 }
142 TProfile * drawProf(TString xexpr, int bins, double xlow, double xhigh, TString yexpr, const char *cut = "", const char *drawopt = "") {
143 TProfile * htemp = new TProfile("htemp", (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr), bins, xlow, xhigh);
144 return drawProf(xexpr,yexpr,cut,drawopt,htemp);
145 }
146
147 TH2 * draw2D(TString xexpr, TString yexpr, const char *cut = "", TString drawopt = "", TH2 *hist = 0) {
148 // prep the machinery
149 helper::ScannerBase scanner(objType);
150 scanner.setIgnoreExceptions(ignoreExceptions_);
151 scanner.addExpression(xexpr);
152 scanner.addExpression(yexpr);
153 if (strlen(cut)) scanner.setCut(cut);
154
155 // make histo, if needed
156 if (hist == 0) {
157 hist = new TH2F("htemp",
158 (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr),
159 gEnv->GetValue("Hist.Binning.2D.x",20), 0, 0,
160 gEnv->GetValue("Hist.Binning.2D.y",20), 0, 0);
161 hist->SetBit(TH1::kCanRebin);
162 }
163
164 // fill histogram
165 for (event_->toBegin(); !event_->atEnd(); ++(*event_)) {
166 handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
167 const std::vector<T> & vals = *handle_;
168 for (size_t j = 0, n = vals.size(); j < n; ++j) {
169 scanner.fill2D(&vals[j], hist);
170 }
171 }
172
173 if (!strlen(hist->GetXaxis()->GetTitle())) hist->GetXaxis()->SetTitle(xexpr);
174 if (!strlen(hist->GetYaxis()->GetTitle())) hist->GetYaxis()->SetTitle(yexpr);
175 if (!TString(drawopt).Contains("goff",TString::kIgnoreCase)) hist->Draw(drawopt);
176 return hist;
177 }
178 TH2 * draw2D(TString xexpr, int xbins, double xlow, double xhigh,
179 TString yexpr, int ybins, double ylow, double yhigh,
180 const char *cut = "", const char *drawopt = "") {
181 TH2 * htemp = new TH2F("htemp", (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr),
182 xbins, xlow, xhigh, ybins,ylow,yhigh);
183 return draw2D(xexpr,yexpr,cut,drawopt,htemp);
184 }
185
186 void setPrintFullEventId(bool printIt=true) { printFullEventId_ = printIt; }
187 void setExpressionSeparator(TString separator) { exprSep_ = separator; }
188 void setIgnoreExceptions(bool ignoreThem) { ignoreExceptions_ = ignoreThem; }
189 private:
190 fwlite::EventBase *event_;
191 std::string label_, instance_, process_;
192 bool printFullEventId_;
193 bool ignoreExceptions_;
194 TString exprSep_;
195 HandleT handle_;
196 helper::Parser helpme;
197 Reflex::Type objType;
198
199 };
200 }