CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
fwlite::Scanner< Collection > Class Template Reference

fwlite::Scanner<C>, a way to inspect or plots elements of a collection C by using the StringParser. More...

#include <Scanner.h>

Public Types

typedef fwlite::Handle< Collection > HandleT
 The type of the Handle to read the Ts from the event. Needed to resolve its Type. More...
 

Public Member Functions

void addEventSelector (fwlite::EventSelector *selector)
 
void clearEventSelector ()
 
size_t count (const char *cut)
 
size_t countEvents ()
 
TH1 * draw (const char *expr, const char *cut, TString drawopt, TH1 *hist)
 
TH1 * draw (const char *expr, const char *cut="", TString drawopt="", const char *hname="htemp", const TH1 *htemplate=0)
 
TH1 * draw (const char *expr, int nbins, double xlow, double xhigh, const char *cut="", const char *drawopt="", const char *hname="htemp")
 
TH1 * draw (const char *expr, int nbins, double *xbins, const char *cut="", const char *drawopt="", const char *hname="htemp")
 
TH2 * draw2D (TString xexpr, TString yexpr, const char *cut, TString drawopt, TH2 *hist)
 Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it's usually (y,x)! More...
 
TH2 * draw2D (TString xexpr, TString yexpr, const char *cut="", TString drawopt="", const char *hname="htemp", TH2 *htemplate=0)
 
TH2 * draw2D (TString xexpr, int xbins, double xlow, double xhigh, TString yexpr, int ybins, double ylow, double yhigh, const char *cut="", const char *drawopt="", const char *hname="htemp")
 Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it's usually (y,x)! More...
 
TGraph * drawGraph (TString xexpr, TString yexpr, const char *cut, TString drawopt, TGraph *graph)
 
TGraph * drawGraph (TString xexpr, TString yexpr, const char *cut="", TString drawopt="AP", const char *gname="htemp")
 
TProfile * drawProf (TString xexpr, TString yexpr, const char *cut, TString drawopt, TProfile *hist)
 Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it's usually (y,x)! More...
 
TProfile * drawProf (TString xexpr, TString yexpr, const char *cut="", TString drawopt="", const char *hname="htemp", TProfile *htemplate=0)
 Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it's usually (y,x)! More...
 
TProfile * drawProf (TString xexpr, int bins, double xlow, double xhigh, TString yexpr, const char *cut="", const char *drawopt="", const char *hname="htemp")
 Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it's usually (y,x)! More...
 
TObjArray & eventSelectors ()
 
RooDataSet * fillDataSet (const char *realvars, const char *boolvars, const char *cut="", const char *name="data")
 
void scan (const char *exprs, const char *cut="", int nmax=-1)
 
 Scanner (fwlite::EventBase *ev, const char *label, const char *instance="", const char *process="")
 
bool selectEvent (const fwlite::EventBase &ev) const
 
void setExpressionSeparator (TString separator)
 
void setIgnoreExceptions (bool ignoreThem)
 
void setMaxEvents (int max)
 
void setMaxLinesToPrint (int lines)
 
void setPrintFullEventId (bool printIt=true)
 

Private Member Functions

TH1 * getSameH1 (const char *hname)
 
TH2 * getSameH2 (const char *hname)
 
TProfile * getSameProf (const char *hname)
 
void htempDelete ()
 
bool wantMore () const
 

Private Attributes

fwlite::EventBaseevent_
 
TObjArray eventSelectors_
 
TString exprSep_
 
HandleT handle_
 
bool ignoreExceptions_
 
std::string instance_
 
std::string label_
 
int maxEvents_
 
int maxLinesToPrint_
 
edm::TypeWithDict objType
 
bool printFullEventId_
 
std::string process_
 

Detailed Description

template<typename Collection>
class fwlite::Scanner< Collection >

fwlite::Scanner<C>, a way to inspect or plots elements of a collection C by using the StringParser.

fwlite::Scanner<C>, a way to inspect or plots elements of a collection C by using the StringParser.

The collection can be something as easy as std::vector<T>, but also some other fancy EDM collections like RefVector, RefToBaseVector and OwnVector (and probably PtrVector, but it was not tested)

If you're using something other than std::vector, you must provide the full typename, including all optional template parameters; e.g. you can't have C = edm::RefVector<reco::MuonCollection>, but you need C = edm::RefVector<vector<reco::Muon>,reco::Muon,edm::refhelper::FindUsingAdvance<vector<reco::Muon>,reco::Muon> > In order to figure out what is the correct full name for a collection in an event, open it in ROOT/FWLite, get the branch name including the trailing ".obj" (hint: Events->GetAlias("label")) usually works), and then do Events->GetBranch("xxx.obj")->GetClassName() to get something like edm::Wrapper<X>. then X is what you want to use to create the fwlite::Scanner. Don't use typedefs, they don't work.

Definition at line 44 of file Scanner.h.

Member Typedef Documentation

template<typename Collection >
typedef fwlite::Handle<Collection> fwlite::Scanner< Collection >::HandleT

The type of the Handle to read the Ts from the event. Needed to resolve its Type.

Definition at line 47 of file Scanner.h.

Constructor & Destructor Documentation

template<typename Collection >
fwlite::Scanner< Collection >::Scanner ( fwlite::EventBase ev,
const char *  label,
const char *  instance = "",
const char *  process = "" 
)
inline

Create a Scanner, passing a fwlite Event and the labels (just like you would in 'getByLabel')

Definition at line 50 of file Scanner.h.

References helper::Parser::elementType(), fwlite::Scanner< Collection >::objType, and edm::Wrapper< T >::typeInfo().

50  :
53  ignoreExceptions_(false),
54  exprSep_(":"),
55  maxEvents_(-1),
57  {
59  }
bool ignoreExceptions_
Definition: Scanner.h:582
static PFTauRenderPlugin instance
bool printFullEventId_
Definition: Scanner.h:581
TString exprSep_
Definition: Scanner.h:583
bool isRealData() const
Definition: EventBase.h:64
std::string label_
Definition: Scanner.h:580
edm::TypeWithDict objType
Definition: Scanner.h:585
fwlite::EventBase * event_
Definition: Scanner.h:579
std::string instance_
Definition: Scanner.h:580
int maxLinesToPrint_
Definition: Scanner.h:591
static edm::TypeWithDict elementType(const edm::TypeWithDict &wrapperType)
Perform the type deduction form edm::Wrapper<C> to C::value_type and resolves typedefs.
static std::type_info const & typeInfo()
Definition: Wrapper.h:37

