File indexing completed on 2023-03-17 11:19:17
0001
0002
0003
0004
0005
0006 #include "DTRecHitProducer.h"
0007
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011 #include "FWCore/Framework/interface/ConsumesCollector.h"
0012
0013 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
0014
0015 #include "Geometry/DTGeometry/interface/DTLayer.h"
0016 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0017 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0018 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
0019 #include "DataFormats/DTRecHit/interface/DTRecHit1DPair.h"
0020
0021 #include "RecoLocalMuon/DTRecHit/interface/DTRecHitBaseAlgo.h"
0022 #include "RecoLocalMuon/DTRecHit/interface/DTRecHitAlgoFactory.h"
0023 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
0024 #include <string>
0025
0026 using namespace edm;
0027 using namespace std;
0028
0029 DTRecHitProducer::DTRecHitProducer(const ParameterSet& config)
0030 :
0031 debug(config.getUntrackedParameter<bool>("debug", false)),
0032 dtGeomToken_(esConsumes()),
0033
0034 theAlgo{DTRecHitAlgoFactory::get()->create(config.getParameter<string>("recAlgo"),
0035 config.getParameter<ParameterSet>("recAlgoConfig"),
0036 consumesCollector())} {
0037 if (debug)
0038 cout << "[DTRecHitProducer] Constructor called" << endl;
0039
0040 produces<DTRecHitCollection>();
0041
0042 DTDigiToken_ = consumes<DTDigiCollection>(config.getParameter<InputTag>("dtDigiLabel"));
0043 }
0044
0045 DTRecHitProducer::~DTRecHitProducer() {
0046 if (debug)
0047 cout << "[DTRecHitProducer] Destructor called" << endl;
0048 }
0049
0050 void DTRecHitProducer::produce(Event& event, const EventSetup& setup) {
0051
0052 const DTGeometry& dtGeom = setup.getData(dtGeomToken_);
0053
0054
0055 Handle<DTDigiCollection> digis;
0056 event.getByToken(DTDigiToken_, digis);
0057
0058
0059 theAlgo->setES(setup);
0060
0061
0062 auto recHitCollection = std::make_unique<DTRecHitCollection>();
0063
0064
0065 DTDigiCollection::DigiRangeIterator dtLayerIt;
0066 for (dtLayerIt = digis->begin(); dtLayerIt != digis->end(); ++dtLayerIt) {
0067
0068 const DTLayerId& layerId = (*dtLayerIt).first;
0069
0070 const DTLayer* layer = dtGeom.layer(layerId);
0071
0072
0073 const DTDigiCollection::Range& range = (*dtLayerIt).second;
0074
0075 OwnVector<DTRecHit1DPair> recHits = theAlgo->reconstruct(layer, layerId, range);
0076
0077 if (debug)
0078 cout << "Number of hits in this layer: " << recHits.size() << endl;
0079 if (!recHits.empty())
0080 recHitCollection->put(layerId, recHits.begin(), recHits.end());
0081 }
0082
0083 event.put(std::move(recHitCollection));
0084 }