CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
EventTimeDistribution Class Reference

#include <DPGAnalysis/SiStripTools/plugins/EventTimeDistribution.cc>

Inheritance diagram for EventTimeDistribution:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 EventTimeDistribution (const edm::ParameterSet &)
 
 ~EventTimeDistribution ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
virtual void beginJob () override
 
virtual void beginRun (const edm::Run &, const edm::EventSetup &) override
 
virtual void endJob () override
 
virtual void endRun (const edm::Run &, const edm::EventSetup &) override
 

Private Attributes

edm::EDGetTokenT
< APVCyclePhaseCollection
_apvphasecollToken
 
TH1F ** _bx
 
TH1F ** _bxincycle
 
TH2F ** _bxincyclevsbx
 
TH1F ** _dbx
 
TH2F ** _dbxvsbx
 
TH2F ** _dbxvsbxincycle
 
edm::EDGetTokenT
< EventWithHistory
_historyProductToken
 
unsigned int _nevents
 
TH1F ** _orbit
 
TH2F ** _orbitvsbxincycle
 
const std::string _phasepart
 
RunHistogramManager _rhm
 
const bool _wantbxincyclevsbx
 
const bool _wantdbxvsbx
 
const bool _wantdbxvsbxincycle
 
const bool _wantorbitvsbxincycle
 
std::vector< TH1F ** > m_dbxhistos
 
std::vector< std::pair
< unsigned int, unsigned int > > 
m_dbxindices
 
TH1F ** m_ewhdepth
 
const bool m_ewhdepthHisto
 
const unsigned int m_LSfrac
 
const unsigned int m_maxLS
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

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

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

Definition at line 54 of file EventTimeDistribution.cc.

Constructor & Destructor Documentation

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

Definition at line 110 of file EventTimeDistribution.cc.

References _bx, _bxincycle, _bxincyclevsbx, _dbx, _dbxvsbx, _dbxvsbxincycle, _orbit, _orbitvsbxincycle, _rhm, _wantbxincyclevsbx, _wantdbxvsbx, _wantdbxvsbxincycle, _wantorbitvsbxincycle, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), LogDebug, m_dbxhistos, m_dbxindices, m_ewhdepth, m_ewhdepthHisto, m_LSfrac, m_maxLS, RunHistogramManager::makeTH1F(), and RunHistogramManager::makeTH2F().

