File indexing completed on 2024-04-06 12:19:46
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDProducer.h"
0025
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 #include "FWCore/Framework/interface/ConsumesCollector.h"
0029
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 #include "FWCore/Utilities/interface/StreamID.h"
0032 #include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"
0033 #include "FWCore/Framework/interface/ConsumesCollector.h"
0034 #include "DataFormats/Common/interface/Handle.h"
0035
0036 #include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"
0037
0038 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0039 #include "FWCore/Framework/interface/EventSetup.h"
0040 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
0041 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
0042 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0043 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0044 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0045
0046 #include "L1Trigger/DTTriggerPhase2/interface/constants.h"
0047
0048 namespace edm {
0049 class ParameterSet;
0050 class EventSetup;
0051 }
0052
0053 using namespace std;
0054 using namespace edm;
0055 using namespace cmsdt;
0056
0057
0058
0059
0060 class CalibratedDigis : public edm::stream::EDProducer<> {
0061 public:
0062 explicit CalibratedDigis(const edm::ParameterSet&);
0063 ~CalibratedDigis() override;
0064
0065 private:
0066 int timeOffset_;
0067 int flat_calib_;
0068 int scenario;
0069
0070 void produce(edm::Event&, const edm::EventSetup&) override;
0071
0072 std::unique_ptr<DTTTrigBaseSync> theSync;
0073
0074
0075 edm::EDGetTokenT<DTDigiCollection> dtDigisToken;
0076 edm::Handle<DTDigiCollection> DTDigiHandle;
0077 edm::InputTag dtDigiTag;
0078
0079 static constexpr float bxspacing = 25.0;
0080 static constexpr float timeshift = 400.0;
0081 static constexpr float flatcalib = 325.0;
0082 };
0083
0084
0085
0086
0087 CalibratedDigis::CalibratedDigis(const edm::ParameterSet& iConfig) {
0088
0089 dtDigiTag = iConfig.getParameter<InputTag>("dtDigiTag");
0090 dtDigisToken = consumes<DTDigiCollection>(dtDigiTag);
0091
0092 theSync = DTTTrigSyncFactory::get()->create(iConfig.getParameter<string>("tTrigMode"),
0093 iConfig.getParameter<ParameterSet>("tTrigModeConfig"),
0094 consumesCollector());
0095
0096 flat_calib_ = iConfig.getParameter<int>("flat_calib");
0097 timeOffset_ = iConfig.getParameter<int>("timeOffset");
0098
0099 scenario = iConfig.getParameter<int>("scenario");
0100
0101 produces<DTDigiCollection>();
0102
0103 }
0104
0105 CalibratedDigis::~CalibratedDigis() {
0106
0107
0108 }
0109
0110
0111
0112
0113
0114
0115 void CalibratedDigis::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0116 using namespace edm;
0117 theSync->setES(iSetup);
0118 iEvent.getByToken(dtDigisToken, DTDigiHandle);
0119 DTDigiCollection mydigis;
0120
0121 for (const auto& dtLayerIt : *DTDigiHandle) {
0122 const DTLayerId& layerId = dtLayerIt.first;
0123 for (DTDigiCollection::const_iterator digiIt = dtLayerIt.second.first; digiIt != dtLayerIt.second.second;
0124 ++digiIt) {
0125 DTWireId wireId(layerId, (*digiIt).wire());
0126 float digiTime = (*digiIt).time();
0127 int wire = (*digiIt).wire();
0128 int number = (*digiIt).number();
0129 float newTime = 0;
0130 if (flat_calib_ != 0)
0131 newTime = digiTime - flatcalib + bxspacing * iEvent.eventAuxiliary().bunchCrossing() + float(timeOffset_);
0132 else {
0133 if (scenario == MC)
0134 newTime = digiTime + bxspacing * timeshift;
0135 else if (scenario == SLICE_TEST)
0136 newTime = digiTime;
0137 else
0138 newTime = digiTime - theSync->offset(wireId) + bxspacing * iEvent.eventAuxiliary().bunchCrossing() +
0139 float(timeOffset_);
0140 }
0141 DTDigi newDigi(wire, newTime, number);
0142 mydigis.insertDigi(layerId, newDigi);
0143 }
0144 }
0145 auto CorrectedDTDigiCollection = std::make_unique<DTDigiCollection>(mydigis);
0146 iEvent.put(std::move(CorrectedDTDigiCollection));
0147 }
0148
0149
0150 DEFINE_FWK_MODULE(CalibratedDigis);