CMS 3D CMS Logo

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