CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EvtPlaneProducer Class Reference

#include <RecoHI/EvtPlaneProducer/src/EvtPlaneProducer.cc>

Inheritance diagram for EvtPlaneProducer:
edm::stream::EDProducer<>

Public Member Functions

 EvtPlaneProducer (const edm::ParameterSet &)
 
 ~EvtPlaneProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

double caloCentRef_
 
double caloCentRefWidth_
 
edm::Handle< CaloTowerCollectioncaloCollection_
 
edm::InputTag caloTag_
 
edm::EDGetTokenT< CaloTowerCollectioncaloToken
 
edm::Handle< std::vector< reco::CastorTower > > castorCollection_
 
edm::InputTag castorTag_
 
edm::EDGetTokenT< std::vector< reco::CastorTower > > castorToken
 
int CentBinCompression_
 
edm::InputTag centralityBinTag_
 
edm::EDGetTokenT< int > centralityBinToken
 
std::string centralityLabel_
 
std::string centralityMC_
 
std::string centralityVariable_
 
double chi2_
 
double dzerr_
 
HiEvtPlaneFlattenflat [NumEPNames]
 
int FlatOrder_
 
edm::ESWatcher< HeavyIonRPRcdhirpWatcher
 
edm::ESWatcher< HeavyIonRcdhiWatcher
 
bool loadDB_
 
double maxet_
 
double maxpt_
 
double maxvtx_
 
double minet_
 
double minpt_
 
double minvtx_
 
double nCentBins_
 
int NumFlatBins_
 
GenPlanerp [NumEPNames]
 
edm::Handle< reco::TrackCollectiontrackCollection_
 
edm::InputTag trackTag_
 
edm::EDGetTokenT< reco::TrackCollectiontrackToken
 
edm::Handle< std::vector< reco::Vertex > > vertex_
 
edm::InputTag vertexTag_
 
edm::EDGetTokenT< std::vector< reco::Vertex > > vertexToken
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 167 of file EvtPlaneProducer.cc.

Constructor & Destructor Documentation

EvtPlaneProducer::EvtPlaneProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 223 of file EvtPlaneProducer.cc.

References caloTag_, caloToken, castorTag_, castorToken, centralityBinTag_, centralityBinToken, centralityLabel_, centralityMC_, centralityVariable_, hi::EPEtaMax1, hi::EPEtaMax2, hi::EPEtaMin1, hi::EPEtaMin2, hi::EPNames, hi::EPOrder, edm::ParameterSet::exists(), flat, FlatOrder_, edm::ParameterSet::getParameter(), mps_fire::i, HiEvtPlaneFlatten::init(), nCentBins_, hi::NumEPNames, NumFlatBins_, rp, AlCaHLTBitMon_QueryRunRegistry::string, trackTag_, trackToken, vertexTag_, and vertexToken.

