ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/GPetrucc/interface/fwlite/Scanner.h
Revision: 1.4
Committed: Sun Jan 31 10:20:21 2010 UTC (15 years, 3 months ago) by gpetrucc
Content type: text/plain
Branch: MAIN
Changes since 1.3: +111 -3 lines
Log Message:
First implementation of draw

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