Member Function Documentation

template<typename Collection >
void fwlite::Scanner< Collection >::addEventSelector ( fwlite::EventSelector selector)
inline

Definition at line 567 of file Scanner.h.

References fwlite::Scanner< Collection >::eventSelectors_.

567 { eventSelectors_.Add(selector); }
TObjArray eventSelectors_
Definition: Scanner.h:587
template<typename Collection >
void fwlite::Scanner< Collection >::clearEventSelector ( )
inline

Definition at line 568 of file Scanner.h.

References fwlite::Scanner< Collection >::eventSelectors_.

568 { eventSelectors_.Clear(); }
TObjArray eventSelectors_
Definition: Scanner.h:587
template<typename Collection >
size_t fwlite::Scanner< Collection >::count ( const char *  cut)
inline

Count the number of entries that pass a given cut. See setMaxEvents() to specify how many events to loop on when counting. Events can be further selected by using addEventSelector().

Definition at line 145 of file Scanner.h.

References fwlite::EventBase::atEnd(), fwlite::Scanner< Collection >::event_, fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::handle_, fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, fwlite::Scanner< Collection >::maxEvents_, gen::n, fwlite::Scanner< Collection >::objType, fwlite::Scanner< Collection >::process_, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), helper::ScannerBase::test(), fwlite::EventBase::toBegin(), and create_public_pileup_plots::vals.

145  {
146  helper::ScannerBase scanner(objType);
147  scanner.setIgnoreExceptions(ignoreExceptions_);
148 
149  scanner.setCut(cut);
150 
151  size_t npass = 0;
152  int iev = 0;
153  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
154  if (maxEvents_ > -1 && iev > maxEvents_) break;
155  if (!selectEvent(*event_)) continue;
156  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
157  const Collection & vals = *handle_;
158  for (size_t j = 0, n = vals.size(); j < n; ++j) {
159  if (scanner.test(&vals[j])) npass++;
160  }
161  }
162  return npass;
163  }
bool ignoreExceptions_
Definition: Scanner.h:582
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
std::string instance_
Definition: Scanner.h:580
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
template<typename Collection >
size_t fwlite::Scanner< Collection >::countEvents ( )
inline

Count the number of events, taking into account setMaxEvents() and the event selectors

Definition at line 166 of file Scanner.h.

References fwlite::EventBase::atEnd(), fwlite::Scanner< Collection >::event_, fwlite::Scanner< Collection >::maxEvents_, fwlite::Scanner< Collection >::selectEvent(), and fwlite::EventBase::toBegin().

166  {
167  size_t npass = 0;
168  int iev = 0;
169  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
170  if (maxEvents_ > -1 && iev > maxEvents_) break;
171  if (selectEvent(*event_)) npass++;
172  }
173  return npass;
174  }
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
template<typename Collection >
TH1* fwlite::Scanner< Collection >::draw ( const char *  expr,
const char *  cut,
TString  drawopt,
TH1 *  hist 
)
inline

Plot the expression expr for events passing 'cut, into histogram hist. hist is not reset before filling it, so it will add to the existing content.

If "NORM" is specified in the draw options, the output histogram is normalized If "GOFF" is specified in the draw options, the output histogram is not drawn

See setMaxEvents() to specify how many events to loop on when plotting. Events can be further selected by using addEventSelector().

Definition at line 185 of file Scanner.h.

References helper::ScannerBase::addExpression(), fwlite::EventBase::atEnd(), MessageLogger_cfi::cerr, fwlite::Scanner< Collection >::event_, helper::ScannerBase::fill1D(), fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::handle_, create_public_lumi_plots::hist, fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, fwlite::Scanner< Collection >::maxEvents_, gen::n, fwlite::Scanner< Collection >::objType, fwlite::Scanner< Collection >::process_, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), fwlite::EventBase::toBegin(), and create_public_pileup_plots::vals.

Referenced by fwlite::Scanner< Collection >::draw().

185  {
186  // prep the machinery
187  helper::ScannerBase scanner(objType);
188  scanner.setIgnoreExceptions(ignoreExceptions_);
189  if (!scanner.addExpression(expr)) return 0;
190  if (strlen(cut)) scanner.setCut(cut);
191 
192  // check histo
193  if (hist == 0) {
194  std::cerr << "Method draw(expr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other draw methods instead." << std::endl;
195  return 0;
196  }
197 
198  // fill histogram
199  int iev = 0;
200  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
201  if (maxEvents_ > -1 && iev > maxEvents_) break;
202  if (!selectEvent(*event_)) continue;
203  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
204  const Collection & vals = *handle_;
205  for (size_t j = 0, n = vals.size(); j < n; ++j) {
206  scanner.fill1D(&vals[j], hist);
207  }
208  }
209 
210  if (drawopt.Contains("NORM",TString::kIgnoreCase) && (hist->Integral() != 0)) {
211  hist->Sumw2();
212  hist->Scale(1.0/hist->Integral());
213  // remove the "NORM" because THistPainter doesn't understand it
214  drawopt(TRegexp("[Nn][Oo][Rr][Mm]")) = "";
215  }
216 
217  if (!drawopt.Contains("GOFF",TString::kIgnoreCase)) hist->Draw(drawopt);
218  return hist;
219  }
bool ignoreExceptions_
Definition: Scanner.h:582
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
std::string instance_
Definition: Scanner.h:580
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
template<typename Collection >
TH1* fwlite::Scanner< Collection >::draw ( const char *  expr,
const char *  cut = "",
TString  drawopt = "",
const char *  hname = "htemp",
const TH1 *  htemplate = 0 
)
inline

