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 46 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 49 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 52 of file Scanner.h.

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

52  :
55  ignoreExceptions_(false),
56  exprSep_(":"),
57  maxEvents_(-1),
59  {
61  }
bool ignoreExceptions_
Definition: Scanner.h:584
static PFTauRenderPlugin instance
bool printFullEventId_
Definition: Scanner.h:583
TString exprSep_
Definition: Scanner.h:585
bool isRealData() const
Definition: EventBase.h:64
std::string label_
Definition: Scanner.h:582
edm::TypeWithDict objType
Definition: Scanner.h:587
fwlite::EventBase * event_
Definition: Scanner.h:581
std::string instance_
Definition: Scanner.h:582
int maxLinesToPrint_
Definition: Scanner.h:593
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:42

Member Function Documentation

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

Definition at line 569 of file Scanner.h.

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

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

Definition at line 570 of file Scanner.h.

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

570 { eventSelectors_.Clear(); }
TObjArray eventSelectors_
Definition: Scanner.h:589
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 147 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.

147  {
148  helper::ScannerBase scanner(objType);
149  scanner.setIgnoreExceptions(ignoreExceptions_);
150 
151  scanner.setCut(cut);
152 
153  size_t npass = 0;
154  int iev = 0;
155  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
156  if (maxEvents_ > -1 && iev > maxEvents_) break;
157  if (!selectEvent(*event_)) continue;
158  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
159  const Collection & vals = *handle_;
160  for (size_t j = 0, n = vals.size(); j < n; ++j) {
161  if (scanner.test(&vals[j])) npass++;
162  }
163  }
164  return npass;
165  }
bool ignoreExceptions_
Definition: Scanner.h:584
HandleT handle_
Definition: Scanner.h:586
std::string label_
Definition: Scanner.h:582
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:587
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
std::string instance_
Definition: Scanner.h:582
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:582
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 168 of file Scanner.h.

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

168  {
169  size_t npass = 0;
170  int iev = 0;
171  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
172  if (maxEvents_ > -1 && iev > maxEvents_) break;
173  if (selectEvent(*event_)) npass++;
174  }
175  return npass;
176  }
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
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 187 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().

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

229  {
230  TH1 *hist = 0;
231  if (htemplate != 0) {
232  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname,htemplate->GetName()) != 0)) htempDelete();
233  hist = (TH1*) hist->Clone(hname);
234  } else if (drawopt.Contains("SAME",TString::kIgnoreCase)) {
235  hist = getSameH1(hname);
236  }
237 
238  // if in the end we found no way to make "hist"
239  if (hist == 0) {
240  if (strcmp(hname, "htemp") == 0) htempDelete();
241  hist = new TH1F(hname, "", gEnv->GetValue("Hist.Binning.1D.x",100), 0, 0);
242  hist->SetCanExtend(TH1::kAllAxes);
243  }
244  hist->SetTitle((strlen(cut) ? TString(expr)+"{"+cut+"}" : TString(expr)));
245  hist->GetXaxis()->SetTitle(expr);
246  return draw(expr, cut, drawopt, hist);
247  }
TH1 * getSameH1(const char *hname)
Definition: Scanner.h:614
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
Definition: Scanner.h:187
void htempDelete()
Definition: Scanner.h:605
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 252 of file Scanner.h.

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

252  {
253  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
254  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
255  TH1 *hsame = getSameH1(hname);
256  return draw(expr, cut, drawopt, hsame);
257  }
258  if (strcmp(hname, "htemp") == 0) htempDelete();
259  TH1 * htemp = new TH1F(hname, expr, nbins, xlow, xhigh);
260  if (strlen(cut)) htemp->SetTitle(TString(expr)+"{"+cut+"}");
261  htemp->GetXaxis()->SetTitle(expr);
262  return draw(expr,cut,drawopt,htemp);
263  }
TH1 * getSameH1(const char *hname)
Definition: Scanner.h:614
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
Definition: Scanner.h:187
void htempDelete()
Definition: Scanner.h:605
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 266 of file Scanner.h.

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

