CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  const edmNew::DetSetVector<FTLCluster>& input = *theInput;
77 
78  LogDebug("MTDTrackingRecHitProducer") << "inputCollection " << input.size();
79  for (const auto& DSVit : input) {
80  unsigned int detid = DSVit.detId();
81  DetId detIdObject(detid);
82  const auto genericDet = geom.idToDetUnit(detIdObject);
83  if (genericDet == nullptr) {
84  throw cms::Exception("MTDTrackingRecHitProducer")
85  << "GeographicalID: " << std::hex << detid << " is invalid!" << std::dec << std::endl;
86  }
87 
88  MTDTrackingDetSetVector::FastFiller recHitsOnDet(theoutputhits, detid);
89 
90  for (const auto& clustIt : DSVit) {
91  LogDebug("MTDTrackingRcHitProducer") << "Cluster: size " << clustIt.size() << " " << clustIt.x() << ","
92  << clustIt.y() << " " << clustIt.energy() << " " << clustIt.time();
93  MTDClusterParameterEstimator::ReturnType tuple = cpe.getParameters(clustIt, *genericDet);
94  LocalPoint lp(std::get<0>(tuple));
95  LocalError le(std::get<1>(tuple));
96 
97  // Create a persistent edm::Ref to the cluster
99  // Make a RecHit and add it to the DetSet
100  MTDTrackingRecHit hit(lp, le, *genericDet, cluster);
101  LogDebug("MTDTrackingRcHitProducer")
102  << "MTD_TRH: " << hit.localPosition().x() << "," << hit.localPosition().y() << " : "
103  << hit.localPositionError().xx() << "," << hit.localPositionError().yy() << " : " << hit.time() << " : "
104  << hit.timeError();
105  // Now save it =================
106  recHitsOnDet.push_back(hit);
107  } // <-- End loop on Clusters
108  } // <-- End loop on DetUnits
109  LogDebug("MTDTrackingRcHitProducer") << "outputCollection " << theoutputhits.size();
110  }
111 
112  evt.put(std::move(outputhits));
113 }
114 
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)
float xx() const
Definition: LocalError.h:22
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
float time() const
const edm::EDGetTokenT< FTLClusterCollection > ftlbClusters_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T y() const
Definition: PV3DBase.h:60
const edm::EDGetTokenT< FTLClusterCollection > ftleClusters_
~MTDTrackingRecHitProducer() override=default
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static std::string const input
Definition: EdmProvDump.cc:47
bool getData(T &iHolder) const
Definition: EventSetup.h:128
float yy() const
Definition: LocalError.h:24
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &es) const override
LocalError localPositionError() 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)
size_type size() const
LocalPoint localPosition() const override
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > mtdgeoToken_
const edm::ESGetToken< MTDClusterParameterEstimator, MTDCPERecord > cpeToken_
float timeError() const
T x() const
Definition: PV3DBase.h:59
std::tuple< LocalPoint, LocalError, TimeValue, TimeValueError > ReturnType
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
#define LogDebug(id)