-
Notifications
You must be signed in to change notification settings - Fork 4.2k
/
APVCyclePhaseProducerFromL1ABC.cc
297 lines (213 loc) · 8.26 KB
/
APVCyclePhaseProducerFromL1ABC.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
// -*- C++ -*-
//
// Package: SiStripTools
// Class: APVCyclePhaseProducerFromL1ABC
//
/**\class APVCyclePhaseProducerFromL1ABC APVCyclePhaseProducerFromL1ABC.cc DPGAnalysis/SiStripTools/plugins/APVCyclePhaseProducerFromL1ABC.cc
Description: EDproducer for APVCyclePhaseCollection which uses the configuration file to assign a phase to the run
Implementation:
<Notes on implementation>
*/
//
// Original Author: Andrea Venturi
// Created: Mon Jan 12 09:05:45 CET 2009
//
//
// system include files
#include <memory>
// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Run.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/Exception.h"
#include <map>
#include <vector>
#include <string>
#include "TH1F.h"
#include "TProfile.h"
#include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
#include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
#include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
//
// class decleration
//
class APVCyclePhaseProducerFromL1ABC : public edm::EDProducer {
public:
explicit APVCyclePhaseProducerFromL1ABC(const edm::ParameterSet&);
~APVCyclePhaseProducerFromL1ABC();
private:
virtual void beginJob() override ;
virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
virtual void endRun(const edm::Run&, const edm::EventSetup&) override;
virtual void produce(edm::Event&, const edm::EventSetup&) override;
virtual void endJob() override ;
// ----------member data ---------------------------
edm::EDGetTokenT<L1AcceptBunchCrossingCollection> _l1abccollectionToken;
const std::vector<std::string> _defpartnames;
const std::vector<int> _defphases;
const int _orbitoffsetSOR;
const bool _wantHistos;
RunHistogramManager m_rhm;
TH1F** _hbx;
TH1F** _hdbx;
TH1F** _hdorbit;
const unsigned int _firstgoodrun;
std::map<edm::EventNumber_t, long long> _offsets;
long long _curroffset;
edm::EventNumber_t _curroffevent;
};
//
// constants, enums and typedefs
//
//
// static data member definitions
//
//
// constructors and destructor
//
APVCyclePhaseProducerFromL1ABC::APVCyclePhaseProducerFromL1ABC(const edm::ParameterSet& iConfig):
_l1abccollectionToken(mayConsume<L1AcceptBunchCrossingCollection>(iConfig.getParameter<edm::InputTag>("l1ABCCollection"))),
_defpartnames(iConfig.getParameter<std::vector<std::string> >("defaultPartitionNames")),
_defphases(iConfig.getParameter<std::vector<int> >("defaultPhases")),
_orbitoffsetSOR(iConfig.getParameter<int>("StartOfRunOrbitOffset")),
_wantHistos(iConfig.getUntrackedParameter<bool>("wantHistos",false)),
m_rhm(consumesCollector()),
_hbx(0),_hdbx(0),_hdorbit(0),_firstgoodrun(110878),
_offsets(), _curroffset(0), _curroffevent(0)
{
produces<APVCyclePhaseCollection,edm::InEvent>();
//now do what ever other initialization is needed
if(_wantHistos) {
_hbx = m_rhm.makeTH1F("l1abcbx","BX number from L1ABC",4096,-0.5,4095.5);
_hdbx = m_rhm.makeTH1F("dbx","BX number difference",4096*2-1,-4095.5,4095.5);
_hdorbit = m_rhm.makeTH1F("dorbit","Orbit Number difference",9999,-4999.5,4999.5);
}
}
APVCyclePhaseProducerFromL1ABC::~APVCyclePhaseProducerFromL1ABC()
{
// do anything here that needs to be done at desctruction time
// (e.g. close files, deallocate resources etc.)
}
//
// member functions
//
// ------------ method called to produce the data ------------
void
APVCyclePhaseProducerFromL1ABC::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
{
// reset offset vector
_offsets.clear();
edm::LogInfo("AbsoluteBXOffsetReset") << "Absolute BX offset map reset";
if(_wantHistos) {
m_rhm.beginRun(iRun);
if(_hbx && *_hbx) {
(*_hbx)->GetXaxis()->SetTitle("BX"); (*_hbx)->GetYaxis()->SetTitle("Events");
}
if(_hdbx && *_hdbx) {
(*_hdbx)->GetXaxis()->SetTitle("#DeltaBX"); (*_hdbx)->GetYaxis()->SetTitle("Events");
}
if(_hdorbit && *_hdorbit) {
(*_hdorbit)->GetXaxis()->SetTitle("#Deltaorbit"); (*_hdorbit)->GetYaxis()->SetTitle("Events");
}
}
if(iRun.run() < _firstgoodrun) {
edm::LogInfo("UnreliableMissingL1AcceptBunchCrossingCollection") <<
"In this run L1AcceptBunchCrossingCollection is missing or unreliable: default phases will be used";
}
}
void
APVCyclePhaseProducerFromL1ABC::endRun(const edm::Run&, const edm::EventSetup&)
{
// summary of absolute bx offset vector
edm::LogInfo("L1AcceptBunchCrossingAbsoluteBXOffsetSummary") << "Absolute BX offset summary:";
for(std::map<edm::EventNumber_t, long long>::const_iterator offset=_offsets.begin();offset!=_offsets.end();++offset) {
edm::LogVerbatim("L1AcceptBunchCrossingAbsoluteBXOffsetSummary") << offset->first << " " << offset->second;
}
}
void
APVCyclePhaseProducerFromL1ABC::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
using namespace edm;
std::unique_ptr<APVCyclePhaseCollection> apvphases(new APVCyclePhaseCollection() );
const std::vector<int>& phases = _defphases;
const std::vector<std::string>& partnames = _defpartnames;
// Look for the L1AcceptBunchCrossingCollection
int phasechange = 0;
if(iEvent.run() >= _firstgoodrun ) {
Handle<L1AcceptBunchCrossingCollection > pIn;
iEvent.getByToken(_l1abccollectionToken,pIn);
// offset computation
long long orbitoffset = _orbitoffsetSOR;
int bxoffset = 0;
for(L1AcceptBunchCrossingCollection::const_iterator l1abc=pIn->begin();l1abc!=pIn->end();++l1abc) {
if(l1abc->l1AcceptOffset()==0) {
if(l1abc->eventType()!=0) {
orbitoffset = (long long)iEvent.orbitNumber() - (long long)l1abc->orbitNumber() ;
bxoffset = iEvent.bunchCrossing() - l1abc->bunchCrossing();
if(_wantHistos) {
if(_hbx && *_hbx) (*_hbx)->Fill(l1abc->bunchCrossing());
if(_hdbx && *_hdbx) (*_hdbx)->Fill(bxoffset);
if(_hdorbit && *_hdorbit) (*_hdorbit)->Fill(orbitoffset);
}
}
else {
edm::LogWarning("L1AcceptBunchCrossingNoType") << "L1AcceptBunchCrossing with no type found: ";
for(L1AcceptBunchCrossingCollection::const_iterator debu=pIn->begin();debu!=pIn->end();++debu) {
edm::LogPrint("L1AcceptBunchCrossingNoType") << *debu;
}
}
}
}
long long absbxoffset = orbitoffset*3564 + bxoffset;
if(orbitoffset != _orbitoffsetSOR) phasechange = (orbitoffset*3564)%70;
if(_offsets.size()==0) {
_curroffset = absbxoffset;
_curroffevent = iEvent.id().event();
_offsets[iEvent.id().event()] = absbxoffset;
}
else {
if(_curroffset != absbxoffset || iEvent.id().event() < _curroffevent ) {
if( _curroffset != absbxoffset) {
edm::LogInfo("L1AcceptBunchCrossingAbsoluteBXOffsetChanged") << "Absolute BX offset changed from "
<< _curroffset << " to "
<< absbxoffset << " at orbit "
<< iEvent.orbitNumber() << " and BX "
<< iEvent.bunchCrossing();
for(L1AcceptBunchCrossingCollection::const_iterator l1abc=pIn->begin();l1abc!=pIn->end();++l1abc) {
edm::LogVerbatim("L1AcceptBunchCrossingAbsoluteBXOffsetChanged") << *l1abc;
}
}
_curroffset = absbxoffset;
_curroffevent = iEvent.id().event();
_offsets[iEvent.id().event()] = absbxoffset;
}
}
}
if(phases.size() < partnames.size() ) {
// throw exception
throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent phases/partitions vector sizes: "
<< phases.size() << " "
<< partnames.size();
}
for(unsigned int ipart=0;ipart<partnames.size();++ipart) {
if(phases[ipart]>=0) {
apvphases->get()[partnames[ipart]] = (phases[ipart]+phasechange)%70;
}
}
iEvent.put(std::move(apvphases));
}
// ------------ method called once each job just before starting event loop ------------
void
APVCyclePhaseProducerFromL1ABC::beginJob()
{
}
// ------------ method called once each job just after ending the event loop ------------
void
APVCyclePhaseProducerFromL1ABC::endJob() {
}
//define this as a plug-in
DEFINE_FWK_MODULE(APVCyclePhaseProducerFromL1ABC);