CMS 3D CMS Logo

MTDTrackingRecHitProducer.cc
Go to the documentation of this file.
6 
9 
16 
19 
21 
24 
26 public:
27  explicit MTDTrackingRecHitProducer(const edm::ParameterSet& ps);
28  ~MTDTrackingRecHitProducer() override = default;
29  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
30 
31  void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& es) const override;
32 
33 private:
34  const edm::EDGetTokenT<FTLClusterCollection> ftlbClusters_; // collection of barrel digis
35  const edm::EDGetTokenT<FTLClusterCollection> ftleClusters_; // collection of endcap digis
36 
39 };
40 
42  : ftlbClusters_(consumes<FTLClusterCollection>(ps.getParameter<edm::InputTag>("barrelClusters"))),
43  ftleClusters_(consumes<FTLClusterCollection>(ps.getParameter<edm::InputTag>("endcapClusters"))),
45  cpeToken_(esConsumes<MTDClusterParameterEstimator, MTDCPERecord>(edm::ESInputTag("", "MTDCPEBase"))) {
46  produces<MTDTrackingDetSetVector>();
47 }
48 
49 // Configuration descriptions
52  desc.add<edm::InputTag>("barrelClusters", edm::InputTag("mtdClusters:FTLBarrel"));
53  desc.add<edm::InputTag>("endcapClusters", edm::InputTag("mtdClusters:FTLEndcap"));
54  descriptions.add("mtdTrackingRecHitProducer", desc);
55 }
56 
58  auto const& geom = es.getData(mtdgeoToken_);
59 
60  auto const& cpe = es.getData(cpeToken_);
61 
63  evt.getByToken(ftlbClusters_, inputBarrel);
64 
66  evt.getByToken(ftleClusters_, inputEndcap);
67 
68  std::array<edm::Handle<FTLClusterCollection>, 2> inputHandle{{inputBarrel, inputEndcap}};
69 
70  auto outputhits = std::make_unique<MTDTrackingDetSetVector>();
71  auto& theoutputhits = *outputhits;
72 
73  //---------------------------------------------------------------------------
76  //---------------------------------------------------------------------------
77 
78  for (auto const& theInput : inputHandle) {
79  if (!theInput.isValid()) {
80  edm::LogWarning("MTDTrackingRecHitProducer") << "MTDTrackingRecHitProducer: Invalid collection";
81  continue;
82  }
83  const edmNew::DetSetVector<FTLCluster>& input = *theInput;
84 
85  LogDebug("MTDTrackingRecHitProducer") << "inputCollection " << input.size();
86  for (const auto& DSVit : input) {
87  unsigned int detid = DSVit.detId();
88  DetId detIdObject(detid);
89  const auto genericDet = geom.idToDetUnit(detIdObject);
90  if (genericDet == nullptr) {
91  throw cms::Exception("MTDTrackingRecHitProducer")
92  << "GeographicalID: " << std::hex << detid << " is invalid!" << std::dec << std::endl;
93  }
94 
95  MTDTrackingDetSetVector::FastFiller recHitsOnDet(theoutputhits, detid);
96 
97  LogDebug("MTDTrackingRecHitProducer") << "MTD cluster DetId " << detid << " # cluster " << DSVit.size();
98 #ifdef EDM_ML_DEBUG
99  const auto& Hit = MTDDetId(detid);
100  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 1)) {
101  const auto& btlHit = BTLDetId(detid);
102  LogDebug("MTDTrackingRecHitProducer") << btlHit;
103  } else if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 2)) {
104  const auto& etlHit = ETLDetId(detid);
105  LogDebug("MTDTrackingRecHitProducer") << etlHit;
106  }
107 #endif
108 
109  for (const auto& clustIt : DSVit) {
110  LogDebug("MTDTrackingRecHitProducer") << "Cluster: size " << clustIt.size() << " " << clustIt.x() << ","
111  << clustIt.y() << " " << clustIt.energy() << " " << clustIt.time();
112  MTDClusterParameterEstimator::ReturnType tuple = cpe.getParameters(clustIt, *genericDet);
113  LocalPoint lp(std::get<0>(tuple));
114  LocalError le(std::get<1>(tuple));
115 
116  // Create a persistent edm::Ref to the cluster
118  // Make a RecHit and add it to the DetSet
119  MTDTrackingRecHit hit(lp, le, *genericDet, cluster);
120  LogDebug("MTDTrackingRecHitProducer")
121  << "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  LogDebug("MTDTrackingRecHitProducer") << "outputCollection " << theoutputhits.size();
129  }
130 
131  evt.put(std::move(outputhits));
132 }
133 
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:528
~MTDTrackingRecHitProducer() override=default
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
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)
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:16
HLT enums.
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:19
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)