Plot the expression 'expr' for events passing 'cut', in a histogram named 'hname'

  • If htemplate is provided it will be cloned,
  • otherwise, if "SAME" is among the draw options, it will clone "htemp" (if it's available)
  • otherwise an automatically binned histogram will be used. in the last case, gEnv->GetValue("Hist.Binning.1D.x", 100) is used for the number of bins See draw(const char *expr, const char *cut, TString drawopt, TH1 *histo) for further documentation

Definition at line 227 of file Scanner.h.

References TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::draw(), fwlite::Scanner< Collection >::getSameH1(), create_public_lumi_plots::hist, and fwlite::Scanner< Collection >::htempDelete().

227  {
228  TH1 *hist = 0;
229  if (htemplate != 0) {
230  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname,htemplate->GetName()) != 0)) htempDelete();
231  hist = (TH1*) hist->Clone(hname);
232  } else if (drawopt.Contains("SAME",TString::kIgnoreCase)) {
233  hist = getSameH1(hname);
234  }
235 
236  // if in the end we found no way to make "hist"
237  if (hist == 0) {
238  if (strcmp(hname, "htemp") == 0) htempDelete();
239  hist = new TH1F(hname, "", gEnv->GetValue("Hist.Binning.1D.x",100), 0, 0);
240  hist->SetCanExtend(TH1::kAllAxes);
241  }
242  hist->SetTitle((strlen(cut) ? TString(expr)+"{"+cut+"}" : TString(expr)));
243  hist->GetXaxis()->SetTitle(expr);
244  return draw(expr, cut, drawopt, hist);
245  }
TH1 * getSameH1(const char *hname)
Definition: Scanner.h:612
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
Definition: Scanner.h:185
void htempDelete()
Definition: Scanner.h:603
template<typename Collection >
TH1* fwlite::Scanner< Collection >::draw ( const char *  expr,
int  nbins,
double  xlow,
double  xhigh,
const char *  cut = "",
const char *  drawopt = "",
const char *  hname = "htemp" 
)
inline

Make a histogram named hname with nbins from xlow to xhigh, and then call draw(). If "SAME" is passed in the draw options, complain and ignore the binning.

Definition at line 250 of file Scanner.h.

References MessageLogger_cfi::cerr, TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::draw(), fwlite::Scanner< Collection >::getSameH1(), and fwlite::Scanner< Collection >::htempDelete().

250  {
251  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
252  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
253  TH1 *hsame = getSameH1(hname);
254  return draw(expr, cut, drawopt, hsame);
255  }
256  if (strcmp(hname, "htemp") == 0) htempDelete();
257  TH1 * htemp = new TH1F(hname, expr, nbins, xlow, xhigh);
258  if (strlen(cut)) htemp->SetTitle(TString(expr)+"{"+cut+"}");
259  htemp->GetXaxis()->SetTitle(expr);
260  return draw(expr,cut,drawopt,htemp);
261  }
TH1 * getSameH1(const char *hname)
Definition: Scanner.h:612
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
Definition: Scanner.h:185
void htempDelete()
Definition: Scanner.h:603
template<typename Collection >
TH1* fwlite::Scanner< Collection >::draw ( const char *  expr,
int  nbins,
double *  xbins,
const char *  cut = "",
const char *  drawopt = "",
const char *  hname = "htemp" 
)
inline

Make a histogram named hname with nbins with boundaries xbins, and then call draw(). If "SAME" is passed in the draw options, complain and ignore the binning.

Definition at line 264 of file Scanner.h.

References MessageLogger_cfi::cerr, TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::draw(), fwlite::Scanner< Collection >::getSameH1(), and fwlite::Scanner< Collection >::htempDelete().

264  {
265  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
266  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
267  TH1 *hsame = getSameH1(hname);
268  return draw(expr, cut, drawopt, hsame);
269  }
270  if (strcmp(hname, "htemp") == 0) htempDelete();
271  TH1 * htemp = new TH1F(hname, expr, nbins, xbins);
272  if (strlen(cut)) htemp->SetTitle(TString(expr)+"{"+cut+"}");
273  htemp->GetXaxis()->SetTitle(expr);
274  return draw(expr,cut,drawopt,htemp);
275  }
const double xbins[]
TH1 * getSameH1(const char *hname)
Definition: Scanner.h:612
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
Definition: Scanner.h:185
void htempDelete()
Definition: Scanner.h:603
template<typename Collection >
TH2* fwlite::Scanner< Collection >::draw2D ( TString  xexpr,
TString  yexpr,
const char *  cut,
TString  drawopt,
TH2 *  hist 
)
inline

Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it's usually (y,x)!

Definition at line 344 of file Scanner.h.

References helper::ScannerBase::addExpression(), fwlite::EventBase::atEnd(), MessageLogger_cfi::cerr, fwlite::Scanner< Collection >::event_, helper::ScannerBase::fill2D(), fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::handle_, create_public_lumi_plots::hist, fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, fwlite::Scanner< Collection >::maxEvents_, gen::n, fwlite::Scanner< Collection >::objType, fwlite::Scanner< Collection >::process_, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), fwlite::EventBase::toBegin(), and create_public_pileup_plots::vals.

Referenced by fwlite::Scanner< Collection >::draw2D().

