File indexing completed on 2023-03-17 11:19:20
0001
0002
0003
0004
0005
0006
0007 #include "RecoLocalMuon/DTSegment/src/DTRecSegment4DProducer.h"
0008
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/Utilities/interface/InputTag.h"
0014
0015 #include "RecoLocalMuon/DTSegment/src/DTSegmentUpdator.h"
0016 #include "RecoLocalMuon/DTSegment/src/DTRecSegment4DAlgoFactory.h"
0017
0018 #include "DataFormats/Common/interface/OwnVector.h"
0019 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
0020
0021 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0022
0023 using namespace edm;
0024 using namespace std;
0025
0026 DTRecSegment4DProducer::DTRecSegment4DProducer(const ParameterSet& pset)
0027 :
0028 the4DAlgo{DTRecSegment4DAlgoFactory::get()->create(pset.getParameter<string>("Reco4DAlgoName"),
0029 pset.getParameter<ParameterSet>("Reco4DAlgoConfig"),
0030 consumesCollector())} {
0031 produces<DTRecSegment4DCollection>();
0032
0033
0034 debug = pset.getUntrackedParameter<bool>("debug", false);
0035
0036 if (debug)
0037 cout << "[DTRecSegment4DProducer] Constructor called" << endl;
0038
0039
0040 recHits1DToken_ = consumes<DTRecHitCollection>(pset.getParameter<InputTag>("recHits1DLabel"));
0041
0042
0043 recHits2DToken_ = consumes<DTRecSegment2DCollection>(pset.getParameter<InputTag>("recHits2DLabel"));
0044
0045 if (debug)
0046 cout << "the Reco4D AlgoName is " << pset.getParameter<string>("Reco4DAlgoName") << endl;
0047 }
0048
0049
0050 DTRecSegment4DProducer::~DTRecSegment4DProducer() {
0051 if (debug)
0052 cout << "[DTRecSegment4DProducer] Destructor called" << endl;
0053 }
0054
0055 void DTRecSegment4DProducer::produce(Event& event, const EventSetup& setup) {
0056
0057 Handle<DTRecHitCollection> all1DHits;
0058 event.getByToken(recHits1DToken_, all1DHits);
0059
0060
0061 Handle<DTRecSegment2DCollection> all2DSegments;
0062 if (the4DAlgo->wants2DSegments())
0063 event.getByToken(recHits2DToken_, all2DSegments);
0064
0065
0066 auto segments4DCollection = std::make_unique<DTRecSegment4DCollection>();
0067
0068
0069 the4DAlgo->setES(setup);
0070
0071
0072 DTRecHitCollection::id_iterator dtLayerIt;
0073
0074 DTChamberId oldChId;
0075
0076 for (dtLayerIt = all1DHits->id_begin(); dtLayerIt != all1DHits->id_end(); ++dtLayerIt) {
0077
0078 const DTChamberId chId = (*dtLayerIt).chamberId();
0079 if (chId == oldChId)
0080 continue;
0081 oldChId = chId;
0082 if (debug)
0083 cout << "ChamberId: " << chId << endl;
0084 the4DAlgo->setChamber(chId);
0085
0086 if (debug)
0087 cout << "Take the DTRecHits1D and set them in the reconstructor" << endl;
0088
0089 the4DAlgo->setDTRecHit1DContainer(all1DHits);
0090
0091 if (debug)
0092 cout << "Take the DTRecSegments2D and set them in the reconstructor" << endl;
0093
0094 the4DAlgo->setDTRecSegment2DContainer(all2DSegments);
0095
0096 if (debug)
0097 cout << "Start 4D-Segments Reco " << endl;
0098
0099 OwnVector<DTRecSegment4D> segments4D = the4DAlgo->reconstruct();
0100
0101 if (debug) {
0102 cout << "Number of reconstructed 4D-segments " << segments4D.size() << endl;
0103 copy(segments4D.begin(), segments4D.end(), ostream_iterator<DTRecSegment4D>(cout, "\n"));
0104 }
0105
0106 if (!segments4D.empty())
0107
0108 segments4DCollection->put(chId, segments4D.begin(), segments4D.end());
0109 }
0110
0111 event.put(std::move(segments4DCollection));
0112 }