Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:09

0001 /** \file
0002  *
0003  * \author Stefano Lacaprara - INFN Legnaro <stefano.lacaprara@pd.infn.it>
0004  */
0005 
0006 /* This Class Header */
0007 #include "RecoLocalMuon/DTSegment/src/DTRecSegment2DExtendedProducer.h"
0008 
0009 /* Collaborating Class Header */
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/ESHandle.h"
0012 #include "FWCore/Framework/interface/ConsumesCollector.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/Utilities/interface/InputTag.h"
0015 using namespace edm;
0016 
0017 #include "Geometry/DTGeometry/interface/DTLayer.h"
0018 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0019 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0020 
0021 #include "DataFormats/DTRecHit/interface/DTRecHit1DPair.h"
0022 #include "DataFormats/DTRecHit/interface/DTRecHit1D.h"
0023 #include "DataFormats/DTRecHit/interface/DTSLRecSegment2D.h"
0024 #include "DataFormats/DTRecHit/interface/DTRecSegment2DCollection.h"
0025 #include "DataFormats/DTRecHit/interface/DTRangeMapAccessor.h"
0026 
0027 #include "RecoLocalMuon/DTSegment/src/DTCombinatorialExtendedPatternReco.h"
0028 
0029 /* C++ Headers */
0030 #include <string>
0031 using namespace std;
0032 
0033 /* ====================================================================== */
0034 
0035 /// Constructor
0036 DTRecSegment2DExtendedProducer::DTRecSegment2DExtendedProducer(const edm::ParameterSet& pset)
0037     : dtGeomToken_(esConsumes()) {
0038   // Set verbose output
0039   debug = pset.getUntrackedParameter<bool>("debug");
0040 
0041   // the name of the 1D rec hits collection
0042   recHits1DToken_ = consumes<DTRecHitCollection>(pset.getParameter<InputTag>("recHits1DLabel"));
0043   recClusToken_ = consumes<DTRecClusterCollection>(pset.getParameter<InputTag>("recClusLabel"));
0044 
0045   if (debug)
0046     cout << "[DTRecSegment2DExtendedProducer] Constructor called" << endl;
0047 
0048   produces<DTRecSegment2DCollection>();
0049 
0050   // Get the concrete reconstruction algo from the factory
0051   theAlgo =
0052       new DTCombinatorialExtendedPatternReco(pset.getParameter<ParameterSet>("Reco2DAlgoConfig"), consumesCollector());
0053 }
0054 
0055 /// Destructor
0056 DTRecSegment2DExtendedProducer::~DTRecSegment2DExtendedProducer() {
0057   if (debug)
0058     cout << "[DTRecSegment2DExtendedProducer] Destructor called" << endl;
0059   delete theAlgo;
0060 }
0061 
0062 /* Operations */
0063 void DTRecSegment2DExtendedProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
0064   if (debug)
0065     cout << "[DTRecSegment2DExtendedProducer] produce called" << endl;
0066   // Get the DT Geometry
0067   const DTGeometry& dtGeom = setup.getData(dtGeomToken_);
0068 
0069   theAlgo->setES(setup);
0070 
0071   // Get the 1D rechits from the event
0072   Handle<DTRecHitCollection> allHits;
0073   event.getByToken(recHits1DToken_, allHits);
0074 
0075   // Get the 1D clusters from the event
0076   Handle<DTRecClusterCollection> dtClusters;
0077   event.getByToken(recClusToken_, dtClusters);
0078   theAlgo->setClusters(vector<DTSLRecCluster>(dtClusters->begin(), dtClusters->end()));
0079 
0080   // Create the pointer to the collection which will store the rechits
0081   auto segments = std::make_unique<DTRecSegment2DCollection>();
0082 
0083   // Iterate through all hit collections ordered by LayerId
0084   DTRecHitCollection::id_iterator dtLayerIt;
0085   DTSuperLayerId oldSlId;
0086   for (dtLayerIt = allHits->id_begin(); dtLayerIt != allHits->id_end(); ++dtLayerIt) {
0087     // The layerId
0088     DTLayerId layerId = (*dtLayerIt);
0089     const DTSuperLayerId SLId = layerId.superlayerId();
0090     if (SLId == oldSlId)
0091       continue;  // I'm on the same SL as before
0092     oldSlId = SLId;
0093 
0094     if (debug)
0095       cout << "Reconstructing the 2D segments in " << SLId << endl;
0096 
0097     const DTSuperLayer* sl = dtGeom.superLayer(SLId);
0098 
0099     // Get all the rec hit in the same superLayer in which layerId relies
0100     DTRecHitCollection::range range = allHits->get(DTRangeMapAccessor::layersBySuperLayer(SLId));
0101 
0102     // Fill the vector with the 1D RecHit
0103     vector<DTRecHit1DPair> pairs(range.first, range.second);
0104 
0105     if (debug)
0106       cout << "Number of 1D-RecHit pairs " << pairs.size() << endl;
0107 
0108     //if(debug) cout << "Start the 2D-segments Reco "<< endl;
0109     OwnVector<DTSLRecSegment2D> segs = theAlgo->reconstruct(sl, pairs);
0110     if (debug) {
0111       cout << "Number of Reconstructed segments: " << segs.size() << endl;
0112       copy(segs.begin(), segs.end(), ostream_iterator<DTSLRecSegment2D>(cout, "\n"));
0113     }
0114 
0115     if (!segs.empty())
0116       segments->put(SLId, segs.begin(), segs.end());
0117   }
0118   event.put(std::move(segments));
0119 }