00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021 #include <memory>
00022
00023
00024 #include "TH1F.h"
00025 #include "TProfile.h"
00026
00027 #include "FWCore/Framework/interface/Frameworkfwd.h"
00028 #include "FWCore/Framework/interface/EDAnalyzer.h"
00029
00030 #include "FWCore/Framework/interface/Event.h"
00031 #include "FWCore/Framework/interface/MakerMacros.h"
00032 #include "FWCore/Framework/interface/Run.h"
00033
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035
00036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00037
00038 #include "FWCore/ServiceRegistry/interface/Service.h"
00039 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00040
00041 #include "FWCore/Utilities/interface/InputTag.h"
00042
00043 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00044
00045 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
00046
00047
00048
00049
00050
00051 class APVCyclePhaseMonitor : public edm::EDAnalyzer {
00052 public:
00053 explicit APVCyclePhaseMonitor(const edm::ParameterSet&);
00054 ~APVCyclePhaseMonitor();
00055
00056
00057 private:
00058 virtual void beginJob() ;
00059 virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
00060 virtual void endRun(const edm::Run&, const edm::EventSetup&) override;
00061 virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
00062 virtual void endJob() ;
00063
00064
00065
00066 const edm::InputTag _apvphasecollection;
00067 std::vector<std::string> _selectedparts;
00068 std::vector<std::string> _selectedvectorparts;
00069 const unsigned int m_maxLS;
00070 const unsigned int m_LSfrac;
00071 RunHistogramManager m_rhm;
00072 std::map<std::string,TH1F*> _hphases;
00073 std::map<std::string,TH1F**> _hselectedphases;
00074 std::map<std::string,TH1F**> _hselectedphasesvector;
00075 std::map<std::string,TH1F**> _hselectedphasessize;
00076 std::map<std::string,TProfile*> _hphasevsorbit;
00077 std::map<std::string,TProfile**> _hselectedphasevsorbit;
00078 std::map<std::string,TProfile**> _hselectedphasevectorvsorbit;
00079 unsigned int _nevents;
00080 };
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093 APVCyclePhaseMonitor::APVCyclePhaseMonitor(const edm::ParameterSet& iConfig):
00094 _apvphasecollection(iConfig.getParameter<edm::InputTag>("apvCyclePhaseCollection")),
00095 _selectedparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedPartitions",std::vector<std::string>())),
00096 _selectedvectorparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedVectorPartitions",std::vector<std::string>())),
00097 m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin",125)),
00098 m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction",16)),
00099 _hphases(),_hselectedphases(),_hselectedphasesvector(),_hselectedphasessize(),
00100 _hphasevsorbit(),_hselectedphasevsorbit(),_hselectedphasevectorvsorbit(),
00101 _nevents(0)
00102 {
00103
00104
00105 edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecollection << " used";
00106
00107 for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
00108
00109 char hname[300];
00110
00111 sprintf(hname,"selected_phase_%s",part->c_str());
00112 edm::LogInfo("SelectedTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00113 _hselectedphases[*part] = m_rhm.makeTH1F(hname,hname,70,-0.5,69.5);
00114
00115 sprintf(hname,"selected_phasevsorbit_%s",part->c_str());
00116 edm::LogInfo("SelectedTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
00117 _hselectedphasevsorbit[*part] = m_rhm.makeTProfile(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
00118 }
00119
00120 for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
00121 part!=_selectedvectorparts.end();++part) {
00122
00123 char hname[300];
00124
00125 sprintf(hname,"selected_phase_vector_%s",part->c_str());
00126 edm::LogInfo("SelectedVectTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00127 _hselectedphasesvector[*part] = m_rhm.makeTH1F(hname,hname,70,-0.5,69.5);
00128
00129 sprintf(hname,"selected_phase_size_%s",part->c_str());
00130 edm::LogInfo("SelectedVectSizeTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00131 _hselectedphasessize[*part] = m_rhm.makeTH1F(hname,hname,10,-0.5,9.5);
00132
00133 sprintf(hname,"selected_phasevectorvsorbit_%s",part->c_str());
00134 edm::LogInfo("SelectedVectTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
00135 _hselectedphasevectorvsorbit[*part] = m_rhm.makeTProfile(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
00136 }
00137
00138
00139 }
00140
00141
00142 APVCyclePhaseMonitor::~APVCyclePhaseMonitor()
00143 {
00144
00145
00146
00147
00148 }
00149
00150
00151
00152
00153
00154
00155
00156 void
00157 APVCyclePhaseMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00158 {
00159 using namespace edm;
00160
00161 _nevents++;
00162
00163 edm::Handle<APVCyclePhaseCollection> apvphases;
00164 iEvent.getByLabel(_apvphasecollection,apvphases);
00165
00166
00167
00168 edm::Service<TFileService> tfserv;
00169
00170 for(std::map<std::string,int>::const_iterator phase = apvphases->get().begin(); phase != apvphases->get().end(); ++phase) {
00171
00172 if(_hphases.find(phase->first)==_hphases.end()) {
00173 char dirname[300];
00174 sprintf(dirname,"run_%d",iEvent.run());
00175 TFileDirectory subrun = tfserv->mkdir(dirname);
00176
00177 char hname[300];
00178
00179 sprintf(hname,"phase_%s",phase->first.c_str());
00180 edm::LogInfo("TH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00181 _hphases[phase->first] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
00182 _hphases[phase->first]->GetXaxis()->SetTitle("BX mod 70"); _hphases[phase->first]->GetYaxis()->SetTitle("Events");
00183
00184 sprintf(hname,"phasevsorbit_%s",phase->first.c_str());
00185 edm::LogInfo("TProfileBeingBooked") << "TProfile " << hname << " being booked" ;
00186 _hphasevsorbit[phase->first] = subrun.make<TProfile>(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
00187 _hphasevsorbit[phase->first]->SetBit(TH1::kCanRebin);
00188 _hphasevsorbit[phase->first]->GetXaxis()->SetTitle("time [orbit#]"); _hphasevsorbit[phase->first]->GetYaxis()->SetTitle("Phase");
00189
00190 }
00191 _hphases[phase->first]->Fill(phase->second);
00192 _hphasevsorbit[phase->first]->Fill(iEvent.orbitNumber(),phase->second);
00193 }
00194
00195
00196
00197 for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
00198 if(_hselectedphases.find(*part)!=_hselectedphases.end() && _hselectedphases[*part] && *_hselectedphases[*part]) {
00199 (*_hselectedphases[*part])->Fill(apvphases->getPhase(*part));
00200 }
00201 if(_hselectedphasevsorbit.find(*part)!=_hselectedphasevsorbit.end() && _hselectedphasevsorbit[*part] && *_hselectedphasevsorbit[*part]) {
00202 (*_hselectedphasevsorbit[*part])->Fill(iEvent.orbitNumber(),apvphases->getPhase(*part));
00203 }
00204 }
00205
00206 for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
00207 part!=_selectedvectorparts.end();++part) {
00208
00209 const std::vector<int> phases = apvphases->getPhases(*part);
00210
00211 if(_hselectedphasessize.find(*part)!=_hselectedphasessize.end() && _hselectedphasessize[*part] && *_hselectedphasessize[*part]) {
00212 (*_hselectedphasessize[*part])->Fill(phases.size());
00213 }
00214
00215 for(std::vector<int>::const_iterator phase=phases.begin();phase!=phases.end();++phase) {
00216 if(_hselectedphasesvector.find(*part)!=_hselectedphasesvector.end() && _hselectedphasesvector[*part] && *_hselectedphasesvector[*part]) {
00217 (*_hselectedphasesvector[*part])->Fill(*phase);
00218 }
00219 if(_hselectedphasevectorvsorbit.find(*part)!=_hselectedphasevectorvsorbit.end() &&
00220 _hselectedphasevectorvsorbit[*part] && *_hselectedphasevectorvsorbit[*part]) {
00221 (*_hselectedphasevectorvsorbit[*part])->Fill(iEvent.orbitNumber(),*phase);
00222 }
00223 }
00224
00225 }
00226 }
00227
00228 void
00229 APVCyclePhaseMonitor::beginRun(const edm::Run& iRun, const edm::EventSetup&)
00230 {
00231
00232 _hphases.clear();
00233 _hphasevsorbit.clear();
00234
00235 m_rhm.beginRun(iRun);
00236
00237 for(std::map<std::string,TH1F**>::const_iterator hist=_hselectedphases.begin();hist!=_hselectedphases.end();++hist) {
00238 if(*(hist->second)) {
00239 (*(hist->second))->GetXaxis()->SetTitle("BX mod 70"); (*(hist->second))->GetYaxis()->SetTitle("Events");
00240 }
00241 }
00242 for(std::map<std::string,TProfile**>::const_iterator prof=_hselectedphasevsorbit.begin();prof!=_hselectedphasevsorbit.end();++prof) {
00243 if(*(prof->second)) {
00244 (*(prof->second))->SetBit(TH1::kCanRebin);
00245 (*(prof->second))->GetXaxis()->SetTitle("time [orbit#]");
00246 (*(prof->second))->GetYaxis()->SetTitle("Phase");
00247 }
00248 }
00249 for(std::map<std::string,TH1F**>::const_iterator hist=_hselectedphasesvector.begin();hist!=_hselectedphasesvector.end();++hist) {
00250 if(*(hist->second)) {
00251 (*(hist->second))->GetXaxis()->SetTitle("BX mod 70"); (*(hist->second))->GetYaxis()->SetTitle("Events");
00252 }
00253 }
00254 for(std::map<std::string,TH1F**>::const_iterator hist=_hselectedphasessize.begin();hist!=_hselectedphasessize.end();++hist) {
00255 if(*(hist->second)) {
00256 (*(hist->second))->GetXaxis()->SetTitle("Number of Phases"); (*(hist->second))->GetYaxis()->SetTitle("Events");
00257 }
00258 }
00259 for(std::map<std::string,TProfile**>::const_iterator prof=_hselectedphasevectorvsorbit.begin();prof!=_hselectedphasevectorvsorbit.end();++prof) {
00260 if(*(prof->second)) {
00261 (*(prof->second))->SetBit(TH1::kCanRebin);
00262 (*(prof->second))->GetXaxis()->SetTitle("time [orbit#]");
00263 (*(prof->second))->GetYaxis()->SetTitle("Phase");
00264 }
00265 }
00266
00267 }
00268
00269 void
00270 APVCyclePhaseMonitor::endRun(const edm::Run& iRun, const edm::EventSetup&)
00271 {
00272 }
00273
00274
00275
00276 void
00277 APVCyclePhaseMonitor::beginJob()
00278 {
00279
00280 }
00281
00282
00283 void
00284 APVCyclePhaseMonitor::endJob() {
00285
00286 edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
00287
00288 }
00289
00290
00291
00292 DEFINE_FWK_MODULE(APVCyclePhaseMonitor);