266  {
267  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
268  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
269  TH1 *hsame = getSameH1(hname);
270  return draw(expr, cut, drawopt, hsame);
271  }
272  if (strcmp(hname, "htemp") == 0) htempDelete();
273  TH1 * htemp = new TH1F(hname, expr, nbins, xbins);
274  if (strlen(cut)) htemp->SetTitle(TString(expr)+"{"+cut+"}");
275  htemp->GetXaxis()->SetTitle(expr);
276  return draw(expr,cut,drawopt,htemp);
277  }
const double xbins[]
TH1 * getSameH1(const char *hname)
Definition: Scanner.h:614
TH1 * draw(const char *expr, const char *cut, TString drawopt, TH1 *hist)
Definition: Scanner.h:187
void htempDelete()
Definition: Scanner.h:605
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 346 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().

346  {
347  // prep the machinery
348  helper::ScannerBase scanner(objType);
349  scanner.setIgnoreExceptions(ignoreExceptions_);
350  if (!scanner.addExpression((const char *)xexpr)) return 0;
351  if (!scanner.addExpression((const char *)yexpr)) return 0;
352  if (strlen(cut)) scanner.setCut(cut);
353 
354  // check histo
355  if (hist == 0) {
356  std::cerr << "Method draw2D(xexpr, yexpr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other draw methods instead." << std::endl;
357  return 0;
358  }
359 
360  // fill histogram
361  int iev = 0;
362  for (event_->toBegin(), iev = 0; !event_->atEnd(); ++(*event_), ++iev) {
363  if (maxEvents_ > -1 && iev > maxEvents_) break;
364  if (!selectEvent(*event_)) continue;
365  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
366  const Collection & vals = *handle_;
367  for (size_t j = 0, n = vals.size(); j < n; ++j) {
368  scanner.fill2D(&vals[j], hist);
369  }
370  }
371 
372  if (!strlen(hist->GetTitle())) hist->SetTitle((strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
373  if (!strlen(hist->GetXaxis()->GetTitle())) hist->GetXaxis()->SetTitle(xexpr);
374  if (!strlen(hist->GetYaxis()->GetTitle())) hist->GetYaxis()->SetTitle(yexpr);
375  if (!TString(drawopt).Contains("GOFF",TString::kIgnoreCase)) hist->Draw(drawopt);
376  return hist;
377  }
bool ignoreExceptions_
Definition: Scanner.h:584
HandleT handle_
Definition: Scanner.h:586
std::string label_
Definition: Scanner.h:582
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:587
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
std::string instance_
Definition: Scanner.h:582
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:582
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 380 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.

380  {
381  TH2 *hist = 0;
382  if (htemplate != 0) {
383  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname,htemplate->GetName()) != 0)) htempDelete();
384  hist = (TH2*) hist->Clone(hname);
385  } else if (drawopt.Contains("SAME",TString::kIgnoreCase)) {
386  hist = getSameH2(hname);
387  }
388 
389  // if in the end we found no way to make "hist"
390  if (hist == 0) {
391  // prep the machinery
392  helper::ScannerBase scanner(objType);
393  scanner.setIgnoreExceptions(ignoreExceptions_);
394  if (!scanner.addExpression((const char *)xexpr)) return 0;
395  if (!scanner.addExpression((const char *)yexpr)) return 0;
396  if (strlen(cut)) scanner.setCut(cut);
397 
398  if (strcmp(hname, "htemp") == 0) htempDelete();
399  // ok this is much more a hack than for the 1D case
400  double xmin = 0, xmax = -1, ymin = 0, ymax = -1; int iev;
401  for (event_->toBegin(), iev = 0; !event_->atEnd(); ++(*event_), ++iev) {
402  if (maxEvents_ > -1 && iev > maxEvents_) break;
403  if (!selectEvent(*event_)) continue;
404  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
405  const Collection & vals = *handle_;
406  for (size_t j = 0, n = vals.size(); j < n; ++j) {
407  if (!scanner.test(&vals[j])) continue;
408  double x = scanner.eval(&vals[j],0);
409  double y = scanner.eval(&vals[j],1);
410  if ((xmax == -1) || (x >= xmax)) xmax = x;
411  if ((xmin == 0) || (x <= xmin)) xmin = x;
412  if ((ymax == -1) || (y >= ymax)) ymax = y;
413  if ((ymin == 0) || (y <= ymin)) ymin = y;
414  }
415  }
416  hist = new TH2F(hname, "",
417  gEnv->GetValue("Hist.Binning.2D.x",20), xmin, xmax,
418  gEnv->GetValue("Hist.Binning.2D.y",20), ymin, ymax);
419  }
420  return draw2D(xexpr, yexpr, cut, drawopt, hist);
421  }
bool ignoreExceptions_
Definition: Scanner.h:584
TH2 * getSameH2(const char *hname)
Definition: Scanner.h:630
HandleT handle_
Definition: Scanner.h:586
std::string label_
Definition: Scanner.h:582
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:587
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
std::string instance_
Definition: Scanner.h:582
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:346
void htempDelete()
Definition: Scanner.h:605
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:582
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 424 of file Scanner.h.

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

