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