344  {
345  // prep the machinery
346  helper::ScannerBase scanner(objType);
347  scanner.setIgnoreExceptions(ignoreExceptions_);
348  if (!scanner.addExpression((const char *)xexpr)) return 0;
349  if (!scanner.addExpression((const char *)yexpr)) return 0;
350  if (strlen(cut)) scanner.setCut(cut);
351 
352  // check histo
353  if (hist == 0) {
354  std::cerr << "Method draw2D(xexpr, yexpr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other draw methods instead." << std::endl;
355  return 0;
356  }
357 
358  // fill histogram
359  int iev = 0;
360  for (event_->toBegin(), iev = 0; !event_->atEnd(); ++(*event_), ++iev) {
361  if (maxEvents_ > -1 && iev > maxEvents_) break;
362  if (!selectEvent(*event_)) continue;
363  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
364  const Collection & vals = *handle_;
365  for (size_t j = 0, n = vals.size(); j < n; ++j) {
366  scanner.fill2D(&vals[j], hist);
367  }
368  }
369 
370  if (!strlen(hist->GetTitle())) hist->SetTitle((strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
371  if (!strlen(hist->GetXaxis()->GetTitle())) hist->GetXaxis()->SetTitle(xexpr);
372  if (!strlen(hist->GetYaxis()->GetTitle())) hist->GetYaxis()->SetTitle(yexpr);
373  if (!TString(drawopt).Contains("GOFF",TString::kIgnoreCase)) hist->Draw(drawopt);
374  return hist;
375  }
bool ignoreExceptions_
Definition: Scanner.h:582
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
std::string instance_
Definition: Scanner.h:580
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
template<typename Collection >
TH2* fwlite::Scanner< Collection >::draw2D ( TString  xexpr,
TString  yexpr,
const char *  cut = "",
TString  drawopt = "",
const char *  hname = "htemp",
TH2 *  htemplate = 0 
)
inline

Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it's usually (y,x)! Note that automatical binning for TH2s is more expensive, as it requires to loop on the events twice!

Definition at line 378 of file Scanner.h.

References helper::ScannerBase::addExpression(), fwlite::EventBase::atEnd(), TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::draw2D(), helper::ScannerBase::eval(), fwlite::Scanner< Collection >::event_, fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::getSameH2(), fwlite::Scanner< Collection >::handle_, create_public_lumi_plots::hist, fwlite::Scanner< Collection >::htempDelete(), fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, fwlite::Scanner< Collection >::maxEvents_, gen::n, fwlite::Scanner< Collection >::objType, fwlite::Scanner< Collection >::process_, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), helper::ScannerBase::test(), fwlite::EventBase::toBegin(), create_public_pileup_plots::vals, TrackerOfflineValidation_Dqm_cff::xmax, TrackerOfflineValidation_Dqm_cff::xmin, Phase2TrackerMonitorDigi_cff::ymax, and Phase2TrackerMonitorDigi_cff::ymin.

378  {
379  TH2 *hist = 0;
380  if (htemplate != 0) {
381  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname,htemplate->GetName()) != 0)) htempDelete();
382  hist = (TH2*) hist->Clone(hname);
383  } else if (drawopt.Contains("SAME",TString::kIgnoreCase)) {
384  hist = getSameH2(hname);
385  }
386 
387  // if in the end we found no way to make "hist"
388  if (hist == 0) {
389  // prep the machinery
390  helper::ScannerBase scanner(objType);
391  scanner.setIgnoreExceptions(ignoreExceptions_);
392  if (!scanner.addExpression((const char *)xexpr)) return 0;
393  if (!scanner.addExpression((const char *)yexpr)) return 0;
394  if (strlen(cut)) scanner.setCut(cut);
395 
396  if (strcmp(hname, "htemp") == 0) htempDelete();
397  // ok this is much more a hack than for the 1D case
398  double xmin = 0, xmax = -1, ymin = 0, ymax = -1; int iev;
399  for (event_->toBegin(), iev = 0; !event_->atEnd(); ++(*event_), ++iev) {
400  if (maxEvents_ > -1 && iev > maxEvents_) break;
401  if (!selectEvent(*event_)) continue;
402  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
403  const Collection & vals = *handle_;
404  for (size_t j = 0, n = vals.size(); j < n; ++j) {
405  if (!scanner.test(&vals[j])) continue;
406  double x = scanner.eval(&vals[j],0);
407  double y = scanner.eval(&vals[j],1);
408  if ((xmax == -1) || (x >= xmax)) xmax = x;
409  if ((xmin == 0) || (x <= xmin)) xmin = x;
410  if ((ymax == -1) || (y >= ymax)) ymax = y;
411  if ((ymin == 0) || (y <= ymin)) ymin = y;
412  }
413  }
414  hist = new TH2F(hname, "",
415  gEnv->GetValue("Hist.Binning.2D.x",20), xmin, xmax,
416  gEnv->GetValue("Hist.Binning.2D.y",20), ymin, ymax);
417  }
418  return draw2D(xexpr, yexpr, cut, drawopt, hist);
419  }
bool ignoreExceptions_
Definition: Scanner.h:582
TH2 * getSameH2(const char *hname)
Definition: Scanner.h:628
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
std::string instance_
Definition: Scanner.h:580
TH2 * draw2D(TString xexpr, TString yexpr, const char *cut, TString drawopt, TH2 *hist)
Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it&#39;s usually (y...
Definition: Scanner.h:344
void htempDelete()
Definition: Scanner.h:603
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
template<typename Collection >
TH2* fwlite::Scanner< Collection >::draw2D ( TString  xexpr,
int  xbins,
double  xlow,
double  xhigh,
TString  yexpr,
int  ybins,
double  ylow,
double  yhigh,
const char *  cut = "",
const char *  drawopt = "",
const char *  hname = "htemp" 
)
inline

Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it's usually (y,x)!

Definition at line 422 of file Scanner.h.

References MessageLogger_cfi::cerr, TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::draw2D(), fwlite::Scanner< Collection >::getSameH2(), and fwlite::Scanner< Collection >::htempDelete().

424  {
425  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
426  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
427  TH2 *hsame = getSameH2(hname);
428  return draw2D(xexpr, yexpr, cut, drawopt, hsame);
429  }
430  if (strcmp(hname, "htemp") == 0) htempDelete();
431  TH2 * htemp = new TH2F("htemp", "", xbins, xlow, xhigh, ybins,ylow,yhigh);
432  return draw2D(xexpr,yexpr,cut,drawopt,htemp);
433  }
const double xbins[]
TH2 * getSameH2(const char *hname)
Definition: Scanner.h:628
TH2 * draw2D(TString xexpr, TString yexpr, const char *cut, TString drawopt, TH2 *hist)
Just like draw() except that it uses TH2. Note that the order is (x,y) while in ROOT it&#39;s usually (y...
Definition: Scanner.h:344
void htempDelete()
Definition: Scanner.h:603
template<typename Collection >
TGraph* fwlite::Scanner< Collection >::drawGraph ( TString  xexpr,
TString  yexpr,
const char *  cut,
TString  drawopt,
TGraph *  graph 
)
inline

Draw a scatter plot of x vs y for events passing the cut.

Definition at line 436 of file Scanner.h.

References helper::ScannerBase::addExpression(), fwlite::EventBase::atEnd(), fwlite::Scanner< Collection >::event_, helper::ScannerBase::fillGraph(), fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::handle_, fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, fwlite::Scanner< Collection >::maxEvents_, gen::n, fwlite::Scanner< Collection >::objType, fwlite::Scanner< Collection >::process_, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), fwlite::EventBase::toBegin(), and create_public_pileup_plots::vals.

Referenced by fwlite::Scanner< Collection >::drawGraph().

