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) do { std::cout << x << std::endl; } while (0)
26 #else
27 #define DEBUG(x)
28 #endif
29 
31 
32  public:
33  explicit MTDTrackingRecHitProducer(const edm::ParameterSet& ps);
34  ~MTDTrackingRecHitProducer() override = default;
35  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
36 
37  void produce(edm::Event& evt, const edm::EventSetup& es) override;
40 
41  private:
42  const edm::EDGetTokenT<FTLClusterCollection> ftlbClusters_; // collection of barrel digis
43  const edm::EDGetTokenT<FTLClusterCollection> ftleClusters_; // collection of endcap digis
44 
47 };
48 
50  ftlbClusters_( consumes<FTLClusterCollection>( ps.getParameter<edm::InputTag>("barrelClusters") ) ),
51  ftleClusters_( consumes<FTLClusterCollection>( ps.getParameter<edm::InputTag>("endcapClusters") ) )
52 {
53  produces< MTDTrackingDetSetVector >();
54 }
55 
56 // Configuration descriptions
57 void
60  desc.add<edm::InputTag>("barrelClusters", edm::InputTag("mtdClusters:FTLBarrel"));
61  desc.add<edm::InputTag>("endcapClusters", edm::InputTag("mtdClusters:FTLEndcap"));
62  descriptions.add("mtdTrackingRecHitProducer", desc);
63 }
64 
65 void
67 
69  es.get<MTDDigiGeometryRecord>().get(geom);
70  geom_ = geom.product();
71 
73  es.get<MTDCPERecord>().get("MTDCPEBase",cpe);
74  cpe_ = cpe.product();
75 
77  evt.getByToken( ftlbClusters_, inputBarrel);
78 
80  evt.getByToken( ftleClusters_, inputEndcap);
81 
82  auto outputhits = std::make_unique<MTDTrackingDetSetVector>();
83  auto& theoutputhits = *outputhits;
84 
85  run(inputBarrel,theoutputhits);
86  run(inputEndcap,theoutputhits);
87 
88  evt.put(std::move(outputhits));
89 }
90 
91 //---------------------------------------------------------------------------
94 //---------------------------------------------------------------------------
97 {
98  const edmNew::DetSetVector<FTLCluster>& input = *inputHandle;
100 
101  DEBUG("inputCollection " << input.size());
102  for ( ; DSViter != input.end() ; DSViter++) {
103  unsigned int detid = DSViter->detId();
104  DetId detIdObject( detid );
105  const auto& genericDet = geom_->idToDetUnit(detIdObject);
106  if( genericDet == nullptr ) {
107  throw cms::Exception("MTDTrackingRecHitProducer") << "GeographicalID: " << std::hex
108  << detid
109  << " is invalid!" << std::dec
110  << std::endl;
111  }
112 
113  MTDTrackingDetSetVector::FastFiller recHitsOnDet(output,detid);
114 
115  edmNew::DetSet<FTLCluster>::const_iterator clustIt = DSViter->begin(), clustEnd = DSViter->end();
116 
117  for ( ; clustIt != clustEnd; clustIt++) {
118  DEBUG("Cluster: size " << clustIt->size() << " " << clustIt->x() << "," << clustIt->y() << " " << clustIt->energy() << " " << clustIt->time());
119  MTDClusterParameterEstimator::ReturnType tuple = cpe_->getParameters( *clustIt, *genericDet );
120  LocalPoint lp( std::get<0>(tuple) );
121  LocalError le( std::get<1>(tuple) );
122 
123 
124  // Create a persistent edm::Ref to the cluster
125  edm::Ref< edmNew::DetSetVector<FTLCluster>, FTLCluster > cluster = edmNew::makeRefTo( inputHandle, clustIt);
126  // Make a RecHit and add it to the DetSet
127  MTDTrackingRecHit hit( lp, le, *genericDet, cluster);
128  DEBUG("MTD_TRH: " << hit.localPosition().x() << "," << hit.localPosition().y()
129  << " : " << hit.localPositionError().xx() << "," << hit.localPositionError().yy()
130  << " : " << hit.time() << " : " << hit.timeError());
131  // Now save it =================
132  recHitsOnDet.push_back(hit);
133  } // <-- End loop on Clusters
134  } // <-- End loop on DetUnits
135  DEBUG("outputCollection " << output.size());
136 }
137 
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:24
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
LocalError localPositionError() const final
float time() const
const edm::EDGetTokenT< FTLClusterCollection > ftlbClusters_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
T y() const
Definition: PV3DBase.h:63
const edm::EDGetTokenT< FTLClusterCollection > ftleClusters_
#define DEBUG(x)
~MTDTrackingRecHitProducer() override=default
data_type const * const_iterator
Definition: DetSetNew.h:30
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:26
const MTDClusterParameterEstimator * cpe_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::tuple< LocalPoint, LocalError, TimeValue, TimeValueError > ReturnType
Definition: DetId.h:18
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:71
LocalPoint localPosition() const final
size_type size() const
Definition: DetSetNew.h:87
float timeError() const
T x() const
Definition: PV3DBase.h:62
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:174