CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl Class Reference

#include <MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl.h>

Inheritance diagram for MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl:
reco::MtdRecoClusterToSimLayerClusterAssociatorBaseImpl

Public Member Functions

reco::RecoToSimCollectionMtd associateRecoToSim (const edm::Handle< FTLClusterCollection > &btlRecoClusH, const edm::Handle< FTLClusterCollection > &etlRecoClusH, const edm::Handle< MtdSimLayerClusterCollection > &simClusH) const override
 Associate a MtdRecoCluster to MtdSimLayerClusters. More...
 
reco::SimToRecoCollectionMtd associateSimToReco (const edm::Handle< FTLClusterCollection > &btlRecoClusH, const edm::Handle< FTLClusterCollection > &etlRecoClusH, const edm::Handle< MtdSimLayerClusterCollection > &simClusH) const override
 Associate a MtdSimLayerClusters to MtdRecoClusters. More...
 
 MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl (edm::EDProductGetter const &, double, double, mtd::MTDGeomUtil &)
 
- Public Member Functions inherited from reco::MtdRecoClusterToSimLayerClusterAssociatorBaseImpl
 MtdRecoClusterToSimLayerClusterAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~MtdRecoClusterToSimLayerClusterAssociatorBaseImpl ()
 Destructor. More...
 

Private Attributes

const double energyCut_
 
mtd::MTDGeomUtil geomTools_
 
edm::EDProductGetter const * productGetter_
 
const double timeCut_
 

Detailed Description

Definition at line 11 of file MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl.h.

Constructor & Destructor Documentation

◆ MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl()

MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl::MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl ( edm::EDProductGetter const &  productGetter,
double  energyCut,
double  timeCut,
mtd::MTDGeomUtil geomTools 
)
explicit

Member Function Documentation

◆ associateRecoToSim()

reco::RecoToSimCollectionMtd MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl::associateRecoToSim ( const edm::Handle< FTLClusterCollection > &  btlRecoClusH,
const edm::Handle< FTLClusterCollection > &  etlRecoClusH,
const edm::Handle< MtdSimLayerClusterCollection > &  simClusH 
) const
overridevirtual

Associate a MtdRecoCluster to MtdSimLayerClusters.

Reimplemented from reco::MtdRecoClusterToSimLayerClusterAssociatorBaseImpl.

Definition at line 19 of file MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl.cc.

References funct::abs(), MtdRecoClusterToSimLayerClusterAssociationMap::emplace_back(), energyCut_, dqmdumpme::first, geomTools_, mps_fire::i, l1ctLayer2EG_cff::id, LogDebug, edmNew::makeRefTo(), MtdRecoClusterToSimLayerClusterAssociationMap::post_insert(), edm::Handle< T >::product(), DetId::rawId(), mtd::MTDGeomUtil::sensorModuleId(), DBoxMetadataHelper::set_intersection(), timeCut_, and HcalDetIdTransform::transform().