436  {
437  // prep the machinery
438  helper::ScannerBase scanner(objType);
439  scanner.setIgnoreExceptions(ignoreExceptions_);
440  if (!scanner.addExpression((const char *)xexpr)) return 0;
441  if (!scanner.addExpression((const char *)yexpr)) return 0;
442  if (strlen(cut)) scanner.setCut(cut);
443 
444  // make graph, if needed
445  if (graph == 0) {
446  graph = new TGraph();
447  graph->SetNameTitle("htemp", (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
448  }
449 
450  // fill graph
451  int iev = 0;
452  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
453  if (maxEvents_ > -1 && iev > maxEvents_) break;
454  if (!selectEvent(*event_)) continue;
455  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
456  const Collection & vals = *handle_;
457  for (size_t j = 0, n = vals.size(); j < n; ++j) {
458  scanner.fillGraph(&vals[j], graph);
459  }
460  }
461 
462  if (!strlen(graph->GetTitle())) graph->SetTitle((strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
463  if (!strlen(graph->GetXaxis()->GetTitle())) graph->GetXaxis()->SetTitle(xexpr);
464  if (!strlen(graph->GetYaxis()->GetTitle())) graph->GetYaxis()->SetTitle(yexpr);
465  if (!TString(drawopt).Contains("GOFF",TString::kIgnoreCase)) graph->Draw(drawopt);
466  return graph;
467  }
bool ignoreExceptions_
Definition: Scanner.h:582
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
std::string instance_
Definition: Scanner.h:580
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
template<typename Collection >
TGraph* fwlite::Scanner< Collection >::drawGraph ( TString  xexpr,
TString  yexpr,
const char *  cut = "",
TString  drawopt = "AP",
const char *  gname = "htemp" 
)
inline

Draw a scatter plot of x vs y for events passing the cut.

Definition at line 470 of file Scanner.h.

References TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::drawGraph(), and fwlite::Scanner< Collection >::htempDelete().

470  {
471  if (strcmp(gname, "htemp") == 0) htempDelete();
472  TGraph *graph = new TGraph();
473  graph->SetNameTitle(gname, (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
474  return drawGraph(xexpr,yexpr,cut,drawopt,graph);
475  }
TGraph * drawGraph(TString xexpr, TString yexpr, const char *cut, TString drawopt, TGraph *graph)
Definition: Scanner.h:436
void htempDelete()
Definition: Scanner.h:603
template<typename Collection >
TProfile* fwlite::Scanner< Collection >::drawProf ( TString  xexpr,
TString  yexpr,
const char *  cut,
TString  drawopt,
TProfile *  hist 
)
inline

Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it's usually (y,x)!

Definition at line 279 of file Scanner.h.

References helper::ScannerBase::addExpression(), fwlite::EventBase::atEnd(), MessageLogger_cfi::cerr, fwlite::Scanner< Collection >::event_, helper::ScannerBase::fillProf(), fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::handle_, create_public_lumi_plots::hist, fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, fwlite::Scanner< Collection >::maxEvents_, gen::n, fwlite::Scanner< Collection >::objType, fwlite::Scanner< Collection >::process_, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), fwlite::EventBase::toBegin(), and create_public_pileup_plots::vals.

Referenced by fwlite::Scanner< Collection >::drawProf().

279  {
280  // prep the machinery
281  helper::ScannerBase scanner(objType);
282  scanner.setIgnoreExceptions(ignoreExceptions_);
283  if (!scanner.addExpression(xexpr.Data())) return 0;
284  if (!scanner.addExpression(yexpr.Data())) return 0;
285  if (strlen(cut)) scanner.setCut(cut);
286 
287  // check histo
288  if (hist == 0) {
289  std::cerr << "Method drawProf(xexpr, yexpr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other draw methods instead." << std::endl;
290  return 0;
291  }
292 
293  // fill histogram
294  int iev = 0;
295  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
296  if (maxEvents_ > -1 && iev > maxEvents_) break;
297  if (!selectEvent(*event_)) continue;
298  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
299  const Collection & vals = *handle_;
300  for (size_t j = 0, n = vals.size(); j < n; ++j) {
301  scanner.fillProf(&vals[j], hist);
302  }
303  }
304 
305  if (!strlen(hist->GetTitle())) hist->SetTitle((strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
306  if (!strlen(hist->GetXaxis()->GetTitle())) hist->GetXaxis()->SetTitle(xexpr);
307  if (!strlen(hist->GetYaxis()->GetTitle())) hist->GetYaxis()->SetTitle(yexpr);
308  if (!TString(drawopt).Contains("GOFF",TString::kIgnoreCase)) hist->Draw(drawopt);
309  return hist;
310  }
bool ignoreExceptions_
Definition: Scanner.h:582
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
std::string instance_
Definition: Scanner.h:580
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
template<typename Collection >
TProfile* fwlite::Scanner< Collection >::drawProf ( TString  xexpr,
TString  yexpr,
const char *  cut = "",
TString  drawopt = "",
const char *  hname = "htemp",
TProfile *  htemplate = 0 
)
inline

Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it's usually (y,x)!

Definition at line 312 of file Scanner.h.

References TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::drawProf(), fwlite::Scanner< Collection >::getSameProf(), create_public_lumi_plots::hist, and fwlite::Scanner< Collection >::htempDelete().

312  {
313  TProfile *hist = 0;
314  if (htemplate != 0) {
315  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname,htemplate->GetName() )!= 0)) htempDelete();
316  hist = (TProfile*) hist->Clone(hname);
317  } else if (drawopt.Contains("SAME",TString::kIgnoreCase)) {
318  hist = getSameProf(hname);
319  }
320 
321  // if in the end we found no way to make "hist"
322  if (hist == 0) {
323  if (strcmp(hname, "htemp") == 0) htempDelete();
324  hist = new TProfile(hname, "", gEnv->GetValue("Hist.Binning.1D.x",100), 0., 0.);
325  hist->SetCanExtend(TH1::kAllAxes);
326  }
327  return drawProf(xexpr, yexpr, cut, drawopt, hist);
328  }
TProfile * drawProf(TString xexpr, TString yexpr, const char *cut, TString drawopt, TProfile *hist)
Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it&#39;s usuall...
Definition: Scanner.h:279
TProfile * getSameProf(const char *hname)
Definition: Scanner.h:644
void htempDelete()
Definition: Scanner.h:603
template<typename Collection >
TProfile* fwlite::Scanner< Collection >::drawProf ( TString  xexpr,
int  bins,
double  xlow,
double  xhigh,
TString  yexpr,
const char *  cut = "",
const char *  drawopt = "",
const char *  hname = "htemp" 
)
inline

Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it's usually (y,x)!

Definition at line 331 of file Scanner.h.

References MessageLogger_cfi::cerr, TkAlMuonSelectors_cfi::cut, fwlite::Scanner< Collection >::drawProf(), fwlite::Scanner< Collection >::getSameProf(), and fwlite::Scanner< Collection >::htempDelete().

331  {
332  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
333  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
334  TProfile *hsame = getSameProf(hname);
335  return drawProf(xexpr, yexpr, cut, drawopt, hsame);
336  }
337  if (strcmp(hname, "htemp") == 0) htempDelete();
338  TProfile * htemp = new TProfile(hname, "", bins, xlow, xhigh);
339  return drawProf(xexpr,yexpr,cut,drawopt,htemp);
340  }
TProfile * drawProf(TString xexpr, TString yexpr, const char *cut, TString drawopt, TProfile *hist)
Just like draw() except that it uses TProfile. Note that the order is (x,y) while in ROOT it&#39;s usuall...
Definition: Scanner.h:279
TProfile * getSameProf(const char *hname)
Definition: Scanner.h:644
void htempDelete()
Definition: Scanner.h:603
template<typename Collection >
TObjArray& fwlite::Scanner< Collection >::eventSelectors ( )
inline

Definition at line 569 of file Scanner.h.

References fwlite::Scanner< Collection >::eventSelectors_.

569 { return eventSelectors_; }
TObjArray eventSelectors_
Definition: Scanner.h:587
template<typename Collection >
RooDataSet* fwlite::Scanner< Collection >::fillDataSet ( const char *  realvars,
const char *  boolvars,
const char *  cut = "",
const char *  name = "data" 
)
inline

Fill a RooDataSet.

  • Real variables are defined just like in the scan() command; a list separated by ":" (see also setExpressionSeparator());
  • Boolean variables are defined just like cuts, and are created as RooCategory with two states: pass(1) and fail(0). For each variable, the name is taken from the expression itself, or can be specified manuall by using the notation "@name=expr" Note: the dataset contains one entry per item, irrespectively of how entries are distributed among events.

Definition at line 485 of file Scanner.h.

References helper::ScannerBase::addExpression(), helper::ScannerBase::addExtraCut(), fwlite::EventBase::atEnd(), eostools::cat(), MessageLogger_cfi::cerr, TkAlMuonSelectors_cfi::cut, helper::ScannerBase::eval(), fwlite::Scanner< Collection >::event_, fwlite::Scanner< Collection >::exprSep_, fwlite::Handle< T >::failedToGet(), fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::handle_, mps_fire::i, fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, fwlite::Scanner< Collection >::maxEvents_, gen::n, dataset::name, fwlite::Scanner< Collection >::objType, fwlite::Scanner< Collection >::process_, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), harvestTrackValidationPlots::str, AlCaHLTBitMon_QueryRunRegistry::string, helper::ScannerBase::test(), fwlite::EventBase::toBegin(), create_public_pileup_plots::vals, and JetChargeProducer_cfi::var.