110  :
111  _historyProductToken(consumes<EventWithHistory>(iConfig.getParameter<edm::InputTag>("historyProduct"))),
112  _apvphasecollToken(consumes<APVCyclePhaseCollection>(iConfig.getParameter<edm::InputTag>("apvPhaseCollection"))),
113  _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition","None")),
114  _wantdbxvsbxincycle(iConfig.getUntrackedParameter<bool>("wantDBXvsBXincycle",false)),
115  _wantdbxvsbx(iConfig.getUntrackedParameter<bool>("wantDBXvsBX",false)),
116  _wantbxincyclevsbx(iConfig.getUntrackedParameter<bool>("wantBXincyclevsBX",false)),
117  _wantorbitvsbxincycle(iConfig.getUntrackedParameter<bool>("wantOrbitvsBXincycle",false)),
118  _nevents(0),
119  m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin",100)),
120  m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction",4)),
121  m_ewhdepthHisto(iConfig.getUntrackedParameter<bool>("wantEWHDepthHisto",false)),
124 {
125  //now do what ever initialization is needed
126 
127  std::vector<edm::ParameterSet> dbxhistoparams(iConfig.getUntrackedParameter<std::vector<edm::ParameterSet> >("dbxHistosParams",std::vector<edm::ParameterSet>()));
128 
129  for(std::vector<edm::ParameterSet>::const_iterator params=dbxhistoparams.begin();params!=dbxhistoparams.end();++params) {
130  m_dbxindices.push_back(std::pair<unsigned int,unsigned int>(params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent")));
131  char hname[300];
132  sprintf(hname,"dbx_%d_%d",params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent"));
133  char htitle[300];
134  sprintf(htitle,"dbx(%d,%d)",params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent"));
135 
136  m_dbxhistos.push_back(_rhm.makeTH1F(hname,htitle,params->getParameter<int>("nbins"),params->getParameter<double>("min"),
137  params->getParameter<double>("max")));
138  LogDebug("DBXHistoPreBooking") << "Booked DBX histo named " << hname << " untitled " << htitle;
139  }
140 
141 
142 
143  _dbx = _rhm.makeTH1F("dbx","dbx",1000,-0.5,999.5);
144  _bx = _rhm.makeTH1F("bx","BX number",3564,-0.5,3563.5);
145  _bxincycle = _rhm.makeTH1F("bxcycle","bxcycle",70,-0.5,69.5);
146  _orbit = _rhm.makeTH1F("orbit","orbit",m_LSfrac*m_maxLS,0,m_maxLS*262144);
147  if(_wantdbxvsbxincycle) _dbxvsbxincycle = _rhm.makeTH2F("dbxvsbxincycle","dbxvsbxincycle",70,-0.5,69.5,1000,-0.5,999.5);
148  if(_wantdbxvsbx) _dbxvsbx = _rhm.makeTH2F("dbxvsbx","dbxvsbx",3564,-0.5,3563.5,1000,-0.5,999.5);
149  if(_wantbxincyclevsbx) _bxincyclevsbx = _rhm.makeTH2F("bxincyclevsbx","bxincyclevsbx",3564,-0.5,3563.5,70,-0.5,69.5);
150  if(_wantorbitvsbxincycle) _orbitvsbxincycle = _rhm.makeTH2F("orbitvsbxincycle","orbitvsbxincycle",70,-0.5,69.5,m_maxLS,0,m_maxLS*262144);
151  if(m_ewhdepthHisto) m_ewhdepth = _rhm.makeTH1F("ewhdepth","EventWithHistory Depth",11,-0.5,10.5);
152 
153  edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << iConfig.getParameter<edm::InputTag>("apvPhaseCollection") << " used";
154 
155 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< APVCyclePhaseCollection > _apvphasecollToken
std::vector< std::pair< unsigned int, unsigned int > > m_dbxindices
edm::EDGetTokenT< EventWithHistory > _historyProductToken
TH1F ** makeTH1F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
const unsigned int m_LSfrac
RunHistogramManager _rhm
const std::string _phasepart
const unsigned int m_maxLS
TH2F ** makeTH2F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax, const unsigned int nbiny, const double ymin, const double ymax)
std::vector< TH1F ** > m_dbxhistos
EventTimeDistribution::~EventTimeDistribution ( )

Definition at line 158 of file EventTimeDistribution.cc.

159 {
160 
161  // do anything here that needs to be done at desctruction time
162  // (e.g. close files, deallocate resources etc.)
163 
164 }

Member Function Documentation

void EventTimeDistribution::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 173 of file EventTimeDistribution.cc.

References _apvphasecollToken, _bxincyclevsbx, _dbxvsbx, _dbxvsbxincycle, _historyProductToken, _nevents, _orbitvsbxincycle, _phasepart, edm::EventBase::bunchCrossing(), HcalObjRepresent::Fill(), edm::Event::getByToken(), APVCyclePhaseCollection::invalid, LogDebug, m_dbxhistos, m_dbxindices, m_ewhdepth, APVCyclePhaseCollection::multiphase, APVCyclePhaseCollection::nopartition, and edm::EventBase::orbitNumber().

174 {
175  using namespace edm;
176 
177  _nevents++;
178 
180  iEvent.getByToken(_historyProductToken,he);
181 
182  // improve the matchin between default and actual partitions
183 
184  (*_dbx)->Fill(he->deltaBX());
185  std::vector<std::pair<unsigned int,unsigned int> >::const_iterator indices=m_dbxindices.begin();
186  for(std::vector<TH1F**>::const_iterator dbxhist=m_dbxhistos.begin();dbxhist!=m_dbxhistos.end();++dbxhist,++indices) {
187  (*(*dbxhist))->Fill(he->deltaBX(indices->first,indices->second));
188  }
189 
190  (*_bx)->Fill(iEvent.bunchCrossing());
191  (*_orbit)->Fill(iEvent.orbitNumber());
192  if(_dbxvsbx && *_dbxvsbx) (*_dbxvsbx)->Fill(iEvent.bunchCrossing(),he->deltaBX());
193  if(m_ewhdepth && *m_ewhdepth) (*m_ewhdepth)->Fill(he->depth());
194 
196  iEvent.getByToken(_apvphasecollToken,apvphase);
197 
198  long long tbx = he->absoluteBX();
199  if(apvphase.isValid() && !apvphase.failedToGet()) {
200  const int thephase = apvphase->getPhase(_phasepart);
201  if(thephase!=APVCyclePhaseCollection::invalid &&
204 
205  tbx -= thephase;
206  (*_bxincycle)->Fill(tbx%70);
207  if(_dbxvsbxincycle && *_dbxvsbxincycle) (*_dbxvsbxincycle)->Fill(tbx%70,he->deltaBX());
208  if(_bxincyclevsbx && *_bxincyclevsbx) (*_bxincyclevsbx)->Fill(iEvent.bunchCrossing(),tbx%70);
209  if(_orbitvsbxincycle && *_orbitvsbxincycle) (*_orbitvsbxincycle)->Fill(tbx%70,iEvent.orbitNumber());
210 
211  }
212  else {
213  LogDebug("InvalidPhase") << "Invalid APVCyclePhase value : " << _phasepart << " " << thephase;
214  }
215  }
216 }
#define LogDebug(id)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
int bunchCrossing() const
Definition: EventBase.h:62
edm::EDGetTokenT< APVCyclePhaseCollection > _apvphasecollToken
std::vector< std::pair< unsigned int, unsigned int > > m_dbxindices
edm::EDGetTokenT< EventWithHistory > _historyProductToken
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int orbitNumber() const
Definition: EventBase.h:63
const std::string _phasepart
std::vector< TH1F ** > m_dbxhistos
void EventTimeDistribution::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 273 of file EventTimeDistribution.cc.

274 {
275 
276 }
void EventTimeDistribution::beginRun ( const edm::Run iRun,
const edm::EventSetup  
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 219 of file EventTimeDistribution.cc.

References _bx, _bxincycle, _bxincyclevsbx, _dbx, _dbxvsbx, _dbxvsbxincycle, _orbit, _orbitvsbxincycle, _rhm, RunHistogramManager::beginRun(), LogDebug, m_dbxhistos, and m_ewhdepth.

220 {
221 
222  _rhm.beginRun(iRun);
223 
224  if(*_dbx) { (*_dbx)->GetXaxis()->SetTitle("#DeltaBX"); }
225 
226  LogDebug("NomberOfHistos") << m_dbxhistos.size();
227  for(std::vector<TH1F**>::const_iterator dbxhist=m_dbxhistos.begin();dbxhist!=m_dbxhistos.end();++dbxhist) {
228  LogDebug("HistoPointer") << *dbxhist;
229  if(*(*dbxhist)) { (*(*dbxhist))->GetXaxis()->SetTitle("#DeltaBX"); }
230  }
231  LogDebug("LabelDone") << "all labels set";
232 
233  if(*_bx) { (*_bx)->GetXaxis()->SetTitle("BX"); }
234 
235  if(*_bxincycle) { (*_bxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); }
236 
237  if(*_orbit) {
238  (*_orbit)->SetBit(TH1::kCanRebin);
239  (*_orbit)->GetXaxis()->SetTitle("time [Orb#]");
240  }
241 
242  LogDebug("StdPlotsDone") << "all labels in std plots set";
243 
245  (*_dbxvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_dbxvsbxincycle)->GetYaxis()->SetTitle("#DeltaBX");
246  }
247 
248  if(_dbxvsbx && *_dbxvsbx) { (*_dbxvsbx)->GetXaxis()->SetTitle("BX"); (*_dbxvsbx)->GetYaxis()->SetTitle("#DeltaBX"); }
249 
251  (*_bxincyclevsbx)->GetXaxis()->SetTitle("BX"); (*_bxincyclevsbx)->GetYaxis()->SetTitle("Event BX mod(70)");
252  }
253 
255  (*_orbitvsbxincycle)->SetBit(TH1::kCanRebin);
256  (*_orbitvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_orbitvsbxincycle)->GetYaxis()->SetTitle("time [Orb#]");
257  }
258 
259  if(m_ewhdepth && *m_ewhdepth) {
260  (*m_ewhdepth)->GetXaxis()->SetTitle("Depth");
261  }
262 
263 }
#define LogDebug(id)
RunHistogramManager _rhm
void beginRun(const edm::Run &iRun)
std::vector< TH1F ** > m_dbxhistos
void EventTimeDistribution::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 280 of file EventTimeDistribution.cc.

References _nevents.

280  {
281 
282  edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
283 
284 }
void EventTimeDistribution::endRun ( const edm::Run iRun,
const edm::EventSetup  
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 266 of file EventTimeDistribution.cc.

267 {
268 }

Member Data Documentation

edm::EDGetTokenT<APVCyclePhaseCollection> EventTimeDistribution::_apvphasecollToken
private

Definition at line 70 of file EventTimeDistribution.cc.

Referenced by analyze().

TH1F** EventTimeDistribution::_bx
private

Definition at line 88 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

TH1F** EventTimeDistribution::_bxincycle
private

Definition at line 89 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

TH2F** EventTimeDistribution::_bxincyclevsbx
private

Definition at line 93 of file EventTimeDistribution.cc.

Referenced by analyze(), beginRun(), and EventTimeDistribution().

TH1F** EventTimeDistribution::_dbx
private

Definition at line 85 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

TH2F** EventTimeDistribution::_dbxvsbx
private

Definition at line 92 of file EventTimeDistribution.cc.

Referenced by analyze(), beginRun(), and EventTimeDistribution().

TH2F** EventTimeDistribution::_dbxvsbxincycle
private

Definition at line 91 of file EventTimeDistribution.cc.

Referenced by analyze(), beginRun(), and EventTimeDistribution().

edm::EDGetTokenT<EventWithHistory> EventTimeDistribution::_historyProductToken
private

Definition at line 69 of file EventTimeDistribution.cc.

Referenced by analyze().

unsigned int EventTimeDistribution::_nevents
private

Definition at line 76 of file EventTimeDistribution.cc.

Referenced by analyze(), and endJob().

TH1F** EventTimeDistribution::_orbit
private

Definition at line 90 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

TH2F** EventTimeDistribution::_orbitvsbxincycle
private

Definition at line 94 of file EventTimeDistribution.cc.

Referenced by analyze(), beginRun(), and EventTimeDistribution().

const std::string EventTimeDistribution::_phasepart
private

Definition at line 71 of file EventTimeDistribution.cc.

Referenced by analyze().

RunHistogramManager EventTimeDistribution::_rhm
private

Definition at line 83 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

const bool EventTimeDistribution::_wantbxincyclevsbx
private

Definition at line 74 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

const bool EventTimeDistribution::_wantdbxvsbx
private

Definition at line 73 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

const bool EventTimeDistribution::_wantdbxvsbxincycle
private

Definition at line 72 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

const bool EventTimeDistribution::_wantorbitvsbxincycle
private

Definition at line 75 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

std::vector<TH1F**> EventTimeDistribution::m_dbxhistos
private

Definition at line 86 of file EventTimeDistribution.cc.

Referenced by analyze(), beginRun(), and EventTimeDistribution().

std::vector<std::pair<unsigned int,unsigned int> > EventTimeDistribution::m_dbxindices
private

Definition at line 87 of file EventTimeDistribution.cc.

Referenced by analyze(), and EventTimeDistribution().

TH1F** EventTimeDistribution::m_ewhdepth
private

Definition at line 95 of file EventTimeDistribution.cc.

Referenced by analyze(), beginRun(), and EventTimeDistribution().

const bool EventTimeDistribution::m_ewhdepthHisto
private

Definition at line 79 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

const unsigned int EventTimeDistribution::m_LSfrac
private

Definition at line 78 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

const unsigned int EventTimeDistribution::m_maxLS
private

Definition at line 77 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().