00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021 #include <memory>
00022
00023
00024 #include <string>
00025 #include <vector>
00026
00027 #include "TH1F.h"
00028 #include "TH2F.h"
00029
00030 #include "FWCore/Framework/interface/Frameworkfwd.h"
00031 #include "FWCore/Framework/interface/EDAnalyzer.h"
00032
00033 #include "FWCore/Framework/interface/Event.h"
00034 #include "FWCore/Framework/interface/MakerMacros.h"
00035 #include "FWCore/Framework/interface/Run.h"
00036
00037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00038
00039 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00040
00041 #include "FWCore/ServiceRegistry/interface/Service.h"
00042 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00043
00044 #include "FWCore/Utilities/interface/InputTag.h"
00045
00046 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
00047 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00048
00049 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
00050
00051
00052
00053
00054 class EventTimeDistribution : public edm::EDAnalyzer {
00055 public:
00056 explicit EventTimeDistribution(const edm::ParameterSet&);
00057 ~EventTimeDistribution();
00058
00059
00060 private:
00061 virtual void beginJob() ;
00062 virtual void beginRun(const edm::Run&, const edm::EventSetup&) ;
00063 virtual void endRun(const edm::Run&, const edm::EventSetup&) ;
00064 virtual void analyze(const edm::Event&, const edm::EventSetup&);
00065 virtual void endJob() ;
00066
00067
00068
00069 const edm::InputTag _historyProduct;
00070 const edm::InputTag _apvphasecoll;
00071 const std::string _phasepart;
00072 const bool _wantdbxvsbxincycle;
00073 const bool _wantdbxvsbx;
00074 const bool _wantbxincyclevsbx;
00075 const bool _wantorbitvsbxincycle;
00076 unsigned int _nevents;
00077 const unsigned int m_maxLS;
00078 const unsigned int m_LSfrac;
00079 const bool m_ewhdepthHisto;
00080
00081
00082
00083 RunHistogramManager _rhm;
00084
00085 TH1F** _dbx;
00086 std::vector<TH1F**> m_dbxhistos;
00087 std::vector<std::pair<unsigned int,unsigned int> > m_dbxindices;
00088 TH1F** _bx;
00089 TH1F** _bxincycle;
00090 TH1F** _orbit;
00091 TH2F** _dbxvsbxincycle;
00092 TH2F** _dbxvsbx;
00093 TH2F** _bxincyclevsbx;
00094 TH2F** _orbitvsbxincycle;
00095 TH1F** m_ewhdepth;
00096
00097 };
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110 EventTimeDistribution::EventTimeDistribution(const edm::ParameterSet& iConfig):
00111 _historyProduct(iConfig.getParameter<edm::InputTag>("historyProduct")),
00112 _apvphasecoll(iConfig.getParameter<edm::InputTag>("apvPhaseCollection")),
00113 _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition","None")),
00114 _wantdbxvsbxincycle(iConfig.getUntrackedParameter<bool>("wantDBXvsBXincycle",false)),
00115 _wantdbxvsbx(iConfig.getUntrackedParameter<bool>("wantDBXvsBX",false)),
00116 _wantbxincyclevsbx(iConfig.getUntrackedParameter<bool>("wantBXincyclevsBX",false)),
00117 _wantorbitvsbxincycle(iConfig.getUntrackedParameter<bool>("wantOrbitvsBXincycle",false)),
00118 _nevents(0),
00119 m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin",100)),
00120 m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction",4)),
00121 m_ewhdepthHisto(iConfig.getUntrackedParameter<bool>("wantEWHDepthHisto",false)),
00122 _rhm(),
00123 _dbxvsbxincycle(0), _dbxvsbx(0), _bxincyclevsbx(0), _orbitvsbxincycle(0), m_ewhdepth(0)
00124 {
00125
00126
00127 std::vector<edm::ParameterSet> dbxhistoparams(iConfig.getUntrackedParameter<std::vector<edm::ParameterSet> >("dbxHistosParams",std::vector<edm::ParameterSet>()));
00128
00129 for(std::vector<edm::ParameterSet>::const_iterator params=dbxhistoparams.begin();params!=dbxhistoparams.end();++params) {
00130 m_dbxindices.push_back(std::pair<unsigned int,unsigned int>(params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent")));
00131 char hname[300];
00132 sprintf(hname,"dbx_%d_%d",params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent"));
00133 char htitle[300];
00134 sprintf(htitle,"dbx(%d,%d)",params->getParameter<unsigned int>("firstEvent"),params->getParameter<unsigned int>("secondEvent"));
00135
00136 m_dbxhistos.push_back(_rhm.makeTH1F(hname,htitle,params->getParameter<int>("nbins"),params->getParameter<double>("min"),
00137 params->getParameter<double>("max")));
00138 LogDebug("DBXHistoPreBooking") << "Booked DBX histo named " << hname << " untitled " << htitle;
00139 }
00140
00141
00142
00143 _dbx = _rhm.makeTH1F("dbx","dbx",1000,-0.5,999.5);
00144 _bx = _rhm.makeTH1F("bx","BX number",3564,-0.5,3563.5);
00145 _bxincycle = _rhm.makeTH1F("bxcycle","bxcycle",70,-0.5,69.5);
00146 _orbit = _rhm.makeTH1F("orbit","orbit",m_LSfrac*m_maxLS,0,m_maxLS*262144);
00147 if(_wantdbxvsbxincycle) _dbxvsbxincycle = _rhm.makeTH2F("dbxvsbxincycle","dbxvsbxincycle",70,-0.5,69.5,1000,-0.5,999.5);
00148 if(_wantdbxvsbx) _dbxvsbx = _rhm.makeTH2F("dbxvsbx","dbxvsbx",3564,-0.5,3563.5,1000,-0.5,999.5);
00149 if(_wantbxincyclevsbx) _bxincyclevsbx = _rhm.makeTH2F("bxincyclevsbx","bxincyclevsbx",3564,-0.5,3563.5,70,-0.5,69.5);
00150 if(_wantorbitvsbxincycle) _orbitvsbxincycle = _rhm.makeTH2F("orbitvsbxincycle","orbitvsbxincycle",70,-0.5,69.5,m_maxLS,0,m_maxLS*262144);
00151 if(m_ewhdepthHisto) m_ewhdepth = _rhm.makeTH1F("ewhdepth","EventWithHistory Depth",11,-0.5,10.5);
00152
00153 edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecoll << " used";
00154
00155 }
00156
00157
00158 EventTimeDistribution::~EventTimeDistribution()
00159 {
00160
00161
00162
00163
00164 }
00165
00166
00167
00168
00169
00170
00171
00172 void
00173 EventTimeDistribution::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00174 {
00175 using namespace edm;
00176
00177 _nevents++;
00178
00179 edm::Handle<EventWithHistory> he;
00180 iEvent.getByLabel(_historyProduct,he);
00181
00182
00183
00184 (*_dbx)->Fill(he->deltaBX());
00185 std::vector<std::pair<unsigned int,unsigned int> >::const_iterator indices=m_dbxindices.begin();
00186 for(std::vector<TH1F**>::const_iterator dbxhist=m_dbxhistos.begin();dbxhist!=m_dbxhistos.end();++dbxhist,++indices) {
00187 (*(*dbxhist))->Fill(he->deltaBX(indices->first,indices->second));
00188 }
00189
00190 (*_bx)->Fill(iEvent.bunchCrossing());
00191 (*_orbit)->Fill(iEvent.orbitNumber());
00192 if(_dbxvsbx && *_dbxvsbx) (*_dbxvsbx)->Fill(iEvent.bunchCrossing(),he->deltaBX());
00193 if(m_ewhdepth && *m_ewhdepth) (*m_ewhdepth)->Fill(he->depth());
00194
00195 edm::Handle<APVCyclePhaseCollection> apvphase;
00196 iEvent.getByLabel(_apvphasecoll,apvphase);
00197
00198 long long tbx = he->absoluteBX();
00199 if(apvphase.isValid() && !apvphase.failedToGet()) {
00200 const int thephase = apvphase->getPhase(_phasepart);
00201 if(thephase!=APVCyclePhaseCollection::invalid &&
00202 thephase!=APVCyclePhaseCollection::multiphase &&
00203 thephase!=APVCyclePhaseCollection::nopartition) {
00204
00205 tbx -= thephase;
00206 (*_bxincycle)->Fill(tbx%70);
00207 if(_dbxvsbxincycle && *_dbxvsbxincycle) (*_dbxvsbxincycle)->Fill(tbx%70,he->deltaBX());
00208 if(_bxincyclevsbx && *_bxincyclevsbx) (*_bxincyclevsbx)->Fill(iEvent.bunchCrossing(),tbx%70);
00209 if(_orbitvsbxincycle && *_orbitvsbxincycle) (*_orbitvsbxincycle)->Fill(tbx%70,iEvent.orbitNumber());
00210
00211 }
00212 else {
00213 LogDebug("InvalidPhase") << "Invalid APVCyclePhase value : " << _phasepart << " " << thephase;
00214 }
00215 }
00216 }
00217
00218 void
00219 EventTimeDistribution::beginRun(const edm::Run& iRun, const edm::EventSetup&)
00220 {
00221
00222 _rhm.beginRun(iRun);
00223
00224 if(*_dbx) { (*_dbx)->GetXaxis()->SetTitle("#DeltaBX"); }
00225
00226 LogDebug("NomberOfHistos") << m_dbxhistos.size();
00227 for(std::vector<TH1F**>::const_iterator dbxhist=m_dbxhistos.begin();dbxhist!=m_dbxhistos.end();++dbxhist) {
00228 LogDebug("HistoPointer") << *dbxhist;
00229 if(*(*dbxhist)) { (*(*dbxhist))->GetXaxis()->SetTitle("#DeltaBX"); }
00230 }
00231 LogDebug("LabelDone") << "all labels set";
00232
00233 if(*_bx) { (*_bx)->GetXaxis()->SetTitle("BX"); }
00234
00235 if(*_bxincycle) { (*_bxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); }
00236
00237 if(*_orbit) {
00238 (*_orbit)->SetBit(TH1::kCanRebin);
00239 (*_orbit)->GetXaxis()->SetTitle("time [Orb#]");
00240 }
00241
00242 LogDebug("StdPlotsDone") << "all labels in std plots set";
00243
00244 if(_dbxvsbxincycle && *_dbxvsbxincycle) {
00245 (*_dbxvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_dbxvsbxincycle)->GetYaxis()->SetTitle("#DeltaBX");
00246 }
00247
00248 if(_dbxvsbx && *_dbxvsbx) { (*_dbxvsbx)->GetXaxis()->SetTitle("BX"); (*_dbxvsbx)->GetYaxis()->SetTitle("#DeltaBX"); }
00249
00250 if(_bxincyclevsbx && *_bxincyclevsbx) {
00251 (*_bxincyclevsbx)->GetXaxis()->SetTitle("BX"); (*_bxincyclevsbx)->GetYaxis()->SetTitle("Event BX mod(70)");
00252 }
00253
00254 if(_orbitvsbxincycle && *_orbitvsbxincycle) {
00255 (*_orbitvsbxincycle)->SetBit(TH1::kCanRebin);
00256 (*_orbitvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_orbitvsbxincycle)->GetYaxis()->SetTitle("time [Orb#]");
00257 }
00258
00259 if(m_ewhdepth && *m_ewhdepth) {
00260 (*m_ewhdepth)->GetXaxis()->SetTitle("Depth");
00261 }
00262
00263 }
00264
00265 void
00266 EventTimeDistribution::endRun(const edm::Run& iRun, const edm::EventSetup&)
00267 {
00268 }
00269
00270
00271
00272 void
00273 EventTimeDistribution::beginJob()
00274 {
00275
00276 }
00277
00278
00279 void
00280 EventTimeDistribution::endJob() {
00281
00282 edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
00283
00284 }
00285
00286
00287
00288 DEFINE_FWK_MODULE(EventTimeDistribution);