485  {
486  helper::ScannerBase scanner(objType);
487  scanner.setIgnoreExceptions(ignoreExceptions_);
488 
489  RooArgList vars;
490  TObjArray *exprArray = TString(realvars).Tokenize(exprSep_);
491  TObjArray *catArray = TString(boolvars).Tokenize(exprSep_);
492  int nreals = exprArray->GetEntries();
493  int nbools = catArray->GetEntries();
494  for (int i = 0; i < nreals; ++i) {
495  TString str = ((TObjString *)(*exprArray)[i])->GetString();
496  std::string lb = str.Data();
497  std::string ex = str.Data();
498  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
499  lb = lb.substr(1,ex.find('=')-1);
500  ex = ex.substr(ex.find('=')+1);
501  }
502  if (!scanner.addExpression(ex.c_str())) {
503  std::cerr << "Filed to define real variable '" << lb << "', expr = '" << ex << "'" << std::endl;
504  return 0;
505  }
506  // FIXME: I have to leave it dangling on the HEAP otherwise ROOT segfaults...
507  RooRealVar *var = new RooRealVar(lb.c_str(),lb.c_str(), 0.0);
508  vars.add(*var);
509  }
510  for (int i = 0; i < nbools; ++i) {
511  TString str = ((TObjString *)(*catArray)[i])->GetString();
512  std::string lb = str.Data();
513  std::string ex = str.Data();
514  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
515  lb = lb.substr(1,ex.find('=')-1);
516  ex = ex.substr(ex.find('=')+1);
517  }
518  if (!scanner.addExtraCut(ex.c_str())) {
519  std::cerr << "Filed to define bool variable '" << lb << "', cut = '" << ex << "'" << std::endl;
520  return 0;
521  }
522  RooCategory *cat = new RooCategory(lb.c_str(), lb.c_str());
523  cat->defineType("fail",0);
524  cat->defineType("pass",1);
525  vars.add(*cat);
526  }
527 
528  RooDataSet *ds = new RooDataSet(name, name, vars);
529 
530  if (strlen(cut)) scanner.setCut(cut);
531  int iev = 0;
532  for (event_->toBegin(); !event_->atEnd(); ++iev, ++(*event_)) {
533  if (maxEvents_ > -1 && iev > maxEvents_) break;
534  if (!selectEvent(*event_)) continue;
535  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
536  if (handle_.failedToGet()) {
537  if (ignoreExceptions_) continue;
538  }
539  const Collection & vals = *handle_;
540  for (size_t j = 0, n = vals.size(); j < n; ++j) {
541  if (!scanner.test(&vals[j])) continue;
542  for (int i = 0; i < nreals; ++i) {
543  RooRealVar *var = (RooRealVar *)vars.at(i);
544  var->setVal(scanner.eval(&vals[j], i));
545  }
546  for (int i = 0; i < nbools; ++i) {
547  RooCategory *cat = (RooCategory*) vars.at(i+nreals);
548  cat->setIndex(int(scanner.test(&vals[j], i+1))); // 0 is the event selection cut
549  }
550  ds->add(vars);
551  }
552  }
553 
554  delete exprArray;
555  delete catArray;
556 
557  return ds;
558  }
bool ignoreExceptions_
Definition: Scanner.h:582
TString exprSep_
Definition: Scanner.h:583
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
def cat(path)
Definition: eostools.py:400
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
bool failedToGet() const
Returns true only if Handle was used in a &#39;get&#39; call and the data could not be found.
Definition: Handle.h:64
std::string instance_
Definition: Scanner.h:580
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
template<typename Collection >
TH1* fwlite::Scanner< Collection >::getSameH1 ( const char *  hname)
inlineprivate