22  {
23  RecoToSimCollectionMtd outputCollection;
24 
25  // -- get the collections
26  std::array<edm::Handle<FTLClusterCollection>, 2> inputRecoClusH{{btlRecoClusH, etlRecoClusH}};
27 
28  const auto& simClusters = *simClusH.product();
29 
30  // -- create temporary map DetId, SimClusterRef
31  std::map<uint32_t, std::vector<MtdSimLayerClusterRef>> simClusIdsMap;
32  for (auto simClusIt = simClusters.begin(); simClusIt != simClusters.end(); simClusIt++) {
33  const auto& simClus = *simClusIt;
34 
36  edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIt - simClusters.begin());
37 
38  std::vector<std::pair<uint32_t, std::pair<uint8_t, uint8_t>>> detIdsAndRows = simClus.detIds_and_rows();
39  std::vector<uint32_t> detIds(detIdsAndRows.size());
40  std::transform(detIdsAndRows.begin(),
41  detIdsAndRows.end(),
42  detIds.begin(),
43  [](const std::pair<uint32_t, std::pair<uint8_t, uint8_t>>& pair) { return pair.first; });
44 
45  // -- get the sensor module id from the first hit id of the MtdSimCluster
46  DetId id(detIds[0]);
47  uint32_t modId = geomTools_.sensorModuleId(id);
48  simClusIdsMap[modId].push_back(simClusterRef);
49 
50  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl") << "Sim cluster detId = " << modId << std::endl;
51  }
52 
53  for (auto const& recoClusH : inputRecoClusH) {
54  // -- loop over detSetVec
55  for (const auto& detSet : *recoClusH) {
56  // -- loop over reco clusters
57  for (const auto& recoClus : detSet) {
58  MTDDetId clusId = recoClus.id();
59 
60  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl") << "Reco cluster : " << clusId;
61 
62  std::vector<uint64_t> recoClusHitIds;
63 
64  // -- loop over hits in the reco cluster and find their unique ids
65  for (int ihit = 0; ihit < recoClus.size(); ++ihit) {
66  int hit_row = recoClus.minHitRow() + recoClus.hitOffset()[ihit * 2];
67  int hit_col = recoClus.minHitCol() + recoClus.hitOffset()[ihit * 2 + 1];
68 
69  // -- Get an unique id from sensor module detId , row, column
70  uint64_t uniqueId = static_cast<uint64_t>(clusId.rawId()) << 32;
71  uniqueId |= hit_row << 16;
72  uniqueId |= hit_col;
73  recoClusHitIds.push_back(uniqueId);
74 
75  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
76  << " ======= cluster raw Id : " << clusId.rawId() << " row : " << hit_row << " column: " << hit_col
77  << " recHit uniqueId : " << uniqueId;
78  }
79 
80  // -- loop over sim clusters and check if this reco clus shares some hits
81  std::vector<MtdSimLayerClusterRef> simClusterRefs;
82 
83  for (const auto& simClusterRef : simClusIdsMap[clusId.rawId()]) {
84  const auto& simClus = *simClusterRef;
85 
86  // get the hit ids of hits in the SimCluster
87  std::vector<std::pair<uint32_t, std::pair<uint8_t, uint8_t>>> detIdsAndRows = simClus.detIds_and_rows();
88  std::vector<uint64_t> simClusHitIds(detIdsAndRows.size());
89  for (unsigned int i = 0; i < detIdsAndRows.size(); i++) {
90  DetId id(detIdsAndRows[i].first);
91  uint32_t modId = geomTools_.sensorModuleId(id);
92  // build the unique id from sensor module detId , row, column
93  uint64_t uniqueId = static_cast<uint64_t>(modId) << 32;
94  uniqueId |= detIdsAndRows[i].second.first << 16;
95  uniqueId |= detIdsAndRows[i].second.second;
96  simClusHitIds.push_back(uniqueId);
97  }
98 
99  // -- Get shared hits
100  std::vector<uint64_t> sharedHitIds;
101  std::set_intersection(recoClusHitIds.begin(),
102  recoClusHitIds.end(),
103  simClusHitIds.begin(),
104  simClusHitIds.end(),
105  std::back_inserter(sharedHitIds));
106 
107  if (sharedHitIds.empty())
108  continue;
109 
110  float dE = recoClus.energy() * 0.001 / simClus.simLCEnergy(); // reco cluster energy is in MeV!
111  float dtSig = std::abs((recoClus.time() - simClus.simLCTime()) / recoClus.timeError());
112 
113  // -- If the sim and reco clusters have common hits, fill the std:vector of sim clusters refs
114  if (!sharedHitIds.empty() && dE < energyCut_ &&
115  dtSig < timeCut_) { // at least one hit in common + requirement on energy and time compatibility
116 
117  simClusterRefs.push_back(simClusterRef);
118 
119  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
120  << "RecoToSim --> Found " << sharedHitIds.size() << " shared hits";
121  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
122  << "E_recoClus = " << recoClus.energy() << " E_simClus = " << simClus.simLCEnergy()
123  << " E_recoClus/E_simClus = " << dE;
124  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
125  << "(t_recoClus-t_simClus)/sigma_t = " << dtSig;
126  }
127  } // -- end loop over sim clus refs
128 
129  // -- Now fill the output collection
130  edm::Ref<edmNew::DetSetVector<FTLCluster>, FTLCluster> recoClusterRef = edmNew::makeRefTo(recoClusH, &recoClus);
131  outputCollection.emplace_back(recoClusterRef, simClusterRefs);
132 
133  } // -- end loop over reco clus
134  } // -- end loop over detsetclus
135  }
136 
137  outputCollection.post_insert();
138  return outputCollection;
139 }
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)
T const * product() const
Definition: Handle.h:70
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
uint32_t sensorModuleId(const DetId &id) const
Definition: MTDGeomUtil.cc:165
Definition: DetId.h:17
unsigned long long uint64_t
Definition: Time.h:13
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
#define LogDebug(id)
void emplace_back(const FTLClusterRef &recoClus, std::vector< MtdSimLayerClusterRef > &simClusVect)
std::vector< std::string > set_intersection(std::vector< std::string > const &v1, std::vector< std::string > const &v2)
unsigned transform(const HcalDetId &id, unsigned transformCode)

