CMS 3D CMS Logo

MTDTrackingRecHitProducer.cc
Go to the documentation of this file.
6 
9 
16 
19 
21 
23 public:
24  explicit MTDTrackingRecHitProducer(const edm::ParameterSet& ps);
25  ~MTDTrackingRecHitProducer() override = default;
26  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
27 
28  void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& es) const override;
29 
30 private:
31  const edm::EDGetTokenT<FTLClusterCollection> ftlbClusters_; // collection of barrel digis
32  const edm::EDGetTokenT<FTLClusterCollection> ftleClusters_; // collection of endcap digis
33 
36 };
37 
39  : ftlbClusters_(consumes<FTLClusterCollection>(ps.getParameter<edm::InputTag>("barrelClusters"))),
40  ftleClusters_(consumes<FTLClusterCollection>(ps.getParameter<edm::InputTag>("endcapClusters"))),
42  cpeToken_(esConsumes<MTDClusterParameterEstimator, MTDCPERecord>(edm::ESInputTag("", "MTDCPEBase"))) {
43  produces<MTDTrackingDetSetVector>();
44 }
45 
46 // Configuration descriptions
49  desc.add<edm::InputTag>("barrelClusters", edm::InputTag("mtdClusters:FTLBarrel"));
50  desc.add<edm::InputTag>("endcapClusters", edm::InputTag("mtdClusters:FTLEndcap"));
51  descriptions.add("mtdTrackingRecHitProducer", desc);
52 }
53 
55  auto const& geom = es.getData(mtdgeoToken_);
56 
57  auto const& cpe = es.getData(cpeToken_);
58 
60  evt.getByToken(ftlbClusters_, inputBarrel);
61 
63  evt.getByToken(ftleClusters_, inputEndcap);
64 
65  std::array<edm::Handle<FTLClusterCollection>, 2> inputHandle{{inputBarrel, inputEndcap}};
66 
67  auto outputhits = std::make_unique<MTDTrackingDetSetVector>();
68  auto& theoutputhits = *outputhits;
69 
70  //---------------------------------------------------------------------------
73  //---------------------------------------------------------------------------
74 
75  for (auto const& theInput : inputHandle) {
76  if (!theInput.isValid()) {
77  edm::LogWarning("MTDTrackingRecHitProducer") << "MTDTrackingRecHitProducer: Invalid collection";
78  continue;
79  }
80  const edmNew::DetSetVector<FTLCluster>& input = *theInput;
81 
82  LogDebug("MTDTrackingRecHitProducer") << "inputCollection " << input.size();
83  for (const auto& DSVit : input) {
84  unsigned int detid = DSVit.detId();
85  DetId detIdObject(detid);
86  const auto genericDet = geom.idToDetUnit(detIdObject);
87  if (genericDet == nullptr) {
88  throw cms::Exception("MTDTrackingRecHitProducer")
89  << "GeographicalID: " << std::hex << detid << " is invalid!" << std::dec << std::endl;
90  }
91 
92  MTDTrackingDetSetVector::FastFiller recHitsOnDet(theoutputhits, detid);
93 
94  LogDebug("MTDTrackingRecHitProducer") << "MTD cluster DetId " << detid << " # cluster " << DSVit.size();
95 
96  for (const auto& clustIt : DSVit) {
97  LogDebug("MTDTrackingRecHitProducer") << "Cluster: size " << clustIt.size() << " " << clustIt.x() << ","
98  << clustIt.y() << " " << clustIt.energy() << " " << clustIt.time();
99  MTDClusterParameterEstimator::ReturnType tuple = cpe.getParameters(clustIt, *genericDet);
100  LocalPoint lp(std::get<0>(tuple));
101  LocalError le(std::get<1>(tuple));
102 
103  // Create a persistent edm::Ref to the cluster
105  // Make a RecHit and add it to the DetSet
106  MTDTrackingRecHit hit(lp, le, *genericDet, cluster);
107  LogDebug("MTDTrackingRecHitProducer")
108  << "MTD_TRH: " << hit.localPosition().x() << "," << hit.localPosition().y() << " : "
109  << hit.localPositionError().xx() << "," << hit.localPositionError().yy() << " : " << hit.time() << " : "
110  << hit.timeError();
111  // Now save it =================
112  recHitsOnDet.push_back(hit);
113  } // <-- End loop on Clusters
114  } // <-- End loop on DetUnits
115  LogDebug("MTDTrackingRecHitProducer") << "outputCollection " << theoutputhits.size();
116  }
117 
118  evt.put(std::move(outputhits));
119 }
120 
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
void push_back(data_type const &d)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::EDGetTokenT< FTLClusterCollection > ftlbClusters_
const edm::EDGetTokenT< FTLClusterCollection > ftleClusters_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
~MTDTrackingRecHitProducer() override=default
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static std::string const input
Definition: EdmProvDump.cc:50
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &es) const override
Definition: DetId.h:17
MTDTrackingRecHitProducer(const edm::ParameterSet &ps)
A 2D TrackerRecHit with time and time error information.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > mtdgeoToken_
const edm::ESGetToken< MTDClusterParameterEstimator, MTDCPERecord > cpeToken_
Log< level::Warning, false > LogWarning
std::tuple< LocalPoint, LocalError, TimeValue, TimeValueError > ReturnType
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)