Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:06:35

0001 // -*- C++ -*-
0002 //
0003 // Package:    MultiplicityTimeCorrelations
0004 // Class:      MultiplicityTimeCorrelations
0005 //
0006 /**\class MultiplicityTimeCorrelations MultiplicityTimeCorrelations.cc DPGAnalysis/SiStripTools/src/MultiplicityTimeCorrelations.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Andrea Venturi
0015 //         Created:  Mon Oct 27 17:37:53 CET 2008
0016 // $Id: MultiplicityTimeCorrelations.cc,v 1.1 2011/03/10 16:15:13 venturia Exp $
0017 //
0018 //
0019 
0020 // system include files
0021 #include <memory>
0022 
0023 // user include files
0024 
0025 #include <vector>
0026 #include <map>
0027 #include <climits>
0028 
0029 #include "FWCore/Framework/interface/Frameworkfwd.h"
0030 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0031 
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/Run.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035 
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 
0038 #include "TH1F.h"
0039 
0040 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0041 
0042 #include "FWCore/ServiceRegistry/interface/Service.h"
0043 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0044 
0045 #include "FWCore/Utilities/interface/InputTag.h"
0046 
0047 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
0048 
0049 #include "DPGAnalysis/SiStripTools/interface/SiStripTKNumbers.h"
0050 #include "DPGAnalysis/SiStripTools/interface/DigiBXCorrHistogramMaker.h"
0051 
0052 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
0053 #include "DPGAnalysis/SiStripTools/interface/EventWithHistoryFilter.h"
0054 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
0055 
0056 //
0057 // class decleration
0058 //
0059 
0060 class MultiplicityTimeCorrelations : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0061 public:
0062   explicit MultiplicityTimeCorrelations(const edm::ParameterSet&);
0063   ~MultiplicityTimeCorrelations() override;
0064 
0065 private:
0066   void beginJob() override;
0067   void analyze(const edm::Event&, const edm::EventSetup&) override;
0068   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0069   void endRun(const edm::Run&, const edm::EventSetup&) override {}
0070   void endJob() override;
0071 
0072   // ----------member data ---------------------------
0073 
0074   DigiBXCorrHistogramMaker<EventWithHistory> _digibxcorrhmevent;
0075   EventWithHistoryFilter _evfilter;
0076 
0077   std::map<int, std::map<unsigned int, TH1F*> > _dbxhistos;
0078   /*
0079   std::map<int,TH1F*> _dbxtkhistos;
0080   std::map<int,TH1F*> _dbxtibhistos;
0081   std::map<int,TH1F*> _dbxtidhistos;
0082   std::map<int,TH1F*> _dbxtobhistos;
0083   std::map<int,TH1F*> _dbxtechistos;
0084   */
0085 
0086   edm::InputTag _hecollection;
0087   edm::EDGetTokenT<EventWithHistory> _hecollectionToken;
0088   edm::EDGetTokenT<APVCyclePhaseCollection> _apvphasecollToken;
0089   edm::EDGetTokenT<std::map<unsigned int, int> > _multiplicityMapToken;
0090   std::map<unsigned int, std::string> _subdets;
0091   std::map<unsigned int, int> _binmax;
0092 
0093   int _loworbit;
0094   int _highorbit;
0095 
0096   int _mindbx;
0097   int _mintrpltdbx;
0098 
0099   SiStripTKNumbers _trnumb;
0100 
0101   std::vector<int> _dbxbins;
0102 };
0103 
0104 //
0105 // constants, enums and typedefs
0106 //
0107 
0108 //
0109 // static data member definitions
0110 //
0111 
0112 //
0113 // constructors and destructor
0114 //
0115 MultiplicityTimeCorrelations::MultiplicityTimeCorrelations(const edm::ParameterSet& iConfig)
0116     : _digibxcorrhmevent(iConfig, consumesCollector()),
0117       _evfilter(),
0118       _hecollection(iConfig.getParameter<edm::InputTag>("historyProduct")),
0119       _hecollectionToken(consumes<EventWithHistory>(_hecollection)),
0120       _apvphasecollToken(consumes<APVCyclePhaseCollection>(iConfig.getParameter<edm::InputTag>("apvPhaseCollection"))),
0121       _multiplicityMapToken(
0122           mayConsume<std::map<unsigned int, int> >(iConfig.getParameter<edm::InputTag>("multiplicityMap"))),
0123       _subdets(),
0124       _binmax(),
0125       _loworbit(iConfig.getUntrackedParameter<int>("lowedgeOrbit")),
0126       _highorbit(iConfig.getUntrackedParameter<int>("highedgeOrbit")),
0127       _mindbx(iConfig.getUntrackedParameter<int>("minDBX")),
0128       _mintrpltdbx(iConfig.getUntrackedParameter<int>("minTripletDBX")),
0129       _trnumb(),
0130       _dbxbins(iConfig.getUntrackedParameter<std::vector<int> >("dbxBins")) {
0131   //now do what ever initialization is needed
0132 
0133   usesResource(TFileService::kSharedResource);
0134 
0135   // configure the filter
0136 
0137   edm::ParameterSet filterConfig;
0138   filterConfig.addUntrackedParameter<edm::InputTag>("historyProduct", _hecollection);
0139   if (_mindbx > 0) {
0140     std::vector<int> dbxrange;
0141     dbxrange.push_back(_mindbx + 1);
0142     dbxrange.push_back(-1);
0143     filterConfig.addUntrackedParameter<std::vector<int> >("dbxRange", dbxrange);
0144   }
0145   if (_mintrpltdbx > 0) {
0146     std::vector<int> dbxtrpltrange;
0147     dbxtrpltrange.push_back(_mintrpltdbx + 1);
0148     dbxtrpltrange.push_back(-1);
0149     filterConfig.addUntrackedParameter<std::vector<int> >("dbxTripletRange", dbxtrpltrange);
0150   }
0151 
0152   _evfilter.set(filterConfig, consumesCollector());
0153 
0154   //
0155 
0156   edm::Service<TFileService> tfserv;
0157 
0158   // create map of labels
0159 
0160   std::vector<edm::ParameterSet> wantedsubds(
0161       iConfig.getUntrackedParameter<std::vector<edm::ParameterSet> >("wantedSubDets"));
0162 
0163   for (std::vector<edm::ParameterSet>::iterator ps = wantedsubds.begin(); ps != wantedsubds.end(); ++ps) {
0164     _subdets[ps->getParameter<unsigned int>("detSelection")] = ps->getParameter<std::string>("detLabel");
0165     _binmax[ps->getParameter<unsigned int>("detSelection")] = ps->getParameter<int>("binMax");
0166   }
0167   std::map<int, std::string> labels;
0168 
0169   for (std::map<unsigned int, std::string>::const_iterator subd = _subdets.begin(); subd != _subdets.end(); ++subd) {
0170     labels[int(subd->first)] = subd->second;
0171   }
0172 
0173   //
0174 
0175   _digibxcorrhmevent.book("EventProcs", labels);
0176 
0177   TFileDirectory subdbxbin = tfserv->mkdir("DBXDebugging");
0178 
0179   for (std::vector<int>::const_iterator bin = _dbxbins.begin(); bin != _dbxbins.end(); bin++) {
0180     char hname[200];
0181     char htitle[200];
0182 
0183     edm::LogInfo("DBXHistosBinMaxValue") << "Setting bin max values";
0184 
0185     for (std::map<unsigned int, std::string>::const_iterator subd = _subdets.begin(); subd != _subdets.end(); ++subd) {
0186       if (_binmax.find(subd->first) == _binmax.end()) {
0187         edm::LogVerbatim("DBXHistosNotConfiguredBinMax")
0188             << "Bin max for " << subd->second << " not configured: " << _trnumb.nstrips(int(subd->first)) << " used";
0189         _binmax[subd->first] = _trnumb.nstrips(int(subd->first));
0190       }
0191 
0192       edm::LogVerbatim("DBXHistosBinMaxValue") << "Bin max for " << subd->second << " is " << _binmax[subd->first];
0193 
0194       sprintf(hname, "sumn%sdigi_%d", subd->second.c_str(), *bin);
0195       sprintf(htitle, "%s digi multiplicity at DBX = %d", subd->second.c_str(), *bin);
0196       LogDebug("DBXDebug") << "creating histogram " << hname << " " << htitle;
0197       _dbxhistos[*bin][subd->first] =
0198           subdbxbin.make<TH1F>(hname, htitle, 1000, 0., _binmax[subd->first] / (20 * 1000) * 1000);
0199       _dbxhistos[*bin][subd->first]->GetXaxis()->SetTitle("Number of Digis");
0200     }
0201     /*
0202     sprintf(hname,"sumntkdigi_%d",*bin);
0203     sprintf(htitle,"TK digi multiplicity at DBX = %d",*bin);
0204     LogDebug("DBXDebug") << "creating histogram " << hname << " " << htitle;
0205     _dbxtkhistos[*bin]= subdbxbin.make<TH1F>(hname,htitle,1000,0.,_trnumb.nstrips(0)/(20*1000)*1000);
0206     _dbxtkhistos[*bin]->GetXaxis()->SetTitle("Number of Digis");
0207 
0208     sprintf(hname,"sumntibdigi_%d",*bin);
0209     sprintf(htitle,"TIB digi multiplicity at DBX = %d",*bin);
0210     LogDebug("DBXDebug") << "creating histogram " << hname << " " << htitle;
0211     _dbxtibhistos[*bin]= subdbxbin.make<TH1F>(hname,htitle,1000,0.,_trnumb.nstrips(SiStripDetId::TIB)/(20*1000)*1000);
0212     _dbxtibhistos[*bin]->GetXaxis()->SetTitle("Number of Digis");
0213 
0214     sprintf(hname,"sumntiddigi_%d",*bin);
0215     sprintf(htitle,"TID digi multiplicity at DBX = %d",*bin);
0216     LogDebug("DBXDebug") << "creating histogram " << hname << " " << htitle;
0217     _dbxtidhistos[*bin]= subdbxbin.make<TH1F>(hname,htitle,1000,0.,_trnumb.nstrips(SiStripDetId::TID)/(20*1000)*1000);
0218     _dbxtidhistos[*bin]->GetXaxis()->SetTitle("Number of Digis");
0219 
0220     sprintf(hname,"sumntobdigi_%d",*bin);
0221     sprintf(htitle,"TOB digi multiplicity at DBX = %d",*bin);
0222     LogDebug("DBXDebug") << "creating histogram " << hname << " " << htitle;
0223     _dbxtobhistos[*bin]= subdbxbin.make<TH1F>(hname,htitle,1000,0.,_trnumb.nstrips(SiStripDetId::TOB)/(20*1000)*1000);
0224     _dbxtobhistos[*bin]->GetXaxis()->SetTitle("Number of Digis");
0225 
0226     sprintf(hname,"sumntecdigi_%d",*bin);
0227     sprintf(htitle,"TEC digi multiplicity at DBX = %d",*bin);
0228     LogDebug("DBXDebug") << "creating histogram " << hname << " " << htitle;
0229     _dbxtechistos[*bin]= subdbxbin.make<TH1F>(hname,htitle,1000,0.,_trnumb.nstrips(SiStripDetId::TEC)/(20*1000)*1000);
0230     _dbxtechistos[*bin]->GetXaxis()->SetTitle("Number of Digis");
0231     */
0232   }
0233 }
0234 
0235 MultiplicityTimeCorrelations::~MultiplicityTimeCorrelations() {
0236   // do anything here that needs to be done at desctruction time
0237   // (e.g. close files, deallocate resources etc.)
0238 }
0239 
0240 //
0241 // member functions
0242 //
0243 
0244 // ------------ method called to for each event  ------------
0245 void MultiplicityTimeCorrelations::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0246   using namespace edm;
0247 
0248   // get Phase
0249 
0250   Handle<APVCyclePhaseCollection> apvphase;
0251   iEvent.getByToken(_apvphasecollToken, apvphase);
0252 
0253   // get HE
0254 
0255   Handle<EventWithHistory> he;
0256   iEvent.getByToken(_hecollectionToken, he);
0257 
0258   // check if the event is selected
0259 
0260   if ((_loworbit < 0 || iEvent.orbitNumber() >= _loworbit) && (_highorbit < 0 || iEvent.orbitNumber() <= _highorbit)) {
0261     if (_evfilter.selected(iEvent, iSetup)) {
0262       //Compute digi multiplicity
0263       /*
0264       int ntkdigi=0;
0265       int ntibdigi=0;
0266       int ntiddigi=0;
0267       int ntobdigi=0;
0268       int ntecdigi=0;
0269       */
0270       Handle<std::map<unsigned int, int> > mults;
0271       iEvent.getByToken(_multiplicityMapToken, mults);
0272 
0273       // create map of digi multiplicity
0274 
0275       std::map<int, int> digimap;
0276       for (std::map<unsigned int, int>::const_iterator mult = mults->begin(); mult != mults->end(); ++mult) {
0277         if (_subdets.find(mult->first) != _subdets.end())
0278           digimap[int(mult->first)] = mult->second;
0279       }
0280 
0281       _digibxcorrhmevent.fill(*he, digimap, apvphase);
0282 
0283       // fill debug histos
0284 
0285       if (he->depth() != 0) {
0286         long long dbx = he->deltaBX();
0287 
0288         if (_dbxhistos.find(dbx) != _dbxhistos.end()) {
0289           for (std::map<unsigned int, int>::const_iterator ndigi = mults->begin(); ndigi != mults->end(); ++ndigi) {
0290             _dbxhistos[dbx][ndigi->first]->Fill(ndigi->second);
0291           }
0292         }
0293         if (_dbxhistos.find(-1) != _dbxhistos.end()) {
0294           for (std::map<unsigned int, int>::const_iterator ndigi = mults->begin(); ndigi != mults->end(); ++ndigi) {
0295             _dbxhistos[-1][ndigi->first]->Fill(ndigi->second);
0296           }
0297         }
0298         /*
0299     if(_dbxtkhistos.find(dbx)!=_dbxtkhistos.end()) {
0300       _dbxtkhistos[dbx]->Fill(ntkdigi);
0301     }
0302     if(_dbxtkhistos.find(-1)!=_dbxtkhistos.end()) {
0303       _dbxtkhistos[-1]->Fill(ntkdigi);
0304     }
0305 
0306     if(_dbxtibhistos.find(dbx)!=_dbxtibhistos.end()) {
0307       _dbxtibhistos[dbx]->Fill(ntibdigi);
0308     }
0309     if(_dbxtibhistos.find(-1)!=_dbxtibhistos.end()) {
0310       _dbxtibhistos[-1]->Fill(ntibdigi);
0311     }
0312 
0313     if(_dbxtidhistos.find(dbx)!=_dbxtidhistos.end()) {
0314       _dbxtidhistos[dbx]->Fill(ntiddigi);
0315     }
0316     if(_dbxtidhistos.find(-1)!=_dbxtidhistos.end()) {
0317       _dbxtidhistos[-1]->Fill(ntiddigi);
0318     }
0319     if(_dbxtobhistos.find(dbx)!=_dbxtobhistos.end()) {
0320       _dbxtobhistos[dbx]->Fill(ntobdigi);
0321     }
0322     if(_dbxtobhistos.find(-1)!=_dbxtobhistos.end()) {
0323       _dbxtobhistos[-1]->Fill(ntobdigi);
0324     }
0325     if(_dbxtechistos.find(dbx)!=_dbxtechistos.end()) {
0326       _dbxtechistos[dbx]->Fill(ntecdigi);
0327     }
0328     if(_dbxtechistos.find(-1)!=_dbxtechistos.end()) {
0329       _dbxtechistos[-1]->Fill(ntecdigi);
0330     }
0331     */
0332       }
0333     }
0334   }
0335 }
0336 
0337 // ------------ method called once each job just before starting event loop  ------------
0338 void MultiplicityTimeCorrelations::beginJob() {
0339   LogDebug("IntegerDebug") << " int max and min " << INT_MIN << " " << INT_MAX;
0340   LogDebug("IntegerDebug") << " uint max and min " << UINT_MAX;
0341   LogDebug("IntegerDebug") << " long max and min " << LONG_MIN << " " << LONG_MAX;
0342   LogDebug("IntegerDebug") << " ulong max and min " << ULONG_MAX;
0343   LogDebug("IntegerDebug") << " long long max and min " << LLONG_MIN << " " << LLONG_MAX;
0344   LogDebug("IntegerDebug") << " u long long max and min " << ULLONG_MAX;
0345 
0346   edm::LogInfo("MultiplicityTimeCorrelations")
0347       << " Correlation studies performed only in the orbit # range " << _loworbit << " " << _highorbit;
0348 }
0349 
0350 void MultiplicityTimeCorrelations::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0351   _digibxcorrhmevent.beginRun(iRun.run());
0352 }
0353 // ------------ method called once each job just after ending the event loop  ------------
0354 void MultiplicityTimeCorrelations::endJob() {}
0355 //define this as a plug-in
0356 DEFINE_FWK_MODULE(MultiplicityTimeCorrelations);