CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripTools
00004 // Class:      APVShotsAnalyzer
00005 // 
00013 //
00014 // Original Author:  Andrea Venturi
00015 //         Created:  Tue Jul 19 11:56:00 CEST 2009
00016 //
00017 //
00018 
00019 
00020 // system include files
00021 #include <memory>
00022 
00023 // user include files
00024 #include "TH1F.h"
00025 #include "TProfile.h"
00026 #include <vector>
00027 #include <string>
00028 
00029 #include "FWCore/Framework/interface/Frameworkfwd.h"
00030 #include "FWCore/Framework/interface/EDAnalyzer.h"
00031 
00032 #include "FWCore/Framework/interface/Event.h"
00033 #include "FWCore/Framework/interface/MakerMacros.h"
00034 
00035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00036 
00037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00038 
00039 #include "FWCore/ServiceRegistry/interface/Service.h"
00040 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00041 
00042 #include "FWCore/Utilities/interface/InputTag.h"
00043 
00044 #include "DataFormats/Common/interface/DetSetVector.h"
00045 #include "DataFormats/SiStripDigi/interface/SiStripDigi.h"
00046 
00047 #include "DQM/SiStripCommon/interface/APVShotFinder.h"
00048 #include "DQM/SiStripCommon/interface/APVShot.h"
00049 
00050 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
00051 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00052 
00053 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
00054 //******** Single include for the TkMap *************
00055 #include "DQM/SiStripCommon/interface/TkHistoMap.h" 
00056 //***************************************************
00057 
00058 //******** includes for the cabling *************
00059 #include "CalibFormats/SiStripObjects/interface/SiStripDetCabling.h"
00060 #include "CalibTracker/Records/interface/SiStripDetCablingRcd.h"
00061 #include "CondFormats/SiStripObjects/interface/FedChannelConnection.h"
00062 #include "DataFormats/SiStripCommon/interface/SiStripConstants.h"
00063 //***************************************************
00064 
00065 
00066 //
00067 // class decleration
00068 //
00069 
00070 class APVShotsAnalyzer : public edm::EDAnalyzer {
00071 public:
00072   explicit APVShotsAnalyzer(const edm::ParameterSet&);
00073   ~APVShotsAnalyzer();
00074 
00075 
00076 private:
00077   virtual void beginJob() ;
00078   virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
00079   virtual void endRun(const edm::Run&, const edm::EventSetup&) override;
00080   virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
00081   virtual void endJob() ;
00082 
00083   void updateDetCabling( const edm::EventSetup& setup );
00084 
00085       // ----------member data ---------------------------
00086 
00087   const edm::InputTag _digicollection;
00088   const edm::InputTag _historyProduct;
00089   const edm::InputTag _apvphasecoll;
00090   const std::string _phasepart;
00091   bool _zs;
00092   std::string _suffix;
00093   int _nevents;
00094 
00095   TH1F* _nShots;
00096   TH1F* _whichAPV;
00097   TH1F* _stripMult;
00098   TH1F* _median;
00099   TH1F* _subDetector;
00100   TH1F* _fed;
00101   TH2F* _channelvsfed;
00102 
00103   TProfile* _nShotsbxcycle;
00104   TProfile* _nShotsdbx;
00105   TProfile* _nShotsdbxincycle;
00106   TProfile* _nShotsbxcycleprev;
00107   TProfile* _nShotsdbxprev;
00108   TProfile* _nShotsdbxincycleprev;
00109 
00110   TH2F* _medianVsFED;
00111   TH2F* _nShotsVsFED;
00112 
00113   RunHistogramManager _rhm;
00114 
00115   TH1F** _nShotsrun;
00116   TProfile** _nShotsVsTimerun;
00117   TH1F** _whichAPVrun;
00118   TH1F** _stripMultrun;
00119   TH1F** _medianrun;
00120   TH1F** _subDetectorrun;
00121   TH1F** _fedrun;
00122 
00123   TkHistoMap *tkhisto,*tkhisto2; 
00124 
00125   // DetCabling
00126   bool _useCabling;
00127   uint32_t _cacheIdDet;  
00128   const SiStripDetCabling* _detCabling;  
00129   
00130 };
00131 
00132 //
00133 // constants, enums and typedefs
00134 //
00135 
00136 //
00137 // static data member definitions
00138 //
00139 
00140 //
00141 // constructors and destructor
00142 //
00143 APVShotsAnalyzer::APVShotsAnalyzer(const edm::ParameterSet& iConfig):
00144   _digicollection(iConfig.getParameter<edm::InputTag>("digiCollection")),
00145   _historyProduct(iConfig.getParameter<edm::InputTag>("historyProduct")),
00146   _apvphasecoll(iConfig.getParameter<edm::InputTag>("apvPhaseCollection")),
00147   _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition","None")),
00148   _zs(iConfig.getUntrackedParameter<bool>("zeroSuppressed",true)),
00149   _suffix(iConfig.getParameter<std::string>("mapSuffix")),
00150   _nevents(0),
00151   _rhm(),
00152   _useCabling(iConfig.getUntrackedParameter<bool>("useCabling",true)),
00153   _cacheIdDet(0),
00154   _detCabling(0)
00155 {
00156    //now do what ever initialization is needed
00157 
00158   if(!_zs) _suffix += "_notZS";
00159 
00160  edm::Service<TFileService> tfserv;
00161 
00162  _nShots = tfserv->make<TH1F>("nShots","Number of Shots per event",200,-0.5,199.5);
00163  _nShots->GetXaxis()->SetTitle("Shots");  _nShots->GetYaxis()->SetTitle("Events"); 
00164  _nShots->StatOverflows(kTRUE);
00165 
00166  _whichAPV = tfserv->make<TH1F>("whichAPV","APV with shots",6,-0.5,5.5);
00167  _whichAPV->GetXaxis()->SetTitle("APV");  _whichAPV->GetYaxis()->SetTitle("Shots"); 
00168 
00169  _stripMult = tfserv->make<TH1F>("stripMultiplicity","Shot Strip Multiplicity",129,-0.5,128.5);
00170  _stripMult->GetXaxis()->SetTitle("Number of Strips");  _stripMult->GetYaxis()->SetTitle("Shots");
00171 
00172  _median = tfserv->make<TH1F>("median","APV Shot charge median",256,-0.5,255.5);
00173  _median->GetXaxis()->SetTitle("Charge [ADC]");  _median->GetYaxis()->SetTitle("Shots");
00174 
00175  _subDetector = tfserv->make<TH1F>("subDets","SubDetector Shot distribution",10,-0.5,9.5);
00176  _subDetector->GetYaxis()->SetTitle("Shots");
00177 
00178  _nShotsbxcycle = tfserv->make<TProfile>("nShotsBXcycle","Number of shots vs APV cycle bin",70,-0.5,69.5);
00179  _nShotsbxcycle->GetXaxis()->SetTitle("Event BX mod(70)");  _nShotsbxcycle->GetYaxis()->SetTitle("APV shots"); 
00180 
00181  _nShotsdbx = tfserv->make<TProfile>("nShotsDBX","Number of shots vs #Delta(BX)",1000,-0.5,999.5);
00182  _nShotsdbx->GetXaxis()->SetTitle("Event #Delta(BX)");  _nShotsdbx->GetYaxis()->SetTitle("APV shots"); 
00183 
00184  _nShotsdbxincycle = tfserv->make<TProfile>("nShotsDBXincycle","Number of shots vs #Delta(BX) w.r.t. APV cycle",1000,-0.5,999.5);
00185  _nShotsdbxincycle->GetXaxis()->SetTitle("Event #Delta(BX) w.r.t. APV cycle");  _nShotsdbxincycle->GetYaxis()->SetTitle("APV shots"); 
00186 
00187  _nShotsbxcycleprev = tfserv->make<TProfile>("nShotsBXcycleprev","Number of shots vs APV cycle bin of previous L1A",70,-0.5,69.5);
00188  _nShotsbxcycleprev->GetXaxis()->SetTitle("Previous L1A BX mod(70)");  _nShotsbxcycleprev->GetYaxis()->SetTitle("APV shots"); 
00189 
00190  _nShotsdbxprev = tfserv->make<TProfile>("nShotsDBXprev","Number of shots vs #Delta(BX) of previous L1A",1000,-0.5,999.5);
00191  _nShotsdbxprev->GetXaxis()->SetTitle("Previous L1A #Delta(BX)");  _nShotsdbxprev->GetYaxis()->SetTitle("APV shots"); 
00192 
00193  _nShotsdbxincycleprev = tfserv->make<TProfile>("nShotsDBXincycleprev","Number of shots vs #Delta(BX) w.r.t. APV cycle of previous L1A",1000,-0.5,999.5);
00194  _nShotsdbxincycleprev->GetXaxis()->SetTitle("Previous L1A #Delta(BX) w.r.t. APV cycle");  _nShotsdbxincycleprev->GetYaxis()->SetTitle("APV shots"); 
00195 
00196  _nShotsrun = _rhm.makeTH1F("nShotsrun","Number of Shots per event",200,-0.5,199.5);
00197  _nShotsVsTimerun  = _rhm.makeTProfile("nShotsVsTimerun","Mean number of shots vs orbit number",4*500,0,500*262144);
00198  _whichAPVrun = _rhm.makeTH1F("whichAPVrun","APV with shots",6,-0.5,5.5);
00199  _stripMultrun = _rhm.makeTH1F("stripMultiplicityrun","Shot Strip Multiplicity",129,-0.5,128.5);
00200  _medianrun = _rhm.makeTH1F("medianrun","APV Shot charge median",256,-0.5,255.5);
00201  _subDetectorrun = _rhm.makeTH1F("subDetsrun","SubDetector Shot distribution",10,-0.5,9.5);
00202 
00203  if (_useCabling) {
00204    _fed = tfserv->make<TH1F>("fed","FED Shot distribution",440,50,490);
00205    _fed->GetYaxis()->SetTitle("Shots");
00206    _fedrun = _rhm.makeTH1F("fedrun","FED Shot distribution",440,50,490);
00207 
00208    _channelvsfed = tfserv->make<TH2F>("channelvsfed","Channel vs FED Shot distribution",440,50,490,97,-0.5,96.5);
00209    _channelvsfed->GetXaxis()->SetTitle("FED");    _channelvsfed->GetYaxis()->SetTitle("Channel");
00210 
00211 
00212    _nShotsVsFED = tfserv->make<TH2F>("nShotsVsFED","Number of Shots per event vs fedid",440,50,490,200,-0.5,199.5);
00213    _nShotsVsFED->GetXaxis()->SetTitle("fedId");  _nShots->GetYaxis()->SetTitle("Shots");  _nShots->GetZaxis()->SetTitle("Events");
00214    _nShotsVsFED->StatOverflows(kTRUE);
00215 
00216    _medianVsFED = tfserv->make<TH2F>("medianVsFED","APV Shot charge median vs fedid",440,50,490,256,-0.5,255.5);
00217    _medianVsFED->GetXaxis()->SetTitle("fedId");_medianVsFED->GetYaxis()->SetTitle("Charge [ADC]");  _median->GetZaxis()->SetTitle("Shots");
00218  }
00219 
00220  tkhisto      =new TkHistoMap("ShotMultiplicity","ShotMultiplicity",-1); 
00221  tkhisto2      =new TkHistoMap("StripMultiplicity","StripMultiplicity",-1); 
00222 }
00223 
00224 
00225 APVShotsAnalyzer::~APVShotsAnalyzer()
00226 {
00227  
00228    // do anything here that needs to be done at desctruction time
00229    // (e.g. close files, deallocate resources etc.)
00230   if ( _detCabling ) _detCabling = 0;
00231 
00232 }
00233 
00234 
00235 //
00236 // member functions
00237 //
00238 
00239 // ------------ method called to for each event  ------------
00240 void
00241 APVShotsAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00242 {
00243    using namespace edm;
00244 
00245    if (_useCabling){
00246      //retrieve cabling
00247      updateDetCabling( iSetup );
00248    }
00249 
00250    _nevents++;
00251 
00252    edm::Handle<EventWithHistory> he;
00253    iEvent.getByLabel(_historyProduct,he);
00254 
00255    edm::Handle<APVCyclePhaseCollection> apvphase;
00256    iEvent.getByLabel(_apvphasecoll,apvphase);
00257 
00258    int thephase = APVCyclePhaseCollection::invalid;
00259    if(apvphase.isValid() && !apvphase.failedToGet()) {
00260      thephase = apvphase->getPhase(_phasepart);
00261    }
00262    bool isphaseok = (thephase!=APVCyclePhaseCollection::invalid &&
00263                      thephase!=APVCyclePhaseCollection::multiphase &&
00264                      thephase!=APVCyclePhaseCollection::nopartition);
00265 
00266    Handle<edm::DetSetVector<SiStripDigi> > digis;
00267    iEvent.getByLabel(_digicollection,digis);
00268 
00269    // loop on detector with digis
00270 
00271    int nshots=0;
00272    std::vector<int> nshotsperFed;
00273 
00274    const uint16_t lNumFeds = sistrip::FED_ID_MAX-sistrip::FED_ID_MIN+1;
00275    if (_useCabling){
00276      nshotsperFed.resize(lNumFeds,0);
00277    }
00278 
00279    APVShotFinder apvsf(*digis,_zs);
00280    const std::vector<APVShot>& shots = apvsf.getShots();
00281 
00282    for(std::vector<APVShot>::const_iterator shot=shots.begin();shot!=shots.end();++shot) {
00283      if(shot->isGenuine()) {
00284 
00285        //get the fedid from the detid
00286 
00287        uint32_t det=shot->detId();
00288        if (_useCabling){
00289 
00290          int apvPair = shot->apvNumber()/2;
00291          LogDebug("APVPair") << apvPair;
00292 
00293          const FedChannelConnection& theConn = _detCabling->getConnection( det , apvPair);
00294 
00295          int lChannelId = -1;
00296          int thelFEDId = -1;
00297          if(theConn.isConnected()) {
00298            lChannelId = theConn.fedCh();
00299            thelFEDId = theConn.fedId();
00300          }
00301          else {
00302            edm::LogWarning("ConnectionNotFound") << "connection of det " << det << " APV pair " << apvPair << " not found";
00303          }
00304          LogDebug("FED channels") << thelFEDId << " " << lChannelId ;
00305 
00306          const std::vector<const FedChannelConnection *> & conns = _detCabling->getConnections( det );
00307 
00308          if (!(conns.size())) continue;
00309          uint16_t lFedId = 0;
00310          for (uint32_t ch = 0; ch<conns.size(); ch++) {
00311            if(conns[ch] && conns[ch]->isConnected()) {
00312              LogDebug("Dump") << *(conns[ch]);
00313              LogDebug("ReadyForFEDid") << "Ready for FED id " << ch;
00314              lFedId = conns[ch]->fedId();
00315              LogDebug("FEDid") << "obtained FED id " << ch << " " << lFedId;
00316              //uint16_t lFedCh = conns[ch]->fedCh();
00317              
00318              if (lFedId < sistrip::FED_ID_MIN || lFedId > sistrip::FED_ID_MAX){
00319                edm::LogWarning("InvalidFEDid") << lFedId << " for detid " << det << " connection " << ch;
00320                continue;
00321              }
00322              else break;
00323            }
00324          }
00325          if (lFedId < sistrip::FED_ID_MIN || lFedId > sistrip::FED_ID_MAX){
00326            edm::LogWarning("NoValidFEDid") << lFedId <<  "found for detid " << det;
00327            continue;
00328          }
00329 
00330          if(lFedId != thelFEDId) {
00331            edm::LogWarning("FEDidMismatch") << " Mismatch in FED id for det " << det << " APV pair " 
00332                                             << apvPair << " : " << lFedId << " vs " << thelFEDId;
00333          }
00334 
00335          LogDebug("FillingArray") << nshotsperFed.size() << " " << lFedId-sistrip::FED_ID_MIN;  
00336          ++nshotsperFed[lFedId-sistrip::FED_ID_MIN];
00337          
00338          LogDebug("ReadyToBeFilled") << " ready to be filled with " << thelFEDId << " " << lChannelId;
00339          _channelvsfed->Fill(thelFEDId,lChannelId);
00340          LogDebug("Filled") << " filled with " << thelFEDId << " " << lChannelId;
00341 
00342          _fed->Fill(lFedId);
00343 
00344          if(_fedrun && *_fedrun) (*_fedrun)->Fill(lFedId); 
00345          _medianVsFED->Fill(lFedId,shot->median());
00346 
00347 
00348        }
00349 
00350        ++nshots;
00351 
00352 
00353        _whichAPV->Fill(shot->apvNumber());
00354        _median->Fill(shot->median());
00355        _stripMult->Fill(shot->nStrips());
00356        _subDetector->Fill(shot->subDet());
00357 
00358        if(_whichAPVrun && *_whichAPVrun) (*_whichAPVrun)->Fill(shot->apvNumber());
00359        if(_medianrun && *_medianrun) (*_medianrun)->Fill(shot->median());
00360        if(_stripMultrun && *_stripMultrun) (*_stripMultrun)->Fill(shot->nStrips());
00361        if(_subDetectorrun && *_subDetectorrun) (*_subDetectorrun)->Fill(shot->subDet());
00362 
00363        tkhisto2->fill(det,shot->nStrips());;
00364        tkhisto->add(det,1);
00365 
00366        
00367        
00368      }
00369    }
00370 
00371      _nShots->Fill(nshots);
00372      if(_nShotsrun && *_nShotsrun) (*_nShotsrun)->Fill(nshots);
00373 
00374    _nShotsdbx->Fill(he->deltaBX(),nshots);
00375    _nShotsdbxprev->Fill(he->deltaBX(),nshots);
00376    if(isphaseok) {
00377      _nShotsbxcycle->Fill(he->absoluteBXinCycle(thephase)%70,nshots);
00378      _nShotsdbxincycle->Fill(he->deltaBXinCycle(thephase),nshots);
00379      _nShotsbxcycleprev->Fill(he->absoluteBXinCycle(1,thephase)%70,nshots);
00380      _nShotsdbxincycleprev->Fill(he->deltaBXinCycle(1,2,thephase),nshots);
00381    }
00382 
00383    if (_useCabling){
00384      for (uint16_t lFed(0); lFed<lNumFeds; lFed++){
00385        _nShotsVsFED->Fill(lFed+sistrip::FED_ID_MIN,nshotsperFed[lFed]);
00386      }
00387    }
00388 
00389    if(_nShotsVsTimerun && *_nShotsVsTimerun) (*_nShotsVsTimerun)->Fill(iEvent.orbitNumber(),nshots);
00390    
00391 
00392 }
00393 
00394 void 
00395 APVShotsAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup&)
00396 {
00397 
00398 
00399   _rhm.beginRun(iRun);
00400   
00401   if(_nShotsrun && *_nShotsrun) {
00402     (*_nShotsrun)->GetXaxis()->SetTitle("Shots");  (*_nShotsrun)->GetYaxis()->SetTitle("Events"); 
00403     (*_nShotsrun)->StatOverflows(kTRUE);
00404   }
00405   
00406   if(_nShotsVsTimerun && *_nShotsVsTimerun) {
00407     (*_nShotsVsTimerun)->GetXaxis()->SetTitle("Orbit");  (*_nShotsVsTimerun)->GetYaxis()->SetTitle("Number of Shots");
00408     (*_nShotsVsTimerun)->SetBit(TH1::kCanRebin);
00409   }
00410   
00411   if(_whichAPVrun && *_whichAPVrun) {
00412     (*_whichAPVrun)->GetXaxis()->SetTitle("APV");  (*_whichAPVrun)->GetYaxis()->SetTitle("Shots"); 
00413   }    
00414   
00415   if(_stripMultrun && *_stripMultrun) {
00416     (*_stripMultrun)->GetXaxis()->SetTitle("Number of Strips");  (*_stripMultrun)->GetYaxis()->SetTitle("Shots");
00417   }
00418   
00419   if(_medianrun && *_medianrun) {
00420     (*_medianrun)->GetXaxis()->SetTitle("Charge [ADC]");  (*_medianrun)->GetYaxis()->SetTitle("Shots");
00421   }
00422   
00423   if(_subDetectorrun && *_subDetectorrun) {
00424     (*_subDetectorrun)->GetYaxis()->SetTitle("Shots");
00425   }
00426   
00427   if (_useCabling) {
00428     if(_fedrun && *_fedrun) {
00429       (*_fedrun)->GetYaxis()->SetTitle("Shots");
00430     }
00431   }
00432   
00433 }
00434 
00435 void 
00436 APVShotsAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup&)
00437 {
00438 }
00439 
00440 
00441 // ------------ method called once each job just before starting event loop  ------------
00442 void 
00443 APVShotsAnalyzer::beginJob()
00444 {
00445 
00446 }
00447 
00448 // ------------ method called once each job just after ending the event loop  ------------
00449 void 
00450 APVShotsAnalyzer::endJob() {
00451 
00452   edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
00453 
00454 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
00455   TrackerMap tkmap,tkmap2;
00456 
00457   tkmap.setPalette(1);
00458   tkmap2.setPalette(1);
00459   tkhisto->dumpInTkMap(&tkmap);
00460   tkhisto2->dumpInTkMap(&tkmap2);
00461   std::string tkshotmultmapname = "ShotMultiplicity_" + _suffix + ".png";
00462   tkmap.save(true,0,0,tkshotmultmapname);
00463   std::string tkstripmultmapname = "StripMultiplicity_" + _suffix + ".png";
00464   tkmap2.save(true,0,0,tkstripmultmapname);
00465 
00466   std::string rootmapname = "TKMap_"+_suffix+".root";
00467   tkhisto->save(rootmapname);
00468   tkhisto2->save(rootmapname);
00469 }
00470 
00471 
00472 void APVShotsAnalyzer::updateDetCabling( const edm::EventSetup& setup )
00473 {
00474   if (_useCabling){
00475     uint32_t cache_id = setup.get<SiStripDetCablingRcd>().cacheIdentifier();//.get( cabling_ );
00476    
00477     if ( _cacheIdDet != cache_id ) { // If the cache ID has changed since the last update...
00478       // Update the cabling object
00479       edm::ESHandle<SiStripDetCabling> c;
00480       setup.get<SiStripDetCablingRcd>().get( c );
00481       _detCabling = c.product();
00482       _cacheIdDet = cache_id;
00483     } // end of new cache ID check
00484   }
00485 }
00486 
00487 
00488 
00489 
00490 //define this as a plug-in
00491 DEFINE_FWK_MODULE(APVShotsAnalyzer);