CMS 3D CMS Logo

CalibratedDigis.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1Trigger/DTTriggerPhase2
4 // Class: CalibratedDigis
5 //
13 //
14 // Original Author: Luigi Guiducci
15 // Created: Fri, 11 Jan 2019 12:49:12 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
25 
29 
35 
37 
45 
47 
48 namespace edm {
49  class ParameterSet;
50  class EventSetup;
51 } // namespace edm
52 
53 using namespace std;
54 using namespace edm;
55 using namespace cmsdt;
56 //
57 // class declaration
58 //
59 
61 public:
62  explicit CalibratedDigis(const edm::ParameterSet&);
63  ~CalibratedDigis() override;
64 
65 private:
68  int scenario;
69 
70  void produce(edm::Event&, const edm::EventSetup&) override;
71 
72  std::unique_ptr<DTTTrigBaseSync> theSync;
73 
74  // ----------member data ---------------------------
78 
79  static constexpr float bxspacing = 25.0;
80  static constexpr float timeshift = 400.0;
81  static constexpr float flatcalib = 325.0;
82 };
83 
84 //
85 // constructors and destructor
86 //
88  //register your products
89  dtDigiTag = iConfig.getParameter<InputTag>("dtDigiTag");
90  dtDigisToken = consumes<DTDigiCollection>(dtDigiTag);
91 
92  theSync = DTTTrigSyncFactory::get()->create(iConfig.getParameter<string>("tTrigMode"),
93  iConfig.getParameter<ParameterSet>("tTrigModeConfig"),
94  consumesCollector());
95 
96  flat_calib_ = iConfig.getParameter<int>("flat_calib");
97  timeOffset_ = iConfig.getParameter<int>("timeOffset");
98 
99  scenario = iConfig.getParameter<int>("scenario");
100 
101  produces<DTDigiCollection>();
102  //now do what ever other initialization is needed
103 }
104 
106  // do anything here that needs to be done at destruction time
107  // (e.g. close files, deallocate resources etc.)
108 }
109 
110 //
111 // member functions
112 //
113 
114 // ------------ method called to produce the data ------------
116  using namespace edm;
117  theSync->setES(iSetup);
118  iEvent.getByToken(dtDigisToken, DTDigiHandle);
119  DTDigiCollection mydigis;
120 
121  for (const auto& dtLayerIt : *DTDigiHandle) {
122  const DTLayerId& layerId = dtLayerIt.first;
123  for (DTDigiCollection::const_iterator digiIt = dtLayerIt.second.first; digiIt != dtLayerIt.second.second;
124  ++digiIt) {
125  DTWireId wireId(layerId, (*digiIt).wire());
126  float digiTime = (*digiIt).time();
127  int wire = (*digiIt).wire();
128  int number = (*digiIt).number();
129  float newTime = 0;
130  if (flat_calib_ != 0)
131  newTime = digiTime - flatcalib + bxspacing * iEvent.eventAuxiliary().bunchCrossing() + float(timeOffset_);
132  else {
133  if (scenario == MC) //FIX MC
134  newTime = digiTime + bxspacing * timeshift;
135  else if (scenario == SLICE_TEST) //FIX SliceTest
136  newTime = digiTime;
137  else
138  newTime = digiTime - theSync->offset(wireId) + bxspacing * iEvent.eventAuxiliary().bunchCrossing() +
139  float(timeOffset_);
140  }
141  DTDigi newDigi(wire, newTime, number);
142  mydigis.insertDigi(layerId, newDigi);
143  }
144  }
145  auto CorrectedDTDigiCollection = std::make_unique<DTDigiCollection>(mydigis);
146  iEvent.put(std::move(CorrectedDTDigiCollection));
147 }
148 
149 //define this as a plug-in
CalibratedDigis(const edm::ParameterSet &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void produce(edm::Event &, const edm::EventSetup &) override
std::unique_ptr< DTTTrigBaseSync > theSync
edm::EDGetTokenT< DTDigiCollection > dtDigisToken
scenario
Definition: constants.h:173
int iEvent
Definition: GenABIO.cc:224
edm::InputTag dtDigiTag
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::Handle< DTDigiCollection > DTDigiHandle
Definition: DTDigi.h:17
std::vector< DigiType >::const_iterator const_iterator
~CalibratedDigis() override
HLT enums.
#define get
def move(src, dest)
Definition: eostools.py:511