CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/DPGAnalysis/SiStripTools/plugins/APVCyclePhaseProducerFromL1ABC.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripTools
00004 // Class:      APVCyclePhaseProducerFromL1ABC
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/ServiceRegistry/interface/Service.h"
00034 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00035 
00036 #include "FWCore/Utilities/interface/InputTag.h"
00037 
00038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00039 #include "FWCore/Utilities/interface/Exception.h"
00040 
00041 #include <map>
00042 #include <vector>
00043 #include <string>
00044 
00045 #include "TH1F.h"
00046 #include "TProfile.h"
00047 
00048 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
00049 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00050 
00051 //
00052 // class decleration
00053 //
00054 
00055 class APVCyclePhaseProducerFromL1ABC : public edm::EDProducer {
00056    public:
00057       explicit APVCyclePhaseProducerFromL1ABC(const edm::ParameterSet&);
00058       ~APVCyclePhaseProducerFromL1ABC();
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       // ----------member data ---------------------------
00068 
00069   const edm::InputTag _l1abccollection;
00070   const std::vector<std::string> _defpartnames;
00071   const std::vector<int> _defphases;
00072   const int _orbitoffsetSOR; 
00073   const bool _wantHistos;
00074   TH1F* _hbx;
00075   TH1F* _hdbx;
00076   TH1F* _hdorbit;
00077   const unsigned int _firstgoodrun;
00078   std::map<unsigned int, long long> _offsets;
00079   long long _curroffset;
00080   unsigned int _curroffevent;
00081 
00082 };
00083 
00084 //
00085 // constants, enums and typedefs
00086 //
00087 
00088 
00089 //
00090 // static data member definitions
00091 //
00092 
00093 //
00094 // constructors and destructor
00095 //
00096 APVCyclePhaseProducerFromL1ABC::APVCyclePhaseProducerFromL1ABC(const edm::ParameterSet& iConfig):
00097   _l1abccollection(iConfig.getParameter<edm::InputTag>("l1ABCCollection")),
00098   _defpartnames(iConfig.getParameter<std::vector<std::string> >("defaultPartitionNames")),
00099   _defphases(iConfig.getParameter<std::vector<int> >("defaultPhases")),
00100   _orbitoffsetSOR(iConfig.getParameter<int>("StartOfRunOrbitOffset")),
00101   _wantHistos(iConfig.getUntrackedParameter<bool>("wantHistos",false)),
00102   _hbx(0),_hdbx(0),_hdorbit(0),_firstgoodrun(110878),
00103   _offsets(), _curroffset(0), _curroffevent(0)
00104 {
00105 
00106   produces<APVCyclePhaseCollection,edm::InEvent>();
00107 
00108    //now do what ever other initialization is needed
00109 
00110 }
00111 
00112 
00113 APVCyclePhaseProducerFromL1ABC::~APVCyclePhaseProducerFromL1ABC()
00114 {
00115  
00116    // do anything here that needs to be done at desctruction time
00117    // (e.g. close files, deallocate resources etc.)
00118 
00119 }
00120 
00121 
00122 //
00123 // member functions
00124 //
00125 
00126 // ------------ method called to produce the data  ------------
00127 void
00128 APVCyclePhaseProducerFromL1ABC::beginRun(edm::Run& iRun, const edm::EventSetup& iSetup) 
00129 
00130 {
00131 
00132   // reset offset vector
00133 
00134   _offsets.clear();
00135   edm::LogInfo("AbsoluteBXOffsetReset") << "Absolute BX offset map reset";
00136 
00137 
00138   if(_wantHistos) {
00139 
00140     edm::Service<TFileService> tfserv;
00141 
00142     char dirname[300];
00143     sprintf(dirname,"run_%d",iRun.run());
00144     TFileDirectory subrun = tfserv->mkdir(dirname);
00145 
00146     _hbx = subrun.make<TH1F>("l1abcbx","BX number from L1ABC",4096,-0.5,4095.5);
00147     _hbx->GetXaxis()->SetTitle("BX");     _hbx->GetYaxis()->SetTitle("Events"); 
00148 
00149     _hdbx = subrun.make<TH1F>("dbx","BX number difference",4096*2-1,-4095.5,4095.5);
00150     _hdbx->GetXaxis()->SetTitle("#DeltaBX");     _hdbx->GetYaxis()->SetTitle("Events"); 
00151     
00152     _hdorbit = subrun.make<TH1F>("dorbit","Orbit Number difference",9999,-4999.5,4999.5);
00153     _hdorbit->GetXaxis()->SetTitle("#Deltaorbit");     _hdorbit->GetYaxis()->SetTitle("Events"); 
00154     
00155   }
00156 
00157   if(iRun.run() < _firstgoodrun) {
00158     edm::LogInfo("UnreliableMissingL1AcceptBunchCrossingCollection") << 
00159       "In this run L1AcceptBunchCrossingCollection is missing or unreliable: default phases will be used"; 
00160   }
00161 
00162 }
00163 
00164 void 
00165 APVCyclePhaseProducerFromL1ABC::endRun(edm::Run&, const edm::EventSetup&)
00166 {
00167   // summary of absolute bx offset vector
00168 
00169   edm::LogInfo("L1AcceptBunchCrossingAbsoluteBXOffsetSummary") << "Absolute BX offset summary:";
00170   for(std::map<unsigned int, long long>::const_iterator offset=_offsets.begin();offset!=_offsets.end();++offset) {
00171     edm::LogVerbatim("L1AcceptBunchCrossingAbsoluteBXOffsetSummary") << offset->first << " " << offset->second;
00172   }
00173 
00174 }
00175 
00176 
00177 void
00178 APVCyclePhaseProducerFromL1ABC::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00179 
00180   using namespace edm;
00181   
00182   std::auto_ptr<APVCyclePhaseCollection> apvphases(new APVCyclePhaseCollection() );
00183   
00184 
00185   const std::vector<int>& phases = _defphases;
00186   const std::vector<std::string>& partnames = _defpartnames;
00187 
00188   // Look for the L1AcceptBunchCrossingCollection
00189 
00190   int phasechange = 0;
00191 
00192   if(iEvent.run() >= _firstgoodrun ) {
00193     
00194     Handle<L1AcceptBunchCrossingCollection > pIn;
00195     iEvent.getByLabel(_l1abccollection,pIn);
00196     
00197     // offset computation
00198     
00199     long long orbitoffset = _orbitoffsetSOR;
00200     int bxoffset = 0;
00201     
00202     for(L1AcceptBunchCrossingCollection::const_iterator l1abc=pIn->begin();l1abc!=pIn->end();++l1abc) {
00203       if(l1abc->l1AcceptOffset()==0) {
00204         if(l1abc->eventType()!=0) {
00205           orbitoffset = (long long)iEvent.orbitNumber() - (long long)l1abc->orbitNumber() ;
00206           bxoffset = iEvent.bunchCrossing() - l1abc->bunchCrossing();
00207           
00208           if(_wantHistos) {
00209             _hbx->Fill(l1abc->bunchCrossing());
00210             _hdbx->Fill(bxoffset);
00211             _hdorbit->Fill(orbitoffset);
00212           }
00213         }
00214         else {
00215           edm::LogWarning("L1AcceptBunchCrossingNoType") << "L1AcceptBunchCrossing with no type found: ";
00216           for(L1AcceptBunchCrossingCollection::const_iterator debu=pIn->begin();debu!=pIn->end();++debu) {
00217             edm::LogPrint("L1AcceptBunchCrossingNoType") << *debu;
00218           }
00219         }
00220       }
00221     }
00222 
00223     long long absbxoffset = orbitoffset*3564 + bxoffset;
00224 
00225     if(orbitoffset != _orbitoffsetSOR) phasechange = (orbitoffset*3564)%70;
00226 
00227      if(_offsets.size()==0) {
00228        _curroffset = absbxoffset;
00229        _curroffevent = iEvent.id().event();
00230        _offsets[iEvent.id().event()] = absbxoffset;
00231      }
00232      else {
00233        if(_curroffset != absbxoffset || iEvent.id().event() < _curroffevent ) {
00234 
00235          if( _curroffset != absbxoffset) {
00236            edm::LogInfo("L1AcceptBunchCrossingAbsoluteBXOffsetChanged") << "Absolute BX offset changed from " 
00237                                                                         << _curroffset << " to "
00238                                                                         << absbxoffset << " at orbit "
00239                                                                         << iEvent.orbitNumber() << " and BX "
00240                                                                         << iEvent.bunchCrossing();
00241            for(L1AcceptBunchCrossingCollection::const_iterator l1abc=pIn->begin();l1abc!=pIn->end();++l1abc) {
00242              edm::LogVerbatim("L1AcceptBunchCrossingAbsoluteBXOffsetChanged") << *l1abc;
00243            }
00244          }
00245 
00246          _curroffset = absbxoffset;
00247          _curroffevent = iEvent.id().event();
00248          _offsets[iEvent.id().event()] = absbxoffset;
00249        }
00250      }
00251 
00252     
00253   }
00254   
00255 
00256 
00257   if(phases.size() < partnames.size() ) {
00258     // throw exception
00259     throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent phases/partitions vector sizes: " 
00260                                              << phases.size() << " " 
00261                                              << partnames.size();
00262   }
00263 
00264   for(unsigned int ipart=0;ipart<partnames.size();++ipart) {
00265     if(phases[ipart]>=0) {
00266       apvphases->get()[partnames[ipart]] = (phases[ipart]+phasechange)%70;
00267 
00268     }
00269   }
00270 
00271 
00272   iEvent.put(apvphases);
00273 
00274 }
00275 
00276 // ------------ method called once each job just before starting event loop  ------------
00277 void 
00278 APVCyclePhaseProducerFromL1ABC::beginJob()
00279 {
00280 }
00281 
00282 // ------------ method called once each job just after ending the event loop  ------------
00283 void 
00284 APVCyclePhaseProducerFromL1ABC::endJob() {
00285 }
00286 
00287 //define this as a plug-in
00288 DEFINE_FWK_MODULE(APVCyclePhaseProducerFromL1ABC);