Get whatever histogram makes sense for a plot passing "SAME" in drawOpt, and call it hname Currently it won't work if the histogram of which we want to be "SAME" is not called "htemp"

Definition at line 612 of file Scanner.h.

References MessageLogger_cfi::cerr, fftjetcommon_cfi::Class, and create_public_lumi_plots::hist.

Referenced by fwlite::Scanner< Collection >::draw().

612  {
613  if (gDirectory && gDirectory->Get("htemp") != 0 &&
614  gDirectory->Get("htemp")->IsA()->InheritsFrom(TH1::Class())) {
615  TH1 *hist = (TH1*) ((TH1*) gDirectory->Get("htemp"))->Clone(hname);
616  hist->Reset();
617  hist->SetLineColor(kBlack);
618  hist->SetMarkerColor(kBlack);
619  return hist;
620  } else {
621  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
622  return 0;
623  }
624  }
template<typename Collection >
TH2* fwlite::Scanner< Collection >::getSameH2 ( const char *  hname)
inlineprivate

Get whatever histogram makes sense for a plot passing "SAME" in drawOpt, and call it hname Currently it won't work if the histogram of which we want to be "SAME" is not called "htemp"

Definition at line 628 of file Scanner.h.

References MessageLogger_cfi::cerr, fftjetcommon_cfi::Class, and create_public_lumi_plots::hist.

Referenced by fwlite::Scanner< Collection >::draw2D().

628  {
629  if (gDirectory && gDirectory->Get("htemp") != 0 &&
630  gDirectory->Get("htemp")->IsA()->InheritsFrom(TH2::Class())) {
631  TH2 *hist = (TH2*) ((TH2*) gDirectory->Get("htemp"))->Clone(hname);
632  hist->Reset();
633  hist->SetLineColor(kBlack);
634  hist->SetMarkerColor(kBlack);
635  return hist;
636  } else {
637  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
638  return 0;
639  }
640  }
template<typename Collection >
TProfile* fwlite::Scanner< Collection >::getSameProf ( const char *  hname)
inlineprivate

Get whatever histogram makes sense for a plot passing "SAME" in drawOpt, and call it hname Currently it won't work if the histogram of which we want to be "SAME" is not called "htemp"

Definition at line 644 of file Scanner.h.

References MessageLogger_cfi::cerr, fftjetcommon_cfi::Class, and create_public_lumi_plots::hist.

Referenced by fwlite::Scanner< Collection >::drawProf().

644  {
645  if (gDirectory && gDirectory->Get("htemp") != 0 &&
646  gDirectory->Get("htemp")->IsA()->InheritsFrom(TProfile::Class())) {
647  TProfile *hist = (TProfile*) ((TProfile*) gDirectory->Get("htemp"))->Clone(hname);
648  hist->Reset();
649  hist->SetLineColor(kBlack);
650  hist->SetMarkerColor(kBlack);
651  return hist;
652  } else {
653  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
654  return 0;
655  }
656  }
template<typename Collection >
void fwlite::Scanner< Collection >::htempDelete ( )
inlineprivate

Definition at line 603 of file Scanner.h.

References MuonAssociatorByHits_cfi::obj.

Referenced by fwlite::Scanner< Collection >::draw(), fwlite::Scanner< Collection >::draw2D(), fwlite::Scanner< Collection >::drawGraph(), and fwlite::Scanner< Collection >::drawProf().

603  {
604  if (gDirectory) {
605  TObject *obj = gDirectory->Get("htemp");
606  if (obj) obj->Delete();
607  }
608  }
template<typename Collection >
void fwlite::Scanner< Collection >::scan ( const char *  exprs,
const char *  cut = "",
int  nmax = -1 
)
inline

Scan the first nmax entries of the event and print out the values of some expressions.

The cut is applied to the individual entries. To set Event-wide cuts, use addEventSelector().

The different expressions are separated by ":", unless changed using setExpressionSeparator. The title of each column is the text of the expression, unless one specifies it differently by using the notation "@label=expression"

Each row is prefixed by the event id (Run/LS/Event on Data, entry number within the file for MC) and by the index of the object within the collection. The behaviour can be changed through the setPrintFullEventId() method.

The printing will pause by default every 50 lines (see setMaxLinesToPrint() to change this) Scanning will stop after nmax events.

Definition at line 77 of file Scanner.h.

References helper::ScannerBase::addExpression(), fwlite::EventBase::atEnd(), gather_cfg::cout, TkAlMuonSelectors_cfi::cut, event(), fwlite::Scanner< Collection >::event_, edm::EventBase::eventAuxiliary(), fwlite::Scanner< Collection >::exprSep_, fwlite::Handle< T >::failedToGet(), fwlite::Handle< T >::getByLabel(), fwlite::Scanner< Collection >::handle_, mps_fire::i, fwlite::Scanner< Collection >::ignoreExceptions_, fwlite::Scanner< Collection >::instance_, fwlite::Scanner< Collection >::label_, geometryCSVtoXML::line, fwlite::Scanner< Collection >::maxLinesToPrint_, gen::n, fwlite::Scanner< Collection >::objType, helper::ScannerBase::print(), fwlite::Scanner< Collection >::printFullEventId_, fwlite::Scanner< Collection >::process_, findQualityFiles::run, fwlite::Scanner< Collection >::selectEvent(), helper::ScannerBase::setCut(), helper::ScannerBase::setIgnoreExceptions(), harvestTrackValidationPlots::str, AlCaHLTBitMon_QueryRunRegistry::string, helper::ScannerBase::test(), fwlite::EventBase::toBegin(), create_public_pileup_plots::vals, and fwlite::Scanner< Collection >::wantMore().