224  : centralityVariable_(iConfig.getParameter<std::string>("centralityVariable")),
225  centralityBinTag_(iConfig.getParameter<edm::InputTag>("centralityBinTag")),
226  vertexTag_(iConfig.getParameter<edm::InputTag>("vertexTag")),
227  caloTag_(iConfig.getParameter<edm::InputTag>("caloTag")),
228  castorTag_(iConfig.getParameter<edm::InputTag>("castorTag")),
229  trackTag_(iConfig.getParameter<edm::InputTag>("trackTag")),
230  loadDB_(iConfig.getParameter<bool>("loadDB")),
231  minet_(iConfig.getParameter<double>("minet")),
232  maxet_(iConfig.getParameter<double>("maxet")),
233  minpt_(iConfig.getParameter<double>("minpt")),
234  maxpt_(iConfig.getParameter<double>("maxpt")),
235  minvtx_(iConfig.getParameter<double>("minvtx")),
236  maxvtx_(iConfig.getParameter<double>("maxvtx")),
237  dzerr_(iConfig.getParameter<double>("dzerr")),
238  chi2_(iConfig.getParameter<double>("chi2")),
239  FlatOrder_(iConfig.getParameter<int>("FlatOrder")),
240  NumFlatBins_(iConfig.getParameter<int>("NumFlatBins")),
241  caloCentRef_(iConfig.getParameter<double>("caloCentRef")),
242  caloCentRefWidth_(iConfig.getParameter<double>("caloCentRefWidth")),
243  CentBinCompression_(iConfig.getParameter<int>("CentBinCompression")) {
244  nCentBins_ = 200.;
245 
246  if (iConfig.exists("nonDefaultGlauberModel")) {
247  centralityMC_ = iConfig.getParameter<std::string>("nonDefaultGlauberModel");
248  }
250 
251  centralityBinToken = consumes<int>(centralityBinTag_);
252 
253  vertexToken = consumes<std::vector<reco::Vertex>>(vertexTag_);
254 
255  caloToken = consumes<CaloTowerCollection>(caloTag_);
256 
257  castorToken = consumes<std::vector<reco::CastorTower>>(castorTag_);
258 
259  trackToken = consumes<reco::TrackCollection>(trackTag_);
260 
261  produces<reco::EvtPlaneCollection>();
262  for (int i = 0; i < NumEPNames; i++) {
263  rp[i] = new GenPlane(EPNames[i], EPEtaMin1[i], EPEtaMax1[i], EPEtaMin2[i], EPEtaMax2[i], EPOrder[i]);
264  }
265  for (int i = 0; i < NumEPNames; i++) {
266  flat[i] = new HiEvtPlaneFlatten();
268  }
269 }
void init(int order, int nbins, std::string tag, int vord)
edm::InputTag centralityBinTag_
std::string centralityLabel_
T getParameter(std::string const &) const
const int EPOrder[]
const double EPEtaMin2[]
edm::InputTag caloTag_
std::string centralityMC_
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< int > centralityBinToken
edm::EDGetTokenT< std::vector< reco::CastorTower > > castorToken
const std::string EPNames[]
const double EPEtaMin1[]
GenPlane * rp[NumEPNames]
HiEvtPlaneFlatten * flat[NumEPNames]
edm::EDGetTokenT< std::vector< reco::Vertex > > vertexToken
edm::InputTag vertexTag_
edm::EDGetTokenT< reco::TrackCollection > trackToken
std::string centralityVariable_
const double EPEtaMax2[]
static const int NumEPNames
edm::EDGetTokenT< CaloTowerCollection > caloToken
edm::InputTag castorTag_
const double EPEtaMax1[]
edm::InputTag trackTag_
EvtPlaneProducer::~EvtPlaneProducer ( )
override

Definition at line 271 of file EvtPlaneProducer.cc.

References flat, mps_fire::i, and hi::NumEPNames.

271  {
272  // do anything here that needs to be done at desctruction time
273  // (e.g. close files, deallocate resources etc.)
274  for (int i = 0; i < NumEPNames; i++) {
275  delete flat[i];
276  }
277 }
HiEvtPlaneFlatten * flat[NumEPNames]
static const int NumEPNames

Member Function Documentation

void EvtPlaneProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 284 of file EvtPlaneProducer.cc.

References cms::dd::accepted(), hi::GenPlane::addParticle(), newFWLiteAna::bin, caloCentRef_, caloCentRefWidth_, caloCollection_, caloToken, castorCollection_, castorToken, CentBinCompression_, centralityBinToken, centralityLabel_, edm::ESWatcher< T >::check(), chi2_, funct::cos(), cuy::cv, HLTMuonOfflineAnalyzer_cfi::d0, dqmiodatasetharvest::db, DEFINE_FWK_MODULE, PVValHelper::dz, dzerr_, hi::EPDet, hi::EPOrder, flat, edm::EventSetup::get(), hi::GenPlane::getAngle(), edm::Event::getByToken(), HiEvtPlaneFlatten::getEtScale(), HiEvtPlaneFlatten::getW(), CaloTPGTranscoder_cfi::HF, hirpWatcher, hiWatcher, mps_fire::i, fastTrackerRecHitType::isPixel(), LoadEPDB::IsSuccess(), edm::HandleBase::isValid(), dqmiolumiharvest::j, loadDB_, CentralityTable::m_table, HiEvtPlane_cfi::maxet, maxet_, HiEvtPlane_cfi::maxpt, maxpt_, hi::maxTransverse, maxvtx_, HiEvtPlane_cfi::minet, minet_, HiEvtPlane_cfi::minpt, minpt_, hi::minTransverse, hi::MomConsWeight, eostools::move(), nCentBins_, hi::NumEPNames, edm::Handle< T >::product(), edm::Event::put(), reset(), rp, HiEvtPlaneFlatten::setCaloCentRefBins(), funct::sin(), mathSSE::sqrt(), pfDeepBoostedJetPreprocessParams_cfi::sv, trackCollection_, align::Tracker, trackToken, parallelization::uint, vertex_, vertexToken, and w.