426  {
427  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
428  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
429  TH2 *hsame = getSameH2(hname);
430  return draw2D(xexpr, yexpr, cut, drawopt, hsame);
431  }
432  if (strcmp(hname, "htemp") == 0) htempDelete();
433  TH2 * htemp = new TH2F("htemp", "", xbins, xlow, xhigh, ybins,ylow,yhigh);
434  return draw2D(xexpr,yexpr,cut,drawopt,htemp);
435  }
const double xbins[]
TH2 * getSameH2(const char *hname)
Definition: Scanner.h:630
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:346
void htempDelete()
Definition: Scanner.h:605
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 438 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().

438  {
439  // prep the machinery
440  helper::ScannerBase scanner(objType);
441  scanner.setIgnoreExceptions(ignoreExceptions_);
442  if (!scanner.addExpression((const char *)xexpr)) return 0;
443  if (!scanner.addExpression((const char *)yexpr)) return 0;
444  if (strlen(cut)) scanner.setCut(cut);
445 
446  // make graph, if needed
447  if (graph == 0) {
448  graph = new TGraph();
449  graph->SetNameTitle("htemp", (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
450  }
451 
452  // fill graph
453  int iev = 0;
454  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
455  if (maxEvents_ > -1 && iev > maxEvents_) break;
456  if (!selectEvent(*event_)) continue;
457  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
458  const Collection & vals = *handle_;
459  for (size_t j = 0, n = vals.size(); j < n; ++j) {
460  scanner.fillGraph(&vals[j], graph);
461  }
462  }
463 
464  if (!strlen(graph->GetTitle())) graph->SetTitle((strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
465  if (!strlen(graph->GetXaxis()->GetTitle())) graph->GetXaxis()->SetTitle(xexpr);
466  if (!strlen(graph->GetYaxis()->GetTitle())) graph->GetYaxis()->SetTitle(yexpr);
467  if (!TString(drawopt).Contains("GOFF",TString::kIgnoreCase)) graph->Draw(drawopt);
468  return graph;
469  }
bool ignoreExceptions_
Definition: Scanner.h:584
HandleT handle_
Definition: Scanner.h:586
std::string label_
Definition: Scanner.h:582
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:587
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
std::string instance_
Definition: Scanner.h:582
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:582
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 472 of file Scanner.h.

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

472  {
473  if (strcmp(gname, "htemp") == 0) htempDelete();
474  TGraph *graph = new TGraph();
475  graph->SetNameTitle(gname, (strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
476  return drawGraph(xexpr,yexpr,cut,drawopt,graph);
477  }
TGraph * drawGraph(TString xexpr, TString yexpr, const char *cut, TString drawopt, TGraph *graph)
Definition: Scanner.h:438
void htempDelete()
Definition: Scanner.h:605
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 281 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().

281  {
282  // prep the machinery
283  helper::ScannerBase scanner(objType);
284  scanner.setIgnoreExceptions(ignoreExceptions_);
285  if (!scanner.addExpression(xexpr.Data())) return 0;
286  if (!scanner.addExpression(yexpr.Data())) return 0;
287  if (strlen(cut)) scanner.setCut(cut);
288 
289  // check histo
290  if (hist == 0) {
291  std::cerr << "Method drawProf(xexpr, yexpr, cut, drawopt, hist) cannot be called with null 'hist'. Use the other draw methods instead." << std::endl;
292  return 0;
293  }
294 
295  // fill histogram
296  int iev = 0;
297  for (event_->toBegin(); !event_->atEnd(); ++(*event_), ++iev) {
298  if (maxEvents_ > -1 && iev > maxEvents_) break;
299  if (!selectEvent(*event_)) continue;
300  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
301  const Collection & vals = *handle_;
302  for (size_t j = 0, n = vals.size(); j < n; ++j) {
303  scanner.fillProf(&vals[j], hist);
304  }
305  }
306 
307  if (!strlen(hist->GetTitle())) hist->SetTitle((strlen(cut) ? yexpr+":"+xexpr+"{"+cut+"}" : yexpr+":"+xexpr));
308  if (!strlen(hist->GetXaxis()->GetTitle())) hist->GetXaxis()->SetTitle(xexpr);
309  if (!strlen(hist->GetYaxis()->GetTitle())) hist->GetYaxis()->SetTitle(yexpr);
310  if (!TString(drawopt).Contains("GOFF",TString::kIgnoreCase)) hist->Draw(drawopt);
311  return hist;
312  }
bool ignoreExceptions_
Definition: Scanner.h:584
HandleT handle_
Definition: Scanner.h:586
std::string label_
Definition: Scanner.h:582
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:587
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
std::string instance_
Definition: Scanner.h:582
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:582
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 314 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().

314  {
315  TProfile *hist = 0;
316  if (htemplate != 0) {
317  if ((strcmp(hname, "htemp") == 0) && (strcmp(hname,htemplate->GetName() )!= 0)) htempDelete();
318  hist = (TProfile*) hist->Clone(hname);
319  } else if (drawopt.Contains("SAME",TString::kIgnoreCase)) {
320  hist = getSameProf(hname);
321  }
322 
323  // if in the end we found no way to make "hist"
324  if (hist == 0) {
325  if (strcmp(hname, "htemp") == 0) htempDelete();
326  hist = new TProfile(hname, "", gEnv->GetValue("Hist.Binning.1D.x",100), 0., 0.);
327  hist->SetCanExtend(TH1::kAllAxes);
328  }
329  return drawProf(xexpr, yexpr, cut, drawopt, hist);
330  }
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:281
TProfile * getSameProf(const char *hname)
Definition: Scanner.h:646
void htempDelete()
Definition: Scanner.h:605
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 333 of file Scanner.h.

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

333  {
334  if (TString(drawopt).Contains("SAME",TString::kIgnoreCase)) {
335  std::cerr << "Binning is ignored when 'SAME' is specified." << std::endl;
336  TProfile *hsame = getSameProf(hname);
337  return drawProf(xexpr, yexpr, cut, drawopt, hsame);
338  }
339  if (strcmp(hname, "htemp") == 0) htempDelete();
340  TProfile * htemp = new TProfile(hname, "", bins, xlow, xhigh);
341  return drawProf(xexpr,yexpr,cut,drawopt,htemp);
342  }
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:281
TProfile * getSameProf(const char *hname)
Definition: Scanner.h:646
void htempDelete()
Definition: Scanner.h:605
template<typename Collection >
TObjArray& fwlite::Scanner< Collection >::eventSelectors ( )
inline

Definition at line 571 of file Scanner.h.

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

571 { return eventSelectors_; }
TObjArray eventSelectors_
Definition: Scanner.h:589
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 487 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(), str, AlCaHLTBitMon_QueryRunRegistry::string, helper::ScannerBase::test(), fwlite::EventBase::toBegin(), create_public_pileup_plots::vals, and JetChargeProducer_cfi::var.

487  {
488  helper::ScannerBase scanner(objType);
489  scanner.setIgnoreExceptions(ignoreExceptions_);
490 
491  RooArgList vars;
492  TObjArray *exprArray = TString(realvars).Tokenize(exprSep_);
493  TObjArray *catArray = TString(boolvars).Tokenize(exprSep_);
494  int nreals = exprArray->GetEntries();
495  int nbools = catArray->GetEntries();
496  for (int i = 0; i < nreals; ++i) {
497  TString str = ((TObjString *)(*exprArray)[i])->GetString();
498  std::string lb = str.Data();
499  std::string ex = str.Data();
500  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
501  lb = lb.substr(1,ex.find('=')-1);
502  ex = ex.substr(ex.find('=')+1);
503  }
504  if (!scanner.addExpression(ex.c_str())) {
505  std::cerr << "Filed to define real variable '" << lb << "', expr = '" << ex << "'" << std::endl;
506  return 0;
507  }
508  // FIXME: I have to leave it dangling on the HEAP otherwise ROOT segfaults...
509  RooRealVar *var = new RooRealVar(lb.c_str(),lb.c_str(), 0.0);
510  vars.add(*var);
511  }
512  for (int i = 0; i < nbools; ++i) {
513  TString str = ((TObjString *)(*catArray)[i])->GetString();
514  std::string lb = str.Data();
515  std::string ex = str.Data();
516  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
517  lb = lb.substr(1,ex.find('=')-1);
518  ex = ex.substr(ex.find('=')+1);
519  }
520  if (!scanner.addExtraCut(ex.c_str())) {
521  std::cerr << "Filed to define bool variable '" << lb << "', cut = '" << ex << "'" << std::endl;
522  return 0;
523  }
524  RooCategory *cat = new RooCategory(lb.c_str(), lb.c_str());
525  cat->defineType("fail",0);
526  cat->defineType("pass",1);
527  vars.add(*cat);
528  }
529 
530  RooDataSet *ds = new RooDataSet(name, name, vars);
531 
532  if (strlen(cut)) scanner.setCut(cut);
533  int iev = 0;
534  for (event_->toBegin(); !event_->atEnd(); ++iev, ++(*event_)) {
535  if (maxEvents_ > -1 && iev > maxEvents_) break;
536  if (!selectEvent(*event_)) continue;
537  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
538  if (handle_.failedToGet()) {
539  if (ignoreExceptions_) continue;
540  }
541  const Collection & vals = *handle_;
542  for (size_t j = 0, n = vals.size(); j < n; ++j) {
543  if (!scanner.test(&vals[j])) continue;
544  for (int i = 0; i < nreals; ++i) {
545  RooRealVar *var = (RooRealVar *)vars.at(i);
546  var->setVal(scanner.eval(&vals[j], i));
547  }
548  for (int i = 0; i < nbools; ++i) {
549  RooCategory *cat = (RooCategory*) vars.at(i+nreals);
550  cat->setIndex(int(scanner.test(&vals[j], i+1))); // 0 is the event selection cut
551  }
552  ds->add(vars);
553  }
554  }
555 
556  delete exprArray;
557  delete catArray;
558 
559  return ds;
560  }
bool ignoreExceptions_
Definition: Scanner.h:584
TString exprSep_
Definition: Scanner.h:585
HandleT handle_
Definition: Scanner.h:586
std::string label_
Definition: Scanner.h:582
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:587
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
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:582
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
#define str(s)
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:582
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 614 of file Scanner.h.

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

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

614  {
615  if (gDirectory && gDirectory->Get("htemp") != 0 &&
616  gDirectory->Get("htemp")->IsA()->InheritsFrom(TH1::Class())) {
617  TH1 *hist = (TH1*) ((TH1*) gDirectory->Get("htemp"))->Clone(hname);
618  hist->Reset();
619  hist->SetLineColor(kBlack);
620  hist->SetMarkerColor(kBlack);
621  return hist;
622  } else {
623  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
624  return 0;
625  }
626  }
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 630 of file Scanner.h.

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

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

630  {
631  if (gDirectory && gDirectory->Get("htemp") != 0 &&
632  gDirectory->Get("htemp")->IsA()->InheritsFrom(TH2::Class())) {
633  TH2 *hist = (TH2*) ((TH2*) gDirectory->Get("htemp"))->Clone(hname);
634  hist->Reset();
635  hist->SetLineColor(kBlack);
636  hist->SetMarkerColor(kBlack);
637  return hist;
638  } else {
639  std::cerr << "There is no 'htemp' histogram from which to 'SAME'." << std::endl;
640  return 0;
641  }
642  }
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 646 of file Scanner.h.

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

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

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

Definition at line 605 of file Scanner.h.

References GetRecoTauVFromDQM_MC_cff::obj.

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

605  {
606  if (gDirectory) {
607  TObject *obj = gDirectory->Get("htemp");
608  if (obj) obj->Delete();
609  }
610  }
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 79 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(), str, AlCaHLTBitMon_QueryRunRegistry::string, helper::ScannerBase::test(), fwlite::EventBase::toBegin(), create_public_pileup_plots::vals, and fwlite::Scanner< Collection >::wantMore().

79  {
80  helper::ScannerBase scanner(objType);
81  scanner.setIgnoreExceptions(ignoreExceptions_);
82 
83  TObjArray *exprArray = TString(exprs).Tokenize(exprSep_);
84  int rowline = 0;
85  if (printFullEventId_) {
86  printf(" : %9s : %4s : %9s : %3s", "RUN", "LUMI", "EVENT", "#IT");
87  rowline += 3*4+9+4+9+3-1; // -1 as first char remain blank
88  } else {
89  printf(" : %5s : %3s", "EVENT", "#IT");
90  rowline += 3+6+3+3-1; // -1 as first char remain blank
91  }
92  for (int i = 0; i < exprArray->GetEntries(); ++i) {
93  TString str = ((TObjString *)(*exprArray)[i])->GetString();
94  std::string lb = str.Data();
95  std::string ex = str.Data();
96  if ((ex[0] == '@') && (ex.find('=') != std::string::npos)) {
97  lb = lb.substr(1,ex.find('=')-1);
98  ex = ex.substr(ex.find('=')+1);
99  }
100  scanner.addExpression(ex.c_str());
101  printf(" : %8s", (lb.size()>8 ? lb.substr(lb.size()-8) : lb).c_str()); // the rightmost part is usually the more interesting one
102  rowline += 3+8;
103  }
104  std::cout << " :" << std::endl;
105  rowline += 2;
106  delete exprArray;
107 
108  TString rule('-', rowline);
109  std::cout << " " << rule << " " << std::endl;
110 
111  if (strlen(cut)) scanner.setCut(cut);
112 
113  int iev = 0, line = 0;
114  for (event_->toBegin(); (iev != nmax) && !event_->atEnd(); ++iev, ++(*event_)) {
115  if (!selectEvent(*event_)) continue;
116  handle_.getByLabel(*event_, label_.c_str(), instance_.c_str(), process_.c_str());
117  if (handle_.failedToGet()) {
118  if (ignoreExceptions_) continue;
119  }
120  const Collection & vals = *handle_;
121  for (size_t j = 0, n = vals.size(); j < n; ++j) {
122  if (!scanner.test(&vals[j])) continue;
123  if (printFullEventId_) {
125  printf(" : %9u : %4u : %9llu : %3lu", id.run(), id.luminosityBlock(), id.event(), (unsigned long)j);
126  } else {
127  printf(" : %5d : %3lu", iev, (unsigned long)j);
128  }
129  scanner.print(&vals[j]);
130  std::cout << " :" << std::endl;
131  if (++line == maxLinesToPrint_) {
132  line = 0;
133  if (!wantMore()) {
134  iev = nmax-1; // this is to exit the outer loop
135  break; // and this to exit the inner one
136  }
137  }
138  }
139  }
140  std::cout << std::endl;
141  }
bool ignoreExceptions_
Definition: Scanner.h:584
bool printFullEventId_
Definition: Scanner.h:583
TString exprSep_
Definition: Scanner.h:585
HandleT handle_
Definition: Scanner.h:586
std::string label_
Definition: Scanner.h:582
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:587
virtual bool atEnd() const =0
fwlite::EventBase * event_
Definition: Scanner.h:581
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:582
int maxLinesToPrint_
Definition: Scanner.h:593
bool wantMore() const
Definition: Scanner.h:594
bool selectEvent(const fwlite::EventBase &ev) const
Definition: Scanner.h:572
#define str(s)
virtual EventBase const & toBegin()=0
std::string process_
Definition: Scanner.h:582
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 565 of file Scanner.h.

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

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

Definition at line 566 of file Scanner.h.

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

566 { ignoreExceptions_ = ignoreThem; }
bool ignoreExceptions_
Definition: Scanner.h:584
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 567 of file Scanner.h.

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

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

Definition at line 564 of file Scanner.h.

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

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

Definition at line 594 of file Scanner.h.

References submit::answer.

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

594  {
595  // ask if user wants more
596  fprintf(stderr,"Type <CR> to continue or q to quit ==> ");
597  // read first char
598  int readch = getchar(), answer = readch;
599  // poll out remaining chars from buffer
600  while (readch != '\n' && readch != EOF) readch = getchar();
601  // check first char
602  return !(answer == 'q' || answer == 'Q');
603  }
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