CMS 3D CMS Logo

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

List of all members.

Public Member Functions

 EventTimeDistribution (const edm::ParameterSet &)
 ~EventTimeDistribution ()

Private Member Functions

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

Private Attributes

const edm::InputTag _apvphasecoll
TH1F ** _bx
TH1F ** _bxincycle
TH2F ** _bxincyclevsbx
TH1F ** _dbx
TH2F ** _dbxvsbx
TH2F ** _dbxvsbxincycle
const edm::InputTag _historyProduct
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

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 _apvphasecoll, _bx, _bxincycle, _bxincyclevsbx, _dbx, _dbxvsbx, _dbxvsbxincycle, _orbit, _orbitvsbxincycle, _rhm, _wantbxincyclevsbx, _wantdbxvsbx, _wantdbxvsbxincycle, _wantorbitvsbxincycle, edm::ParameterSet::getUntrackedParameter(), LogDebug, m_dbxhistos, m_dbxindices, m_ewhdepth, m_ewhdepthHisto, m_LSfrac, m_maxLS, RunHistogramManager::makeTH1F(), and RunHistogramManager::makeTH2F().

                                                                          :
  _historyProduct(iConfig.getParameter<edm::InputTag>("historyProduct")),
  _apvphasecoll(iConfig.getParameter<edm::InputTag>("apvPhaseCollection")),
  _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition","None")),
  _wantdbxvsbxincycle(iConfig.getUntrackedParameter<bool>("wantDBXvsBXincycle",false)),
  _wantdbxvsbx(iConfig.getUntrackedParameter<bool>("wantDBXvsBX",false)),
  _wantbxincyclevsbx(iConfig.getUntrackedParameter<bool>("wantBXincyclevsBX",false)),
  _wantorbitvsbxincycle(iConfig.getUntrackedParameter<bool>("wantOrbitvsBXincycle",false)),
  _nevents(0),
  m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin",100)),
  m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction",4)),
  m_ewhdepthHisto(iConfig.getUntrackedParameter<bool>("wantEWHDepthHisto",false)),
  _rhm(),
  _dbxvsbxincycle(0),   _dbxvsbx(0),   _bxincyclevsbx(0),   _orbitvsbxincycle(0), m_ewhdepth(0)
{
   //now do what ever initialization is needed

  std::vector<edm::ParameterSet> dbxhistoparams(iConfig.getUntrackedParameter<std::vector<edm::ParameterSet> >("dbxHistosParams",std::vector<edm::ParameterSet>()));

  for(std::vector<edm::ParameterSet>::const_iterator params=dbxhistoparams.begin();params!=dbxhistoparams.end();++params) {
    m_dbxindices.push_back(std::pair<unsigned int,unsigned int>(params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent")));
    char hname[300];
    sprintf(hname,"dbx_%d_%d",params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent"));
    char htitle[300];
    sprintf(htitle,"dbx(%d,%d)",params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent"));

    m_dbxhistos.push_back(_rhm.makeTH1F(hname,htitle,params->getParameter<int>("nbins"),params->getParameter<double>("min"),
                                        params->getParameter<double>("max")));
    LogDebug("DBXHistoPreBooking") << "Booked DBX histo named " << hname << " untitled " << htitle;
  }



  _dbx = _rhm.makeTH1F("dbx","dbx",1000,-0.5,999.5);
  _bx = _rhm.makeTH1F("bx","BX number",3564,-0.5,3563.5);
  _bxincycle = _rhm.makeTH1F("bxcycle","bxcycle",70,-0.5,69.5);
  _orbit = _rhm.makeTH1F("orbit","orbit",m_LSfrac*m_maxLS,0,m_maxLS*262144);
  if(_wantdbxvsbxincycle) _dbxvsbxincycle = _rhm.makeTH2F("dbxvsbxincycle","dbxvsbxincycle",70,-0.5,69.5,1000,-0.5,999.5);
  if(_wantdbxvsbx) _dbxvsbx = _rhm.makeTH2F("dbxvsbx","dbxvsbx",3564,-0.5,3563.5,1000,-0.5,999.5);
  if(_wantbxincyclevsbx) _bxincyclevsbx = _rhm.makeTH2F("bxincyclevsbx","bxincyclevsbx",3564,-0.5,3563.5,70,-0.5,69.5);
  if(_wantorbitvsbxincycle) _orbitvsbxincycle = _rhm.makeTH2F("orbitvsbxincycle","orbitvsbxincycle",70,-0.5,69.5,m_maxLS,0,m_maxLS*262144);
  if(m_ewhdepthHisto) m_ewhdepth = _rhm.makeTH1F("ewhdepth","EventWithHistory Depth",11,-0.5,10.5);

  edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecoll << " used";

}
EventTimeDistribution::~EventTimeDistribution ( )

Definition at line 158 of file EventTimeDistribution.cc.

{
 
   // do anything here that needs to be done at desctruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

void EventTimeDistribution::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [override, private, virtual]

Implements edm::EDAnalyzer.

Definition at line 173 of file EventTimeDistribution.cc.

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

{
   using namespace edm;

   _nevents++;

   edm::Handle<EventWithHistory> he;
   iEvent.getByLabel(_historyProduct,he);

   // improve the matchin between default and actual partitions
   
   (*_dbx)->Fill(he->deltaBX());
   std::vector<std::pair<unsigned int,unsigned int> >::const_iterator indices=m_dbxindices.begin();
   for(std::vector<TH1F**>::const_iterator dbxhist=m_dbxhistos.begin();dbxhist!=m_dbxhistos.end();++dbxhist,++indices) {
     (*(*dbxhist))->Fill(he->deltaBX(indices->first,indices->second));
   }

   (*_bx)->Fill(iEvent.bunchCrossing());
   (*_orbit)->Fill(iEvent.orbitNumber());
   if(_dbxvsbx && *_dbxvsbx) (*_dbxvsbx)->Fill(iEvent.bunchCrossing(),he->deltaBX());
   if(m_ewhdepth && *m_ewhdepth) (*m_ewhdepth)->Fill(he->depth());

   edm::Handle<APVCyclePhaseCollection> apvphase;
   iEvent.getByLabel(_apvphasecoll,apvphase);

   long long tbx = he->absoluteBX();
   if(apvphase.isValid() && !apvphase.failedToGet()) {
     const int thephase = apvphase->getPhase(_phasepart); 
     if(thephase!=APVCyclePhaseCollection::invalid &&
        thephase!=APVCyclePhaseCollection::multiphase &&
        thephase!=APVCyclePhaseCollection::nopartition) {

       tbx -= thephase;
       (*_bxincycle)->Fill(tbx%70);
       if(_dbxvsbxincycle && *_dbxvsbxincycle) (*_dbxvsbxincycle)->Fill(tbx%70,he->deltaBX());
       if(_bxincyclevsbx && *_bxincyclevsbx) (*_bxincyclevsbx)->Fill(iEvent.bunchCrossing(),tbx%70);
       if(_orbitvsbxincycle && *_orbitvsbxincycle) (*_orbitvsbxincycle)->Fill(tbx%70,iEvent.orbitNumber());

     }
     else {
       LogDebug("InvalidPhase") << "Invalid APVCyclePhase value : " << _phasepart << " " << thephase;
     }
   }
}
void EventTimeDistribution::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 273 of file EventTimeDistribution.cc.

{

}
void EventTimeDistribution::beginRun ( const edm::Run iRun,
const edm::EventSetup  
) [override, private, virtual]

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.

{

  _rhm.beginRun(iRun);

  if(*_dbx) {    (*_dbx)->GetXaxis()->SetTitle("#DeltaBX"); }

  LogDebug("NomberOfHistos") << m_dbxhistos.size();
  for(std::vector<TH1F**>::const_iterator dbxhist=m_dbxhistos.begin();dbxhist!=m_dbxhistos.end();++dbxhist) {
    LogDebug("HistoPointer") << *dbxhist;
    if(*(*dbxhist)) { (*(*dbxhist))->GetXaxis()->SetTitle("#DeltaBX"); }
  }
  LogDebug("LabelDone") << "all labels set";

  if(*_bx) { (*_bx)->GetXaxis()->SetTitle("BX");  }

  if(*_bxincycle) {  (*_bxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); }

  if(*_orbit) {
    (*_orbit)->SetBit(TH1::kCanRebin);
    (*_orbit)->GetXaxis()->SetTitle("time [Orb#]"); 
  }

  LogDebug("StdPlotsDone") << "all labels in std plots set";

  if(_dbxvsbxincycle && *_dbxvsbxincycle) {
    (*_dbxvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_dbxvsbxincycle)->GetYaxis()->SetTitle("#DeltaBX"); 
  }

  if(_dbxvsbx && *_dbxvsbx) { (*_dbxvsbx)->GetXaxis()->SetTitle("BX"); (*_dbxvsbx)->GetYaxis()->SetTitle("#DeltaBX"); }

  if(_bxincyclevsbx && *_bxincyclevsbx) {
    (*_bxincyclevsbx)->GetXaxis()->SetTitle("BX"); (*_bxincyclevsbx)->GetYaxis()->SetTitle("Event BX mod(70)");
  }

  if(_orbitvsbxincycle && *_orbitvsbxincycle) {
    (*_orbitvsbxincycle)->SetBit(TH1::kCanRebin);
    (*_orbitvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_orbitvsbxincycle)->GetYaxis()->SetTitle("time [Orb#]"); 
  }

  if(m_ewhdepth && *m_ewhdepth) {
    (*m_ewhdepth)->GetXaxis()->SetTitle("Depth");
  }

}
void EventTimeDistribution::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 280 of file EventTimeDistribution.cc.

References _nevents.

                              {

  edm::LogInfo("EndOfJob") << _nevents << " analyzed events";

}
void EventTimeDistribution::endRun ( const edm::Run iRun,
const edm::EventSetup  
) [override, private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 266 of file EventTimeDistribution.cc.

{
}

Member Data Documentation

Definition at line 70 of file EventTimeDistribution.cc.

Referenced by analyze(), and EventTimeDistribution().

TH1F** EventTimeDistribution::_bx [private]

Definition at line 88 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

Definition at line 89 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

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().

Definition at line 92 of file EventTimeDistribution.cc.

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

Definition at line 91 of file EventTimeDistribution.cc.

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

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().

Definition at line 90 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

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().

Definition at line 83 of file EventTimeDistribution.cc.

Referenced by beginRun(), and EventTimeDistribution().

Definition at line 74 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

Definition at line 73 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

Definition at line 72 of file EventTimeDistribution.cc.

Referenced by EventTimeDistribution().

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().

Definition at line 95 of file EventTimeDistribution.cc.

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

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().