284  {
285  using namespace edm;
286  using namespace std;
287  using namespace reco;
288 
289  if (loadDB_ && (hiWatcher.check(iSetup) || hirpWatcher.check(iSetup))) {
290  //
291  //Get Size of Centrality Table
292  //
294  iSetup.get<HeavyIonRcd>().get(centralityLabel_, centDB_);
295  nCentBins_ = centDB_->m_table.size();
296  for (int i = 0; i < NumEPNames; i++) {
297  if (caloCentRef_ > 0) {
298  int minbin = (caloCentRef_ - caloCentRefWidth_ / 2.) * nCentBins_ / 100.;
299  int maxbin = (caloCentRef_ + caloCentRefWidth_ / 2.) * nCentBins_ / 100.;
300  minbin /= CentBinCompression_;
301  maxbin /= CentBinCompression_;
302  if (minbin > 0 && maxbin >= minbin) {
303  if (EPDet[i] == HF || EPDet[i] == Castor)
304  flat[i]->setCaloCentRefBins(minbin, maxbin);
305  }
306  }
307  }
308 
309  //
310  //Get flattening parameter file.
311  //
312  if (loadDB_) {
313  edm::ESHandle<RPFlatParams> flatparmsDB_;
314  iSetup.get<HeavyIonRPRcd>().get(flatparmsDB_);
315  LoadEPDB db(flatparmsDB_, flat);
316  if (!db.IsSuccess()) {
317  loadDB_ = kFALSE;
318  }
319  }
320 
321  } //rp record change
322 
323  //
324  //Get Centrality
325  //
326  int bin = 0;
327  if (loadDB_) {
328  edm::Handle<int> cbin_;
329  iEvent.getByToken(centralityBinToken, cbin_);
330  int cbin = *cbin_;
331  bin = cbin / CentBinCompression_;
332  }
333  //
334  //Get Vertex
335  //
336  int vs_sell = 0.;
337  float vzr_sell;
338  iEvent.getByToken(vertexToken, vertex_);
339  const reco::VertexCollection *vertices3 = nullptr;
340  if (vertex_.isValid()) {
341  vertices3 = vertex_.product();
342  vs_sell = vertices3->size();
343  }
344  if (vs_sell > 0) {
345  vzr_sell = vertices3->begin()->z();
346  } else
347  vzr_sell = -999.9;
348  //
349  for (int i = 0; i < NumEPNames; i++)
350  rp[i]->reset();
351  if (vzr_sell < minvtx_ or vzr_sell > maxvtx_)
352  return;
353 
354  //calorimetry part
355 
356  double tower_eta, tower_phi;
357  double tower_energyet, tower_energyet_e, tower_energyet_h;
358 
360 
361  if (caloCollection_.isValid()) {
363  tower_eta = j->eta();
364  tower_phi = j->phi();
365  tower_energyet_e = j->emEt();
366  tower_energyet_h = j->hadEt();
367  tower_energyet = tower_energyet_e + tower_energyet_h;
368  double minet = minet_;
369  double maxet = maxet_;
370  for (int i = 0; i < NumEPNames; i++) {
371  if (minet_ < 0)
372  minet = minTransverse[i];
373  if (maxet_ < 0)
374  maxet = maxTransverse[i];
375  if (tower_energyet < minet)
376  continue;
377  if (tower_energyet > maxet)
378  continue;
379  if (EPDet[i] == HF) {
380  double w = tower_energyet;
381  if (loadDB_)
382  w = tower_energyet * flat[i]->getEtScale(vzr_sell, bin);
383  if (EPOrder[i] == 1) {
384  if (MomConsWeight[i][0] == 'y' && loadDB_) {
385  w = flat[i]->getW(tower_energyet, vzr_sell, bin);
386  }
387  if (tower_eta < 0)
388  w = -w;
389  }
390  rp[i]->addParticle(w, tower_energyet, sin(EPOrder[i] * tower_phi), cos(EPOrder[i] * tower_phi), tower_eta);
391  }
392  }
393  }
394  }
395 
396  //Castor part
397 
399 
400  if (castorCollection_.isValid()) {
401  for (std::vector<reco::CastorTower>::const_iterator j = castorCollection_->begin(); j != castorCollection_->end();
402  j++) {
403  tower_eta = j->eta();
404  tower_phi = j->phi();
405  tower_energyet = j->et();
406  double minet = minet_;
407  double maxet = maxet_;
408  for (int i = 0; i < NumEPNames; i++) {
409  if (EPDet[i] == Castor) {
410  if (minet_ < 0)
411  minet = minTransverse[i];
412  if (maxet_ < 0)
413  maxet = maxTransverse[i];
414  if (tower_energyet < minet)
415  continue;
416  if (tower_energyet > maxet)
417  continue;
418  double w = tower_energyet;
419  if (EPOrder[i] == 1) {
420  if (MomConsWeight[i][0] == 'y' && loadDB_) {
421  w = flat[i]->getW(tower_energyet, vzr_sell, bin);
422  }
423  if (tower_eta < 0)
424  w = -w;
425  }
426  rp[i]->addParticle(w, tower_energyet, sin(EPOrder[i] * tower_phi), cos(EPOrder[i] * tower_phi), tower_eta);
427  }
428  }
429  }
430  }
431 
432  //Tracking part
433 
434  double track_eta;
435  double track_phi;
436  double track_pt;
437 
438  double vzErr2 = 0.0, vxyErr = 0.0;
439  math::XYZPoint vtxPoint(0.0, 0.0, 0.0);
440  if (vertex_.isValid() && !vertex_->empty()) {
441  vtxPoint = vertex_->begin()->position();
442  vzErr2 = (vertex_->begin()->zError()) * (vertex_->begin()->zError());
443  vxyErr = vertex_->begin()->xError() * vertex_->begin()->yError();
444  }
445 
447  if (trackCollection_.isValid()) {
448  for (reco::TrackCollection::const_iterator j = trackCollection_->begin(); j != trackCollection_->end(); j++) {
449  bool accepted = true;
450  bool isPixel = false;
451  // determine if the track is a pixel track
452  if (j->numberOfValidHits() < 7)
453  isPixel = true;
454 
455  // determine the vertex significance
456  double d0 = 0.0, dz = 0.0, d0sigma = 0.0, dzsigma = 0.0;
457  d0 = -1. * j->dxy(vtxPoint);
458  dz = j->dz(vtxPoint);
459  d0sigma = sqrt(j->d0Error() * j->d0Error() + vxyErr);
460  dzsigma = sqrt(j->dzError() * j->dzError() + vzErr2);
461 
462  // cuts for pixel tracks
463  if (isPixel) {
464  // dz significance cut
465  if (fabs(dz / dzsigma) > dzerr_)
466  accepted = false;
467  // chi2/ndof cut
468  if (j->normalizedChi2() > chi2_)
469  accepted = false;
470  }
471  // cuts for full tracks
472  if (!isPixel) {
473  // dz and d0 significance cuts
474  if (fabs(dz / dzsigma) > 3)
475  accepted = false;
476  if (fabs(d0 / d0sigma) > 3)
477  accepted = false;
478  // pt resolution cut
479  if (j->ptError() / j->pt() > 0.1)
480  accepted = false;
481  // number of valid hits cut
482  if (j->numberOfValidHits() < 12)
483  accepted = false;
484  }
485  if (accepted) {
486  track_eta = j->eta();
487  track_phi = j->phi();
488  track_pt = j->pt();
489  double minpt = minpt_;
490  double maxpt = maxpt_;
491  for (int i = 0; i < NumEPNames; i++) {
492  if (minpt_ < 0)
493  minpt = minTransverse[i];
494  if (maxpt_ < 0)
495  maxpt = maxTransverse[i];
496  if (track_pt < minpt)
497  continue;
498  if (track_pt > maxpt)
499  continue;
500  if (EPDet[i] == Tracker) {
501  double w = track_pt;
502  if (w > 2.5)
503  w = 2.0; //v2 starts decreasing above ~2.5 GeV/c
504  if (EPOrder[i] == 1) {
505  if (MomConsWeight[i][0] == 'y' && loadDB_) {
506  w = flat[i]->getW(track_pt, vzr_sell, bin);
507  }
508  if (track_eta < 0)
509  w = -w;
510  }
511  rp[i]->addParticle(w, track_pt, sin(EPOrder[i] * track_phi), cos(EPOrder[i] * track_phi), track_eta);
512  }
513  }
514  }
515  } //end for
516  }
517 
518  auto evtplaneOutput = std::make_unique<EvtPlaneCollection>();
519 
520  double ang = -10;
521  double sv = 0;
522  double cv = 0;
523  double svNoWgt = 0;
524  double cvNoWgt = 0;
525 
526  double wv = 0;
527  double wv2 = 0;
528  double pe = 0;
529  double pe2 = 0;
530  uint epmult = 0;
531 
532  for (int i = 0; i < NumEPNames; i++) {
533  rp[i]->getAngle(ang, sv, cv, svNoWgt, cvNoWgt, wv, wv2, pe, pe2, epmult);
534  evtplaneOutput->push_back(EvtPlane(i, 0, ang, sv, cv, wv, wv2, pe, pe2, epmult));
535  evtplaneOutput->back().addLevel(3, 0., svNoWgt, cvNoWgt);
536  }
537 
538  iEvent.put(std::move(evtplaneOutput));
539 }
std::string centralityLabel_
const int EPOrder[]
edm::Handle< CaloTowerCollection > caloCollection_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
std::vector< CBin > m_table
const int EPDet[]
edm::Handle< reco::TrackCollection > trackCollection_
const double w
Definition: UKUtility.cc:23
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::ESWatcher< HeavyIonRcd > hiWatcher
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< CaloTower >::const_iterator const_iterator
void setCaloCentRefBins(const int caloCentRefMinBin, const int caloCentRefMaxBin)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
cv
Definition: cuy.py:364
edm::EDGetTokenT< int > centralityBinToken
edm::EDGetTokenT< std::vector< reco::CastorTower > > castorToken
const double maxTransverse[]
double getEtScale(double vtx, int centbin) const
edm::ESWatcher< HeavyIonRPRcd > hirpWatcher
const std::string MomConsWeight[]
GenPlane * rp[NumEPNames]
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
HiEvtPlaneFlatten * flat[NumEPNames]
void addParticle(double w, double PtOrEt, double s, double c, double eta)
edm::EDGetTokenT< std::vector< reco::Vertex > > vertexToken
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< reco::TrackCollection > trackToken
edm::Handle< std::vector< reco::CastorTower > > castorCollection_
double getW(double pt, double vtx, int centbin) const
T const * product() const
Definition: Handle.h:69
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool accepted(std::vector< std::string_view > const &, std::string_view)
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:73
bool isPixel(HitType hitType)
static const int NumEPNames
edm::EDGetTokenT< CaloTowerCollection > caloToken
double getAngle(double &ang, double &sv, double &cv, double &svNoWgt, double &cvNoWgt, double &w, double &w2, double &PtOrEt, double &PtOrEt2, uint &epmult)
void reset(double vett[256])
Definition: TPedValues.cc:11
def move(src, dest)
Definition: eostools.py:511
edm::Handle< std::vector< reco::Vertex > > vertex_
const double minTransverse[]

