CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/DPGAnalysis/SiStripTools/plugins/APVCyclePhaseProducerFromL1TS.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripTools
00004 // Class:      APVCyclePhaseProducerFromL1TS
00005 // 
00013 //
00014 // Original Author:  Andrea Venturi
00015 //         Created:  Mon Jan 12 09:05:45 CET 2009
00016 //
00017 //
00018 
00019 
00020 // system include files
00021 #include <memory>
00022 
00023 // user include files
00024 #include "FWCore/Framework/interface/Frameworkfwd.h"
00025 #include "FWCore/Framework/interface/EDProducer.h"
00026 
00027 #include "FWCore/Framework/interface/Event.h"
00028 #include "FWCore/Framework/interface/Run.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030 
00031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00032 
00033 #include "FWCore/Utilities/interface/InputTag.h"
00034 
00035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00036 #include "FWCore/Utilities/interface/Exception.h"
00037 
00038 #include <map>
00039 #include <vector>
00040 #include <utility>
00041 #include <string>
00042 
00043 #include "TH1F.h"
00044 #include "TProfile.h"
00045 
00046 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
00047 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00048 
00049 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
00050 
00051 //
00052 // class decleration
00053 //
00054 
00055 class APVCyclePhaseProducerFromL1TS : public edm::EDProducer {
00056    public:
00057       explicit APVCyclePhaseProducerFromL1TS(const edm::ParameterSet&);
00058       ~APVCyclePhaseProducerFromL1TS();
00059 
00060 private:
00061   virtual void beginJob() ;
00062   virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
00063   virtual void produce(edm::Event&, const edm::EventSetup&) override;
00064   virtual void endJob() ;
00065 
00066   bool isBadRun(const unsigned int) const;
00067   
00068       // ----------member data ---------------------------
00069 
00070   const edm::InputTag _l1tscollection;
00071   const std::vector<std::string> _defpartnames;
00072   const std::vector<int> _defphases;
00073   const bool _wantHistos;
00074   const bool _useEC0;
00075   const int _magicOffset;
00076   const unsigned int m_maxLS;
00077   const unsigned int m_LSfrac;
00078 
00079   RunHistogramManager m_rhm;
00080 
00081   TH1F** _hsize;
00082   TH1F** _hlresync;
00083   TH1F** _hlOC0;
00084   TH1F** _hlTE;
00085   TH1F** _hlstart;
00086   TH1F** _hlEC0;
00087   TH1F** _hlHR;
00088 
00089   TH1F** _hdlec0lresync;
00090   TH1F** _hdlresynclHR;
00091 
00092   std::vector<std::pair<unsigned int, unsigned int> > m_badruns;
00093   
00094   long long _lastResync;
00095   long long _lastHardReset;
00096   long long _lastStart;
00097   long long _lastEventCounter0;
00098   long long _lastOrbitCounter0;
00099   long long _lastTestEnable;
00100 
00101 
00102 };
00103 
00104 //
00105 // constants, enums and typedefs
00106 //
00107 
00108 
00109 //
00110 // static data member definitions
00111 //
00112 
00113 //
00114 // constructors and destructor
00115 //
00116 APVCyclePhaseProducerFromL1TS::APVCyclePhaseProducerFromL1TS(const edm::ParameterSet& iConfig):
00117   _l1tscollection(iConfig.getParameter<edm::InputTag>("l1TSCollection")),
00118   _defpartnames(iConfig.getParameter<std::vector<std::string> >("defaultPartitionNames")),
00119   _defphases(iConfig.getParameter<std::vector<int> >("defaultPhases")),
00120   _wantHistos(iConfig.getUntrackedParameter<bool>("wantHistos",false)),
00121   _useEC0(iConfig.getUntrackedParameter<bool>("useEC0",false)),
00122   _magicOffset(iConfig.getUntrackedParameter<int>("magicOffset",8)),
00123   m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin",250)),
00124   m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction",16)),
00125   m_rhm(),
00126   _hsize(0),_hlresync(0),_hlOC0(0),_hlTE(0),_hlstart(0),_hlEC0(0),_hlHR(0),_hdlec0lresync(0),_hdlresynclHR(0),
00127   m_badruns(),
00128   _lastResync(-1),_lastHardReset(-1),_lastStart(-1),
00129   _lastEventCounter0(-1),_lastOrbitCounter0(-1),_lastTestEnable(-1)
00130 {
00131 
00132   produces<APVCyclePhaseCollection,edm::InEvent>();
00133 
00134   m_badruns.push_back(std::pair<unsigned int, unsigned int>(0,131767));
00135   m_badruns.push_back(std::pair<unsigned int, unsigned int>(193150,193733));
00136 
00137    //now do what ever other initialization is needed
00138 
00139   if(_wantHistos) {
00140     _hsize = m_rhm.makeTH1F("size","Level1TriggerScalers Collection size",20,-0.5,19.5);
00141 
00142     _hlresync = m_rhm.makeTH1F("lresync","Orbit of last resync",m_LSfrac*m_maxLS,0,m_maxLS*262144);
00143     _hlOC0 = m_rhm.makeTH1F("lOC0","Orbit of last OC0",m_LSfrac*m_maxLS,0,m_maxLS*262144);
00144     _hlTE = m_rhm.makeTH1F("lTE","Orbit of last TestEnable",m_LSfrac*m_maxLS,0,m_maxLS*262144);
00145     _hlstart = m_rhm.makeTH1F("lstart","Orbit of last Start",m_LSfrac*m_maxLS,0,m_maxLS*262144);
00146     _hlEC0 = m_rhm.makeTH1F("lEC0","Orbit of last EC0",m_LSfrac*m_maxLS,0,m_maxLS*262144);
00147     _hlHR = m_rhm.makeTH1F("lHR","Orbit of last HardReset",m_LSfrac*m_maxLS,0,m_maxLS*262144);
00148     _hdlec0lresync = m_rhm.makeTH1F("dlec0lresync","Orbit difference EC0-Resync",4000,-1999.5,2000.5);
00149     _hdlresynclHR = m_rhm.makeTH1F("dlresynclHR","Orbit difference Resync-HR",4000,-1999.5,2000.5);
00150 
00151   }
00152 
00153 
00154 }
00155 
00156 
00157 APVCyclePhaseProducerFromL1TS::~APVCyclePhaseProducerFromL1TS()
00158 {
00159  
00160    // do anything here that needs to be done at desctruction time
00161    // (e.g. close files, deallocate resources etc.)
00162 
00163 }
00164 
00165 
00166 //
00167 // member functions
00168 //
00169 
00170 // ------------ method called to produce the data  ------------
00171 void
00172 APVCyclePhaseProducerFromL1TS::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) 
00173 
00174 {
00175 
00176   // reset offset vector
00177 
00178   if(_wantHistos) {
00179 
00180     m_rhm.beginRun(iRun);
00181 
00182     if(_hlresync && *_hlresync) {
00183       (*_hlresync)->GetXaxis()->SetTitle("Orbit");     (*_hlresync)->GetYaxis()->SetTitle("Events");
00184       (*_hlresync)->SetBit(TH1::kCanRebin);
00185     }
00186 
00187     if(_hlOC0 && *_hlOC0) {
00188       (*_hlOC0)->GetXaxis()->SetTitle("Orbit");     (*_hlOC0)->GetYaxis()->SetTitle("Events");
00189       (*_hlOC0)->SetBit(TH1::kCanRebin);
00190     }
00191 
00192     if(_hlTE && *_hlTE) {
00193       (*_hlTE)->GetXaxis()->SetTitle("Orbit");     (*_hlTE)->GetYaxis()->SetTitle("Events");
00194       (*_hlTE)->SetBit(TH1::kCanRebin);
00195     }
00196 
00197     if(_hlstart && *_hlstart) {
00198       (*_hlstart)->GetXaxis()->SetTitle("Orbit");     (*_hlstart)->GetYaxis()->SetTitle("Events");
00199       (*_hlstart)->SetBit(TH1::kCanRebin);
00200     }
00201 
00202     if(_hlEC0 && *_hlEC0) {
00203       (*_hlEC0)->GetXaxis()->SetTitle("Orbit");     (*_hlEC0)->GetYaxis()->SetTitle("Events");
00204       (*_hlEC0)->SetBit(TH1::kCanRebin);
00205     }
00206 
00207     if(_hlHR && *_hlHR) {
00208       (*_hlHR)->GetXaxis()->SetTitle("Orbit");     (*_hlHR)->GetYaxis()->SetTitle("Events");
00209       (*_hlHR)->SetBit(TH1::kCanRebin);
00210     }
00211     
00212     if(_hdlec0lresync && *_hdlec0lresync) {
00213       (*_hdlec0lresync)->GetXaxis()->SetTitle("lastEC0-lastResync"); 
00214     }
00215 
00216     if(_hdlresynclHR && *_hdlresynclHR) {
00217       (*_hdlresynclHR)->GetXaxis()->SetTitle("lastEC0-lastResync"); 
00218     }
00219 
00220   }
00221 
00222   if(isBadRun(iRun.run())) {
00223     LogDebug("UnreliableMissingL1TriggerScalers") << 
00224       "In this run L1TriggerScalers is missing or unreliable for phase determination: invlid phase will be returned"; 
00225   }
00226 
00227 }
00228 
00229 
00230 void
00231 APVCyclePhaseProducerFromL1TS::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00232 
00233   using namespace edm;
00234   
00235   std::auto_ptr<APVCyclePhaseCollection> apvphases(new APVCyclePhaseCollection() );
00236   
00237 
00238   std::vector<int> phases(_defphases.size(),APVCyclePhaseCollection::invalid);
00239 
00240   const std::vector<std::string>& partnames = _defpartnames;
00241 
00242   int phasechange = 0;
00243 
00244     
00245   Handle<Level1TriggerScalersCollection> l1ts;
00246   iEvent.getByLabel(_l1tscollection,l1ts);
00247   
00248   if(_wantHistos && _hsize && *_hsize) (*_hsize)->Fill(l1ts->size());
00249   
00250   // offset computation
00251   
00252   long long orbitoffset = 0;
00253   
00254   if(l1ts->size()>0) {
00255     
00256     if((*l1ts)[0].lastResync()!=0) {
00257       orbitoffset = _useEC0 ? (*l1ts)[0].lastEventCounter0() + _magicOffset : (*l1ts)[0].lastResync() + _magicOffset;
00258     }
00259     
00260     if(_wantHistos) {
00261       if(_hlresync && *_hlresync) (*_hlresync)->Fill((*l1ts)[0].lastResync());
00262       if(_hlOC0 && *_hlOC0) (*_hlOC0)->Fill((*l1ts)[0].lastOrbitCounter0());
00263       if(_hlTE && *_hlTE) (*_hlTE)->Fill((*l1ts)[0].lastTestEnable());
00264       if(_hlstart && *_hlstart) (*_hlstart)->Fill((*l1ts)[0].lastStart());
00265       if(_hlEC0 && *_hlEC0) (*_hlEC0)->Fill((*l1ts)[0].lastEventCounter0());
00266       if(_hlHR && *_hlHR) (*_hlHR)->Fill((*l1ts)[0].lastHardReset());
00267     }
00268     
00269     if(_lastResync != (*l1ts)[0].lastResync()) {
00270       _lastResync = (*l1ts)[0].lastResync();
00271       if(_wantHistos && _hdlec0lresync && *_hdlec0lresync) (*_hdlec0lresync)->Fill((*l1ts)[0].lastEventCounter0()-(*l1ts)[0].lastResync());
00272       LogDebug("TTCSignalReceived") << "New Resync at orbit " << _lastResync ;
00273     }
00274     if(_lastHardReset != (*l1ts)[0].lastHardReset()) {
00275       _lastHardReset = (*l1ts)[0].lastHardReset();
00276       if(_wantHistos && _hdlresynclHR && *_hdlresynclHR) (*_hdlresynclHR)->Fill((*l1ts)[0].lastResync()-(*l1ts)[0].lastHardReset());
00277       LogDebug("TTCSignalReceived") << "New HardReset at orbit " << _lastHardReset ;
00278     }
00279     if(_lastTestEnable != (*l1ts)[0].lastTestEnable()) {
00280       _lastTestEnable = (*l1ts)[0].lastTestEnable();
00281       //      LogDebug("TTCSignalReceived") << "New TestEnable at orbit " << _lastTestEnable ;
00282     }
00283     if(_lastOrbitCounter0 != (*l1ts)[0].lastOrbitCounter0()) {
00284       _lastOrbitCounter0 = (*l1ts)[0].lastOrbitCounter0();
00285       LogDebug("TTCSignalReceived") << "New OrbitCounter0 at orbit " << _lastOrbitCounter0 ;
00286     }
00287     if(_lastEventCounter0 != (*l1ts)[0].lastEventCounter0()) {
00288       _lastEventCounter0 = (*l1ts)[0].lastEventCounter0();
00289       LogDebug("TTCSignalReceived") << "New EventCounter0 at orbit " << _lastEventCounter0 ;
00290     }
00291     if(_lastStart != (*l1ts)[0].lastStart()) {
00292       _lastStart = (*l1ts)[0].lastStart();
00293       LogDebug("TTCSignalReceived") << "New Start at orbit " << _lastStart ;
00294     }
00295     
00296     if(!isBadRun(iEvent.run())) {
00297       phasechange = ((long long)(orbitoffset*3564))%70;
00298       
00299       for(unsigned int ipart=0;ipart<phases.size();++ipart) {
00300         phases[ipart] = (_defphases[ipart]+phasechange)%70;
00301       }
00302       
00303     }
00304   }
00305   
00306 
00307   if(phases.size() < partnames.size() ) {
00308     // throw exception
00309     throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent phases/partitions vector sizes: " 
00310                                              << phases.size() << " " 
00311                                              << partnames.size();
00312   }
00313 
00314   for(unsigned int ipart=0;ipart<partnames.size();++ipart) {
00315     //    if(phases[ipart]>=0) {
00316       //      apvphases->get()[partnames[ipart]] = (phases[ipart]+phasechange)%70;
00317     apvphases->get()[partnames[ipart]] = phases[ipart];
00318 
00319       //    }
00320   }
00321 
00322 
00323   iEvent.put(apvphases);
00324 
00325 }
00326 
00327 // ------------ method called once each job just before starting event loop  ------------
00328 void 
00329 APVCyclePhaseProducerFromL1TS::beginJob()
00330 {
00331 }
00332 
00333 // ------------ method called once each job just after ending the event loop  ------------
00334 void 
00335 APVCyclePhaseProducerFromL1TS::endJob() {
00336 }
00337 
00338 bool 
00339 APVCyclePhaseProducerFromL1TS::isBadRun(const unsigned int run) const {
00340 
00341   for(std::vector<std::pair<unsigned int, unsigned int> >::const_iterator runpair = m_badruns.begin();runpair!=m_badruns.end();++runpair) {
00342     if( run >= runpair->first && run <= runpair->second) return true;
00343   }
00344 
00345   return false;
00346 
00347 }
00348 
00349 //define this as a plug-in
00350 DEFINE_FWK_MODULE(APVCyclePhaseProducerFromL1TS);