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/ServiceRegistry/interface/Service.h"
00034 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00035
00036 #include "FWCore/Utilities/interface/InputTag.h"
00037
00038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00039 #include "FWCore/Utilities/interface/Exception.h"
00040
00041 #include <map>
00042 #include <vector>
00043 #include <string>
00044
00045 #include "TH1F.h"
00046 #include "TProfile.h"
00047
00048 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
00049 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00050
00051
00052
00053
00054
00055 class APVCyclePhaseProducerFromL1TS : public edm::EDProducer {
00056 public:
00057 explicit APVCyclePhaseProducerFromL1TS(const edm::ParameterSet&);
00058 ~APVCyclePhaseProducerFromL1TS();
00059
00060 private:
00061 virtual void beginJob() ;
00062 virtual void beginRun(edm::Run&, const edm::EventSetup&);
00063 virtual void endRun(edm::Run&, const edm::EventSetup&);
00064 virtual void produce(edm::Event&, const edm::EventSetup&);
00065 virtual void endJob() ;
00066
00067
00068
00069 const edm::InputTag _l1tscollection;
00070 const std::vector<std::string> _defpartnames;
00071 const std::vector<int> _defphases;
00072 const bool _wantHistos;
00073 const bool _useEC0;
00074 const int _magicOffset;
00075
00076 TH1F* _hsize;
00077 TH1F* _hlresync;
00078 TH1F* _hlOC0;
00079 TH1F* _hlTE;
00080 TH1F* _hlEC0;
00081 TH1F* _hlstart;
00082 TH1F* _hlHR;
00083
00084 TH1F* _hdlec0lresync;
00085 TH1F* _hdlresynclHR;
00086
00087 const unsigned int _firstgoodrun;
00088
00089 long long _lastResync;
00090 long long _lastHardReset;
00091 long long _lastStart;
00092 long long _lastEventCounter0;
00093 long long _lastOrbitCounter0;
00094 long long _lastTestEnable;
00095
00096
00097 };
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111 APVCyclePhaseProducerFromL1TS::APVCyclePhaseProducerFromL1TS(const edm::ParameterSet& iConfig):
00112 _l1tscollection(iConfig.getParameter<edm::InputTag>("l1TSCollection")),
00113 _defpartnames(iConfig.getParameter<std::vector<std::string> >("defaultPartitionNames")),
00114 _defphases(iConfig.getParameter<std::vector<int> >("defaultPhases")),
00115 _wantHistos(iConfig.getUntrackedParameter<bool>("wantHistos",false)),
00116 _useEC0(iConfig.getUntrackedParameter<bool>("useEC0",false)),
00117 _magicOffset(iConfig.getUntrackedParameter<int>("magicOffset",8)),
00118 _firstgoodrun(131768),
00119 _lastResync(-1),_lastHardReset(-1),_lastStart(-1),
00120 _lastEventCounter0(-1),_lastOrbitCounter0(-1),_lastTestEnable(-1)
00121 {
00122
00123 produces<APVCyclePhaseCollection,edm::InEvent>();
00124
00125
00126
00127 }
00128
00129
00130 APVCyclePhaseProducerFromL1TS::~APVCyclePhaseProducerFromL1TS()
00131 {
00132
00133
00134
00135
00136 }
00137
00138
00139
00140
00141
00142
00143
00144 void
00145 APVCyclePhaseProducerFromL1TS::beginRun(edm::Run& iRun, const edm::EventSetup& iSetup)
00146
00147 {
00148
00149
00150
00151 if(_wantHistos) {
00152
00153 edm::Service<TFileService> tfserv;
00154
00155 char dirname[300];
00156 sprintf(dirname,"run_%d",iRun.run());
00157 TFileDirectory subrun = tfserv->mkdir(dirname);
00158
00159 _hsize = subrun.make<TH1F>("size","Level1TriggerScalers Collection size",20,-0.5,19.5);
00160
00161 _hlresync = subrun.make<TH1F>("lresync","Orbit of last resync",3600,0.,3600*11223);
00162 _hlresync->GetXaxis()->SetTitle("Orbit"); _hlresync->GetYaxis()->SetTitle("Events");
00163 _hlresync->SetBit(TH1::kCanRebin);
00164
00165 _hlOC0 = subrun.make<TH1F>("lOC0","Orbit of last OC0",3600,0.,3600*11223);
00166 _hlOC0->GetXaxis()->SetTitle("Orbit"); _hlOC0->GetYaxis()->SetTitle("Events");
00167 _hlOC0->SetBit(TH1::kCanRebin);
00168
00169 _hlTE = subrun.make<TH1F>("lTE","Orbit of last TestEnable",3600,0.,3600*11223);
00170 _hlTE->GetXaxis()->SetTitle("Orbit"); _hlTE->GetYaxis()->SetTitle("Events");
00171 _hlTE->SetBit(TH1::kCanRebin);
00172
00173 _hlstart = subrun.make<TH1F>("lstart","Orbit of last Start",3600,0.,3600*11223);
00174 _hlstart->GetXaxis()->SetTitle("Orbit"); _hlstart->GetYaxis()->SetTitle("Events");
00175 _hlstart->SetBit(TH1::kCanRebin);
00176
00177 _hlEC0 = subrun.make<TH1F>("lEC0","Orbit of last EC0",3600,0.,3600*11223);
00178 _hlEC0->GetXaxis()->SetTitle("Orbit"); _hlEC0->GetYaxis()->SetTitle("Events");
00179 _hlEC0->SetBit(TH1::kCanRebin);
00180
00181 _hlHR = subrun.make<TH1F>("lHR","Orbit of last HardReset",3600,0.,3600*11223);
00182 _hlHR->GetXaxis()->SetTitle("Orbit"); _hlHR->GetYaxis()->SetTitle("Events");
00183 _hlHR->SetBit(TH1::kCanRebin);
00184
00185 _hdlec0lresync = subrun.make<TH1F>("dlec0lresync","Orbit difference EC0-Resync",4000,-1999.5,2000.5);
00186 _hdlec0lresync->GetXaxis()->SetTitle("lastEC0-lastResync");
00187
00188 _hdlresynclHR = subrun.make<TH1F>("dlresynclHR","Orbit difference Resync-HR",4000,-1999.5,2000.5);
00189 _hdlresynclHR->GetXaxis()->SetTitle("lastEC0-lastResync");
00190
00191 }
00192
00193 if(iRun.run() < _firstgoodrun) {
00194 LogDebug("UnreliableMissingL1TriggerScalers") <<
00195 "In this run L1TriggerScalers is missing or unreliable for phase determination: default phases will be used";
00196 }
00197
00198 }
00199
00200 void
00201 APVCyclePhaseProducerFromL1TS::endRun(edm::Run&, const edm::EventSetup&)
00202 {
00203
00204
00205 }
00206
00207
00208 void
00209 APVCyclePhaseProducerFromL1TS::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00210
00211 using namespace edm;
00212
00213 std::auto_ptr<APVCyclePhaseCollection> apvphases(new APVCyclePhaseCollection() );
00214
00215
00216 const std::vector<int>& phases = _defphases;
00217 const std::vector<std::string>& partnames = _defpartnames;
00218
00219 int phasechange = 0;
00220
00221 if(iEvent.run() >= _firstgoodrun ) {
00222
00223 Handle<Level1TriggerScalersCollection> l1ts;
00224 iEvent.getByLabel(_l1tscollection,l1ts);
00225
00226 if(_wantHistos) _hsize->Fill(l1ts->size());
00227
00228
00229
00230 long long orbitoffset = 0;
00231
00232 if(l1ts->size()>0) {
00233
00234 if((*l1ts)[0].lastResync()!=0) {
00235 orbitoffset = _useEC0 ? (*l1ts)[0].lastEventCounter0() + _magicOffset : (*l1ts)[0].lastResync() + _magicOffset;
00236 }
00237
00238 if(_wantHistos) {
00239 _hlresync->Fill((*l1ts)[0].lastResync());
00240 _hlOC0->Fill((*l1ts)[0].lastOrbitCounter0());
00241 _hlTE->Fill((*l1ts)[0].lastTestEnable());
00242 _hlstart->Fill((*l1ts)[0].lastStart());
00243 _hlEC0->Fill((*l1ts)[0].lastEventCounter0());
00244 _hlHR->Fill((*l1ts)[0].lastHardReset());
00245 }
00246
00247 if(_lastResync != (*l1ts)[0].lastResync()) {
00248 _lastResync = (*l1ts)[0].lastResync();
00249 if(_wantHistos) _hdlec0lresync->Fill((*l1ts)[0].lastEventCounter0()-(*l1ts)[0].lastResync());
00250 LogDebug("TTCSignalReceived") << "New Resync at orbit " << _lastResync ;
00251 }
00252 if(_lastHardReset != (*l1ts)[0].lastHardReset()) {
00253 _lastHardReset = (*l1ts)[0].lastHardReset();
00254 if(_wantHistos) _hdlresynclHR->Fill((*l1ts)[0].lastResync()-(*l1ts)[0].lastHardReset());
00255 LogDebug("TTCSignalReceived") << "New HardReset at orbit " << _lastHardReset ;
00256 }
00257 if(_lastTestEnable != (*l1ts)[0].lastTestEnable()) {
00258 _lastTestEnable = (*l1ts)[0].lastTestEnable();
00259
00260 }
00261 if(_lastOrbitCounter0 != (*l1ts)[0].lastOrbitCounter0()) {
00262 _lastOrbitCounter0 = (*l1ts)[0].lastOrbitCounter0();
00263 LogDebug("TTCSignalReceived") << "New OrbitCounter0 at orbit " << _lastOrbitCounter0 ;
00264 }
00265 if(_lastEventCounter0 != (*l1ts)[0].lastEventCounter0()) {
00266 _lastEventCounter0 = (*l1ts)[0].lastEventCounter0();
00267 LogDebug("TTCSignalReceived") << "New EventCounter0 at orbit " << _lastEventCounter0 ;
00268 }
00269 if(_lastStart != (*l1ts)[0].lastStart()) {
00270 _lastStart = (*l1ts)[0].lastStart();
00271 LogDebug("TTCSignalReceived") << "New Start at orbit " << _lastStart ;
00272 }
00273
00274 phasechange = ((long long)(orbitoffset*3564))%70;
00275
00276 }
00277 }
00278
00279
00280 if(phases.size() < partnames.size() ) {
00281
00282 throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent phases/partitions vector sizes: "
00283 << phases.size() << " "
00284 << partnames.size();
00285 }
00286
00287 for(unsigned int ipart=0;ipart<partnames.size();++ipart) {
00288 if(phases[ipart]>=0) {
00289 apvphases->get()[partnames[ipart]] = (phases[ipart]+phasechange)%70;
00290
00291 }
00292 }
00293
00294
00295 iEvent.put(apvphases);
00296
00297 }
00298
00299
00300 void
00301 APVCyclePhaseProducerFromL1TS::beginJob()
00302 {
00303 }
00304
00305
00306 void
00307 APVCyclePhaseProducerFromL1TS::endJob() {
00308 }
00309
00310
00311 DEFINE_FWK_MODULE(APVCyclePhaseProducerFromL1TS);