Member Data Documentation

double EvtPlaneProducer::caloCentRef_
private

Definition at line 217 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::caloCentRefWidth_
private

Definition at line 218 of file EvtPlaneProducer.cc.

Referenced by produce().

edm::Handle<CaloTowerCollection> EvtPlaneProducer::caloCollection_
private

Definition at line 192 of file EvtPlaneProducer.cc.

Referenced by produce().

edm::InputTag EvtPlaneProducer::caloTag_
private

Definition at line 190 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

edm::EDGetTokenT<CaloTowerCollection> EvtPlaneProducer::caloToken
private

Definition at line 191 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().

edm::Handle<std::vector<reco::CastorTower> > EvtPlaneProducer::castorCollection_
private

Definition at line 196 of file EvtPlaneProducer.cc.

Referenced by produce().

edm::InputTag EvtPlaneProducer::castorTag_
private

Definition at line 194 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

edm::EDGetTokenT<std::vector<reco::CastorTower> > EvtPlaneProducer::castorToken
private

Definition at line 195 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().

int EvtPlaneProducer::CentBinCompression_
private

Definition at line 219 of file EvtPlaneProducer.cc.

Referenced by produce().

edm::InputTag EvtPlaneProducer::centralityBinTag_
private

Definition at line 183 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

edm::EDGetTokenT<int> EvtPlaneProducer::centralityBinToken
private