◆ associateSimToReco()

reco::SimToRecoCollectionMtd MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl::associateSimToReco ( const edm::Handle< FTLClusterCollection > &  btlRecoClusH,
const edm::Handle< FTLClusterCollection > &  etlRecoClusH,
const edm::Handle< MtdSimLayerClusterCollection > &  simClusH 
) const
overridevirtual

Associate a MtdSimLayerClusters to MtdRecoClusters.

Reimplemented from reco::MtdRecoClusterToSimLayerClusterAssociatorBaseImpl.

Definition at line 141 of file MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl.cc.

References funct::abs(), MtdSimLayerClusterToRecoClusterAssociationMap::emplace_back(), energyCut_, dqmdumpme::first, geomTools_, mps_fire::i, l1ctLayer2EG_cff::id, LogDebug, edmNew::makeRefTo(), MtdSimLayerClusterToRecoClusterAssociationMap::post_insert(), edm::Handle< T >::product(), DetId::rawId(), mtd::MTDGeomUtil::sensorModuleId(), DBoxMetadataHelper::set_intersection(), and timeCut_.

144  {
145  SimToRecoCollectionMtd outputCollection;
146 
147  // -- get the collections
148  const auto& simClusters = *simClusH.product();
149 
150  std::array<edm::Handle<FTLClusterCollection>, 2> inputH{{btlRecoClusH, etlRecoClusH}};
151 
152  // -- loop over MtdSimLayerClusters
153  for (auto simClusIt = simClusters.begin(); simClusIt != simClusters.end(); simClusIt++) {
154  const auto& simClus = *simClusIt;
155 
156  // get the hit ids of hits in the SimCluster
157  std::vector<std::pair<uint32_t, std::pair<uint8_t, uint8_t>>> detIdsAndRows = simClus.detIds_and_rows();
158  std::vector<uint64_t> simClusHitIds(detIdsAndRows.size());
159  uint32_t modId = 0;
160  for (unsigned int i = 0; i < detIdsAndRows.size(); i++) {
161  DetId id(detIdsAndRows[i].first);
162  modId = geomTools_.sensorModuleId(id);
163  // build the unique id from sensor module detId , row, column
164  uint64_t uniqueId = static_cast<uint64_t>(modId) << 32;
165  uniqueId |= detIdsAndRows[i].second.first << 16;
166  uniqueId |= detIdsAndRows[i].second.second;
167  simClusHitIds.push_back(uniqueId);
168  }
169 
170  DetId simClusId(modId);
171 
172  std::vector<FTLClusterRef> recoClusterRefs;
173 
174  // -- loop over reco clusters
175  for (auto const& recoClusH : inputH) {
176  // -- loop over detSetVec
177  for (const auto& detSet : *recoClusH) {
178  if (detSet.id() != simClusId.rawId())
179  continue;
180 
181  // -- loop over reco clusters
182  for (const auto& recoClus : detSet) {
183  MTDDetId clusId = recoClus.id();
184 
185  std::vector<uint64_t> recoClusHitIds;
186 
187  // -- loop over hits in the reco cluster and find their ids
188  for (int ihit = 0; ihit < recoClus.size(); ++ihit) {
189  int hit_row = recoClus.minHitRow() + recoClus.hitOffset()[ihit * 2];
190  int hit_col = recoClus.minHitCol() + recoClus.hitOffset()[ihit * 2 + 1];
191 
192  // -- Get an unique id from sensor module detId , row, column
193  uint64_t uniqueId = static_cast<uint64_t>(clusId.rawId()) << 32;
194  uniqueId |= hit_row << 16;
195  uniqueId |= hit_col;
196  recoClusHitIds.push_back(uniqueId);
197 
198  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
199  << " ======= cluster raw Id : " << clusId.rawId() << " row : " << hit_row << " column: " << hit_col
200  << " recHit uniqueId : " << uniqueId;
201  }
202 
203  // -- Get shared hits
204  std::vector<uint64_t> sharedHitIds;
205  std::set_intersection(simClusHitIds.begin(),
206  simClusHitIds.end(),
207  recoClusHitIds.begin(),
208  recoClusHitIds.end(),
209  std::back_inserter(sharedHitIds));
210  if (sharedHitIds.empty())
211  continue;
212 
213  float dE = recoClus.energy() * 0.001 / simClus.simLCEnergy(); // reco cluster energy is in MeV
214  float dtSig = std::abs((recoClus.time() - simClus.simLCTime()) / recoClus.timeError());
215 
216  // -- If the sim and reco clusters have common hits, fill the std:vector of reco clusters refs
217  if (!sharedHitIds.empty() && dE < energyCut_ && dtSig < timeCut_) {
218  // Create a persistent edm::Ref to the cluster
220  edmNew::makeRefTo(recoClusH, &recoClus);
221  recoClusterRefs.push_back(recoClusterRef);
222 
223  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
224  << "SimToReco --> Found " << sharedHitIds.size() << " shared hits";
225  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
226  << "E_recoClus = " << recoClus.energy() << " E_simClus = " << simClus.simLCEnergy()
227  << " E_recoClus/E_simClus = " << dE;
228  LogDebug("MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl")
229  << "(t_recoClus-t_simClus)/sigma_t = " << dtSig;
230  }
231 
232  } // end loop ove reco clus
233  } // end loop over detsets
234  }
235 
236  // -- Now fill the output collection
238  edm::Ref<MtdSimLayerClusterCollection>(simClusH, simClusIt - simClusters.begin());
239  outputCollection.emplace_back(simClusterRef, recoClusterRefs);
240 
241  } // -- end loop over sim clusters
242 
243  outputCollection.post_insert();
244  return outputCollection;
245 }
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 emplace_back(const MtdSimLayerClusterRef &simClus, std::vector< FTLClusterRef > &recoClusVect)
T const * product() const
Definition: Handle.h:70
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
uint32_t sensorModuleId(const DetId &id) const
Definition: MTDGeomUtil.cc:165
Definition: DetId.h:17
unsigned long long uint64_t
Definition: Time.h:13
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
#define LogDebug(id)
std::vector< std::string > set_intersection(std::vector< std::string > const &v1, std::vector< std::string > const &v2)

Member Data Documentation

◆ energyCut_

const double MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl::energyCut_
private

◆ geomTools_

mtd::MTDGeomUtil MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl::geomTools_
private

◆ productGetter_

edm::EDProductGetter const* MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl::productGetter_
private

◆ timeCut_

const double MtdRecoClusterToSimLayerClusterAssociatorByHitsImpl::timeCut_
private