77  {
78  helper::ScannerBase scanner(objType);
79  scanner.setIgnoreExceptions(ignoreExceptions_);
80 
81  TObjArray *exprArray = TString(exprs).Tokenize(exprSep_);
82  int rowline = 0;
83  if (printFullEventId_) {
84  printf(" : %9s : %4s : %9s : %3s", "RUN", "LUMI", "EVENT", "#IT");
85  rowline += 3*4+9+4+9+3-1; // -1 as first char remain blank
86  } else {
87  printf(" : %5s : %3s", "EVENT", "#IT");
88  rowline += 3+6+3+3-1; // -1 as first char remain blank
89  }
90  for (int i = 0; i < exprArray->GetEntries(); ++i) {
91  TString str = ((TObjString *)(*exprArray)[i])->GetString();
92  std::string lb = str.Data();
93  std::string ex = str.Data();
94  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
95  lb = lb.substr(1,ex.find('=')-1);
96  ex = ex.substr(ex.find('=')+1);
97  }
98  scanner.addExpression(ex.c_str());
99  printf(" : %8s", (lb.size()>8 ? lb.substr(lb.size()-8) : lb).c_str()); // the rightmost part is usually the more interesting one
100  rowline += 3+8;
101  }
102  std::cout << " :" << std::endl;
103  rowline += 2;
104  delete exprArray;
105 
106  TString rule('-', rowline);
107  std::cout << " " << rule << " " << std::endl;
108 
109  if (strlen(cut)) scanner.setCut(cut);
110 
111  int iev = 0, line = 0;
112  for (event_->toBegin(); (iev != nmax) && !event_->atEnd(); ++iev, ++(*event_)) {
113  if (!selectEvent(*event_)) continue;
114  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
115  if (handle_.failedToGet()) {
116  if (ignoreExceptions_) continue;
117  }
118  const Collection & vals = *handle_;
119  for (size_t j = 0, n = vals.size(); j < n; ++j) {
120  if (!scanner.test(&vals[j])) continue;
121  if (printFullEventId_) {
123  printf(" : %9u : %4u : %9llu : %3lu", id.run(), id.luminosityBlock(), id.event(), (unsigned long)j);
124  } else {
125  printf(" : %5d : %3lu", iev, (unsigned long)j);
126  }
127  scanner.print(&vals[j]);
128  std::cout << " :" << std::endl;
129  if (++line == maxLinesToPrint_) {
130  line = 0;
131  if (!wantMore()) {
132  iev = nmax-1; // this is to exit the outer loop
133  break; // and this to exit the inner one
134  }
135  }
136  }
137  }
138  std::cout << std::endl;
139  }
bool ignoreExceptions_
Definition: Scanner.h:582
bool printFullEventId_
Definition: Scanner.h:581
TString exprSep_
Definition: Scanner.h:583
HandleT handle_
Definition: Scanner.h:584
std::string label_
Definition: Scanner.h:580
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
edm::TypeWithDict objType
Definition: Scanner.h:585
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:579
bool failedToGet() const
Returns true only if Handle was used in a &#39;get&#39; call and the data could not be found.
Definition: Handle.h:64
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
std::string instance_
Definition: Scanner.h:580
int maxLinesToPrint_
Definition: Scanner.h:591
bool wantMore() const
Definition: Scanner.h:592
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:570
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:580
virtual edm::EventAuxiliary const & eventAuxiliary() const =0
template<typename Collection >
bool fwlite::Scanner< Collection >::selectEvent ( const fwlite::EventBase ev) const
inline
template<typename Collection >
void fwlite::Scanner< Collection >::setExpressionSeparator ( TString  separator)
inline

Definition at line 563 of file Scanner.h.

References fwlite::Scanner< Collection >::exprSep_, and mps_merge::separator.

563 { exprSep_ = separator; }
string separator
Definition: mps_merge.py:77
TString exprSep_
Definition: Scanner.h:583
template<typename Collection >
void fwlite::Scanner< Collection >::setIgnoreExceptions ( bool  ignoreThem)
inline

Definition at line 564 of file Scanner.h.

References fwlite::Scanner< Collection >::ignoreExceptions_.

564 { ignoreExceptions_ = ignoreThem; }
bool ignoreExceptions_
Definition: Scanner.h:582
template<typename Collection >
void fwlite::Scanner< Collection >::setMaxEvents ( int  max)
inline
template<typename Collection >
void fwlite::Scanner< Collection >::setMaxLinesToPrint ( int  lines)
inline

Definition at line 565 of file Scanner.h.

References fwlite::Scanner< Collection >::maxLinesToPrint_.

565 { maxLinesToPrint_ = (lines > 0 ? lines : 2147483647); }
int maxLinesToPrint_
Definition: Scanner.h:591
template<typename Collection >
void fwlite::Scanner< Collection >::setPrintFullEventId ( bool  printIt = true)
inline

Definition at line 562 of file Scanner.h.

References fwlite::Scanner< Collection >::printFullEventId_.

562 { printFullEventId_ = printIt; }
bool printFullEventId_
Definition: Scanner.h:581
template<typename Collection >
bool fwlite::Scanner< Collection >::wantMore ( ) const
inlineprivate

Definition at line 592 of file Scanner.h.

References submit::answer.

Referenced by fwlite::Scanner< Collection >::scan().

592  {
593  // ask if user wants more
594  fprintf(stderr,"Type <CR> to continue or q to quit ==> ");
595  // read first char
596  int readch = getchar(), answer = readch;
597  // poll out remaining chars from buffer
598  while (readch != '\n' && readch != EOF) readch = getchar();
599  // check first char
600  return !(answer == 'q' || answer == 'Q');
601  }
answer
Definition: submit.py:44

Member Data Documentation

template<typename Collection >
fwlite::EventBase* fwlite::Scanner< Collection >::event_
private
template<typename Collection >
TObjArray fwlite::Scanner< Collection >::eventSelectors_
private
template<typename Collection >
TString fwlite::Scanner< Collection >::exprSep_
private
template<typename Collection >
HandleT fwlite::Scanner< Collection >::handle_
private
template<typename Collection >
bool fwlite::Scanner< Collection >::ignoreExceptions_
private
template<typename Collection >
std::string fwlite::Scanner< Collection >::instance_
private
template<typename Collection >
std::string fwlite::Scanner< Collection >::label_
private
template<typename Collection >
int fwlite::Scanner< Collection >::maxEvents_
private
template<typename Collection >
int fwlite::Scanner< Collection >::maxLinesToPrint_
private
template<typename Collection >
edm::TypeWithDict fwlite::Scanner< Collection >::objType
private
template<typename Collection >
bool fwlite::Scanner< Collection >::printFullEventId_
private
template<typename Collection >
std::string fwlite::Scanner< Collection >::process_
private