Definition at line 184 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().

std::string EvtPlaneProducer::centralityLabel_
private

Definition at line 180 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().

std::string EvtPlaneProducer::centralityMC_
private

Definition at line 181 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

std::string EvtPlaneProducer::centralityVariable_
private

Definition at line 179 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

double EvtPlaneProducer::chi2_
private

Definition at line 213 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::dzerr_
private

Definition at line 212 of file EvtPlaneProducer.cc.

Referenced by produce().

HiEvtPlaneFlatten* EvtPlaneProducer::flat[NumEPNames]
private

Definition at line 220 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), produce(), and ~EvtPlaneProducer().

int EvtPlaneProducer::FlatOrder_
private

Definition at line 214 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

edm::ESWatcher<HeavyIonRPRcd> EvtPlaneProducer::hirpWatcher
private

Definition at line 203 of file EvtPlaneProducer.cc.

Referenced by produce().

edm::ESWatcher<HeavyIonRcd> EvtPlaneProducer::hiWatcher
private

Definition at line 202 of file EvtPlaneProducer.cc.

Referenced by produce().

bool EvtPlaneProducer::loadDB_
private

Definition at line 205 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::maxet_
private

Definition at line 207 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::maxpt_
private

Definition at line 209 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::maxvtx_
private

Definition at line 211 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::minet_
private

