CMS 3D CMS Logo

MTDTrackingRecHitProducer.cc
Go to the documentation of this file.
6 
9 
16 
20 
22 
23 //#define DEBUG_ENABLED
24 #ifdef DEBUG_ENABLED
25 #define DEBUG(x) \
26  do { \
27  std::cout << x << std::endl; \
28  } while (0)
29 #else
30 #define DEBUG(x)
31 #endif
32 
34 public:
35  explicit MTDTrackingRecHitProducer(const edm::ParameterSet& ps);
36  ~MTDTrackingRecHitProducer() override = default;
37  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
38 
39  void produce(edm::Event& evt, const edm::EventSetup& es) override;
41 
42 private:
43  const edm::EDGetTokenT<FTLClusterCollection> ftlbClusters_; // collection of barrel digis
44  const edm::EDGetTokenT<FTLClusterCollection> ftleClusters_; // collection of endcap digis
45 
48 };
49 
51  : ftlbClusters_(consumes<FTLClusterCollection>(ps.getParameter<edm::InputTag>("barrelClusters"))),
52  ftleClusters_(consumes<FTLClusterCollection>(ps.getParameter<edm::InputTag>("endcapClusters"))) {
53  produces<MTDTrackingDetSetVector>();
54 }
55 
56 // Configuration descriptions
59  desc.add<edm::InputTag>("barrelClusters", edm::InputTag("mtdClusters:FTLBarrel"));
60  desc.add<edm::InputTag>("endcapClusters", edm::InputTag("mtdClusters:FTLEndcap"));
61  descriptions.add("mtdTrackingRecHitProducer", desc);
62 }
63 
66  es.get<MTDDigiGeometryRecord>().get(geom);
67  geom_ = geom.product();
68 
70  es.get<MTDCPERecord>().get("MTDCPEBase", cpe);
71  cpe_ = cpe.product();
72 
74  evt.getByToken(ftlbClusters_, inputBarrel);
75 
77  evt.getByToken(ftleClusters_, inputEndcap);
78 
79  auto outputhits = std::make_unique<MTDTrackingDetSetVector>();
80  auto& theoutputhits = *outputhits;
81 
82  run(inputBarrel, theoutputhits);
83  run(inputEndcap, theoutputhits);
84 
85  evt.put(std::move(outputhits));
86 }
87 
88 //---------------------------------------------------------------------------
91 //---------------------------------------------------------------------------
93  const edmNew::DetSetVector<FTLCluster>& input = *inputHandle;
95 
96  DEBUG("inputCollection " << input.size());
97  for (; DSViter != input.end(); DSViter++) {
98  unsigned int detid = DSViter->detId();
99  DetId detIdObject(detid);
100  const auto& genericDet = geom_->idToDetUnit(detIdObject);
101  if (genericDet == nullptr) {
102  throw cms::Exception("MTDTrackingRecHitProducer")
103  << "GeographicalID: " << std::hex << detid << " is invalid!" << std::dec << std::endl;
104  }
105 
106  MTDTrackingDetSetVector::FastFiller recHitsOnDet(output, detid);
107 
108  edmNew::DetSet<FTLCluster>::const_iterator clustIt = DSViter->begin(), clustEnd = DSViter->end();
109 
110  for (; clustIt != clustEnd; clustIt++) {
111  DEBUG("Cluster: size " << clustIt->size() << " " << clustIt->x() << "," << clustIt->y() << " "
112  << clustIt->energy() << " " << clustIt->time());
113  MTDClusterParameterEstimator::ReturnType tuple = cpe_->getParameters(*clustIt, *genericDet);
114  LocalPoint lp(std::get<0>(tuple));
115  LocalError le(std::get<1>(tuple));
116 
117  // Create a persistent edm::Ref to the cluster
119  // Make a RecHit and add it to the DetSet
120  MTDTrackingRecHit hit(lp, le, *genericDet, cluster);
121  DEBUG("MTD_TRH: " << hit.localPosition().x() << "," << hit.localPosition().y() << " : "
122  << hit.localPositionError().xx() << "," << hit.localPositionError().yy() << " : " << hit.time()
123  << " : " << hit.timeError());
124  // Now save it =================
125  recHitsOnDet.push_back(hit);
126  } // <-- End loop on Clusters
127  } // <-- End loop on DetUnits
128  DEBUG("outputCollection " << output.size());
129 }
130 
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
const_iterator end(bool update=false) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
LocalError localPositionError() const final
float time() const
const edm::EDGetTokenT< FTLClusterCollection > ftlbClusters_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
T y() const
Definition: PV3DBase.h:60
const edm::EDGetTokenT< FTLClusterCollection > ftleClusters_
#define DEBUG(x)
~MTDTrackingRecHitProducer() override=default
data_type const * const_iterator
Definition: DetSetNew.h:31
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static std::string const input
Definition: EdmProvDump.cc:48
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float yy() const
Definition: LocalError.h:24
const MTDClusterParameterEstimator * cpe_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
MTDTrackingRecHitProducer(const edm::ParameterSet &ps)
virtual ReturnType getParameters(const FTLCluster &cl, const GeomDetUnit &det) const =0
A 2D TrackerRecHit with time and time error information.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void run(edm::Handle< edmNew::DetSetVector< FTLCluster > > inputHandle, MTDTrackingDetSetVector &output)
size_type size() const
HLT enums.
T get() const
Definition: EventSetup.h:73
LocalPoint localPosition() const final
size_type size() const
Definition: DetSetNew.h:68
float timeError() const
T x() const
Definition: PV3DBase.h:59
std::tuple< LocalPoint, LocalError, TimeValue, TimeValueError > ReturnType
T const * product() const
Definition: ESHandle.h:86
void produce(edm::Event &evt, const edm::EventSetup &es) override
def move(src, dest)
Definition: eostools.py:511
const_iterator begin(bool update=false) const
const MTDGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: MTDGeometry.cc:152