Definition at line 206 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::minpt_
private

Definition at line 208 of file EvtPlaneProducer.cc.

Referenced by produce().

double EvtPlaneProducer::minvtx_
private

Definition at line 210 of file EvtPlaneProducer.cc.

double EvtPlaneProducer::nCentBins_
private

Definition at line 216 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().

int EvtPlaneProducer::NumFlatBins_
private

Definition at line 215 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

GenPlane* EvtPlaneProducer::rp[NumEPNames]
private

Definition at line 173 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().

edm::Handle<reco::TrackCollection> EvtPlaneProducer::trackCollection_
private

Definition at line 200 of file EvtPlaneProducer.cc.

Referenced by produce().

edm::InputTag EvtPlaneProducer::trackTag_
private

Definition at line 198 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

edm::EDGetTokenT<reco::TrackCollection> EvtPlaneProducer::trackToken
private

Definition at line 199 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().

edm::Handle<std::vector<reco::Vertex> > EvtPlaneProducer::vertex_
private

Definition at line 188 of file EvtPlaneProducer.cc.

Referenced by produce().

edm::InputTag EvtPlaneProducer::vertexTag_
private

Definition at line 186 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer().

edm::EDGetTokenT<std::vector<reco::Vertex> > EvtPlaneProducer::vertexToken
private

Definition at line 187 of file EvtPlaneProducer.cc.

Referenced by EvtPlaneProducer(), and produce().