CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Types | Private Attributes
RealisticHitToClusterAssociator Class Reference

#include <RealisticHitToClusterAssociator.h>

Classes

struct  RealisticHit
 

Public Member Functions

void computeAssociation (float exclusiveFraction, bool useMCFractionsForExclEnergy, unsigned int fhOffset, unsigned int bhOffset)
 
void filterHitsByDistance (float maxDistance)
 
void findAndMergeInvisibleClusters (float invisibleFraction, float exclusiveFraction)
 
void findCentersOfGravity ()
 
void init (std::size_t numberOfHits, std::size_t numberOfSimClusters, std::size_t numberOfLayers)
 
void insertHitEnergy (float energy, unsigned int hitIndex)
 
void insertHitPosition (float x, float y, float z, unsigned int hitIndex)
 
void insertLayerId (unsigned int layerId, unsigned int hitIndex)
 
void insertSimClusterIdAndFraction (unsigned int scIdx, float fraction, unsigned int hitIndex, float associatedEnergy)
 
const std::vector< RealisticCluster > & realisticClusters () const
 
float XYdistanceFromMaxHit (unsigned int hitId, unsigned int clusterId)
 
float XYdistanceFromPointOnSameLayer (unsigned int hitId, const Hit3DPosition &point)
 

Private Types

using Hit3DPosition = std::array< float, 3 >
 

Private Attributes

std::vector< RealisticHitrealisticHits_
 
std::vector< RealisticClusterrealisticSimClusters_
 

Detailed Description

Definition at line 29 of file RealisticHitToClusterAssociator.h.

Member Typedef Documentation

using RealisticHitToClusterAssociator::Hit3DPosition = std::array<float, 3>
private

Definition at line 30 of file RealisticHitToClusterAssociator.h.

Member Function Documentation

void RealisticHitToClusterAssociator::computeAssociation ( float  exclusiveFraction,
bool  useMCFractionsForExclEnergy,
unsigned int  fhOffset,
unsigned int  bhOffset 
)
inline

Definition at line 82 of file RealisticHitToClusterAssociator.h.

References MillePedeFileConverter_cfg::e, JetChargeProducer_cfi::exp, f, or, realisticHits_, realisticSimClusters_, x, and XYdistanceFromMaxHit().

Referenced by RealisticSimClusterMapper::buildClusters().

85  {
86  //if more than exclusiveFraction of a hit's energy belongs to a cluster, that rechit is not counted as shared
87  unsigned int numberOfHits = realisticHits_.size();
88  std::vector<float> partialEnergies;
89  for (unsigned int hitId = 0; hitId < numberOfHits; ++hitId) {
90  partialEnergies.clear();
91  std::vector<unsigned int> removeAssociation;
92  auto& realisticHit = realisticHits_[hitId];
93  unsigned int numberOfClusters = realisticHit.hitToCluster_.size();
94  if (numberOfClusters == 1) {
95  unsigned int simClusterId = realisticHit.hitToCluster_[0].simClusterId_;
96  float assignedFraction = 1.f;
97  realisticHit.hitToCluster_[0].realisticEnergyFraction_ = assignedFraction;
98  float assignedEnergy = realisticHit.totalEnergy_;
99  realisticSimClusters_[simClusterId].increaseEnergy(assignedEnergy);
100  realisticSimClusters_[simClusterId].addHitAndFraction(hitId, assignedFraction);
101  realisticSimClusters_[simClusterId].increaseExclusiveEnergy(assignedEnergy);
102  } else {
103  partialEnergies.resize(numberOfClusters, 0.f);
104  unsigned int layer = realisticHit.layerId_;
105  float sumE = 0.f;
106  float energyDecayLength = getDecayLength(layer, fhOffset, bhOffset);
107  for (unsigned int clId = 0; clId < numberOfClusters; ++clId) {
108  auto simClusterId = realisticHit.hitToCluster_[clId].simClusterId_;
109  realisticHit.hitToCluster_[clId].distanceFromMaxHit_ = XYdistanceFromMaxHit(hitId, simClusterId);
110  // partial energy is computed based on the distance from the maximum energy hit and its energy
111  // partial energy is only needed to compute a fraction and it's not the energy assigned to the cluster
112  auto maxEnergyOnLayer = realisticSimClusters_[simClusterId].getMaxEnergy(layer);
113  if (maxEnergyOnLayer > 0.f) {
114  partialEnergies[clId] =
115  maxEnergyOnLayer * std::exp(-realisticHit.hitToCluster_[clId].distanceFromMaxHit_ / energyDecayLength);
116  }
117  sumE += partialEnergies[clId];
118  }
119  if (sumE > 0.f) {
120  float invSumE = 1.f / sumE;
121  for (unsigned int clId = 0; clId < numberOfClusters; ++clId) {
122  unsigned int simClusterIndex = realisticHit.hitToCluster_[clId].simClusterId_;
123  float assignedFraction = partialEnergies[clId] * invSumE;
124  if (assignedFraction > 1e-3) {
125  realisticHit.hitToCluster_[clId].realisticEnergyFraction_ = assignedFraction;
126  float assignedEnergy = assignedFraction * realisticHit.totalEnergy_;
127  realisticSimClusters_[simClusterIndex].increaseEnergy(assignedEnergy);
128  realisticSimClusters_[simClusterIndex].addHitAndFraction(hitId, assignedFraction);
129  // if the hits energy belongs for more than exclusiveFraction to a cluster, also the cluster's
130  // exclusive energy is increased. The exclusive energy will be needed to evaluate if
131  // a realistic cluster will be invisible, i.e. absorbed by other clusters
133  realisticHit.hitToCluster_[clId].mcEnergyFraction_ > exclusiveFraction) or
134  (!useMCFractionsForExclEnergy and assignedFraction > exclusiveFraction)) {
135  realisticSimClusters_[simClusterIndex].increaseExclusiveEnergy(assignedEnergy);
136  }
137  } else {
138  removeAssociation.push_back(simClusterIndex);
139  }
140  }
141  }
142  }
143 
144  while (!removeAssociation.empty()) {
145  auto clusterToRemove = removeAssociation.back();
146  removeAssociation.pop_back();
147 
148  realisticHit.hitToCluster_.erase(std::remove_if(realisticHit.hitToCluster_.begin(),
149  realisticHit.hitToCluster_.end(),
150  [clusterToRemove](const RealisticHit::HitToCluster& x) {
151  return x.simClusterId_ == clusterToRemove;
152  }),
153  realisticHit.hitToCluster_.end());
154  }
155  }
156  }
std::vector< RealisticCluster > realisticSimClusters_
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
double f[11][100]
float XYdistanceFromMaxHit(unsigned int hitId, unsigned int clusterId)
void RealisticHitToClusterAssociator::filterHitsByDistance ( float  maxDistance)
inline

Definition at line 266 of file RealisticHitToClusterAssociator.h.

References f, mps_fire::i, realisticHits_, realisticSimClusters_, and XYdistanceFromPointOnSameLayer().

Referenced by RealisticSimClusterMapper::buildClusters().

266  {
267  for (auto& cluster : realisticSimClusters_) {
268  if (cluster.isVisible()) {
269  auto& hAndF = cluster.hitsIdsAndFractions();
270  for (unsigned int i = 0; i < hAndF.size(); ++i) {
271  auto hitId = hAndF[i].first;
272  const auto& hit = realisticHits_[hitId];
273  auto layerId = hit.layerId_;
274  if (XYdistanceFromPointOnSameLayer(hitId, cluster.getCenterOfGravity(layerId)) > maxDistance) {
275  cluster.increaseEnergy(-hit.totalEnergy_ * hAndF[i].second);
276  cluster.modifyFractionByIndex(0.f, i);
277  }
278  }
279  }
280  }
281  }
std::vector< RealisticCluster > realisticSimClusters_
double f[11][100]
float XYdistanceFromPointOnSameLayer(unsigned int hitId, const Hit3DPosition &point)
void RealisticHitToClusterAssociator::findAndMergeInvisibleClusters ( float  invisibleFraction,
float  exclusiveFraction 
)
inline

Definition at line 158 of file RealisticHitToClusterAssociator.h.

References pfMETCorrectionType0_cfi::correction, MillePedeFileConverter_cfg::e, f, HLT_2018_cff::fraction, mps_fire::i, realisticHits_, realisticSimClusters_, and egamma::sharedEnergy().

Referenced by RealisticSimClusterMapper::buildClusters().

158  {
159  unsigned int numberOfRealSimClusters = realisticSimClusters_.size();
160 
161  for (unsigned int clId = 0; clId < numberOfRealSimClusters; ++clId) {
162  if (realisticSimClusters_[clId].getExclusiveEnergyFraction() < invisibleFraction) {
163  realisticSimClusters_[clId].setVisible(false);
164  auto& hAndF = realisticSimClusters_[clId].hitsIdsAndFractions();
165  std::unordered_map<unsigned int, float> energyInNeighbors;
166  float totalSharedEnergy = 0.f;
167 
168  for (auto& elt : hAndF) {
169  unsigned int hitId = elt.first;
170  float fraction = elt.second;
171  auto& realisticHit = realisticHits_[hitId];
172 
173  if (realisticHit.hitToCluster_.size() > 1 && fraction < 1.f) {
174  float correction = 1.f - fraction;
175  unsigned int numberOfClusters = realisticHit.hitToCluster_.size();
176  int clusterToRemove = -1;
177  for (unsigned int i = 0; i < numberOfClusters; ++i) {
178  auto simClusterIndex = realisticHit.hitToCluster_[i].simClusterId_;
179  if (simClusterIndex == clId) {
180  clusterToRemove = i;
181  } else if (realisticSimClusters_[simClusterIndex].isVisible()) {
182  float oldFraction = realisticHit.hitToCluster_[i].realisticEnergyFraction_;
183  float newFraction = oldFraction / correction;
184  float oldEnergy = oldFraction * realisticHit.totalEnergy_;
185  float newEnergy = newFraction * realisticHit.totalEnergy_;
186  float sharedEnergy = newEnergy - oldEnergy;
187  energyInNeighbors[simClusterIndex] += sharedEnergy;
188  totalSharedEnergy += sharedEnergy;
189  realisticSimClusters_[simClusterIndex].increaseEnergy(sharedEnergy);
190  realisticSimClusters_[simClusterIndex].modifyFractionForHitId(newFraction, hitId);
191  realisticHit.hitToCluster_[i].realisticEnergyFraction_ = newFraction;
192  if (newFraction > exclusiveFraction) {
193  realisticSimClusters_[simClusterIndex].increaseExclusiveEnergy(sharedEnergy);
194  if (oldFraction <= exclusiveFraction) {
195  realisticSimClusters_[simClusterIndex].increaseExclusiveEnergy(oldEnergy);
196  }
197  }
198  }
199  }
200  realisticSimClusters_[realisticHit.hitToCluster_[clusterToRemove].simClusterId_].modifyFractionForHitId(
201  0.f, hitId);
202  realisticHit.hitToCluster_.erase(realisticHit.hitToCluster_.begin() + clusterToRemove);
203  }
204  }
205 
206  for (auto& elt : hAndF) {
207  unsigned int hitId = elt.first;
208  auto& realisticHit = realisticHits_[hitId];
209  if (realisticHit.hitToCluster_.size() == 1 and realisticHit.hitToCluster_[0].simClusterId_ == clId and
210  totalSharedEnergy > 0.f) {
211  for (auto& pair : energyInNeighbors) {
212  // hits that belonged completely to the absorbed cluster are redistributed
213  // based on the fraction of energy shared in the shared hits
214  float sharedFraction = pair.second / totalSharedEnergy;
215  if (sharedFraction > 1e-6) {
216  float assignedEnergy = realisticHit.totalEnergy_ * sharedFraction;
217  realisticSimClusters_[pair.first].increaseEnergy(assignedEnergy);
218  realisticSimClusters_[pair.first].addHitAndFraction(hitId, sharedFraction);
219  realisticHit.hitToCluster_.emplace_back(
220  RealisticHit::HitToCluster{pair.first, 0.f, -1.f, sharedFraction});
221  if (sharedFraction > exclusiveFraction)
222  realisticSimClusters_[pair.first].increaseExclusiveEnergy(assignedEnergy);
223  }
224  }
225  }
226  }
227  }
228  }
229  }
std::vector< RealisticCluster > realisticSimClusters_
double f[11][100]
float sharedEnergy(reco::CaloCluster const &clu1, reco::CaloCluster const &clu2, EcalRecHitCollection const &barrelRecHits, EcalRecHitCollection const &endcapRecHits)
void RealisticHitToClusterAssociator::findCentersOfGravity ( )
inline

Definition at line 231 of file RealisticHitToClusterAssociator.h.

References f, HLT_2018_cff::fraction, realisticHits_, and realisticSimClusters_.

Referenced by RealisticSimClusterMapper::buildClusters().

231  {
232  for (auto& cluster : realisticSimClusters_) {
233  if (cluster.isVisible()) {
234  unsigned int layersNum = cluster.getLayersNum();
235  std::vector<float> totalEnergyPerLayer(layersNum, 0.f);
236  std::vector<float> xEnergyPerLayer(layersNum, 0.f);
237  std::vector<float> yEnergyPerLayer(layersNum, 0.f);
238  std::vector<float> zPositionPerLayer(layersNum, 0.f);
239  const auto& hAndF = cluster.hitsIdsAndFractions();
240  for (auto& elt : hAndF) {
241  auto hitId = elt.first;
242  auto fraction = elt.second;
243  const auto& hit = realisticHits_[hitId];
244  const auto& hitPos = hit.hitPosition_;
245  auto layerId = hit.layerId_;
246  auto hitEinCluster = hit.totalEnergy_ * fraction;
247  totalEnergyPerLayer[layerId] += hitEinCluster;
248  xEnergyPerLayer[layerId] += hitPos[0] * hitEinCluster;
249  yEnergyPerLayer[layerId] += hitPos[1] * hitEinCluster;
250  zPositionPerLayer[layerId] = hitPos[2];
251  }
252  Hit3DPosition centerOfGravity;
253  for (unsigned int layerId = 0; layerId < layersNum; layerId++) {
254  auto energyOnLayer = totalEnergyPerLayer[layerId];
255  if (energyOnLayer > 0.f) {
256  centerOfGravity = {{xEnergyPerLayer[layerId] / energyOnLayer,
257  yEnergyPerLayer[layerId] / energyOnLayer,
258  zPositionPerLayer[layerId]}};
259  cluster.setCenterOfGravity(layerId, centerOfGravity);
260  }
261  }
262  }
263  }
264  }
std::vector< RealisticCluster > realisticSimClusters_
double f[11][100]
void RealisticHitToClusterAssociator::init ( std::size_t  numberOfHits,
std::size_t  numberOfSimClusters,
std::size_t  numberOfLayers 
)
inline

Definition at line 47 of file RealisticHitToClusterAssociator.h.

References realisticHits_, realisticSimClusters_, and SimDataFormats::CaloAnalysis::sc.

Referenced by RealisticSimClusterMapper::buildClusters().

47  {
48  realisticHits_.resize(numberOfHits);
49  realisticSimClusters_.resize(numberOfSimClusters);
50  for (auto& sc : realisticSimClusters_)
51  sc.setLayersNum(numberOfLayers);
52  }
std::vector< RealisticCluster > realisticSimClusters_
const map< TString, int > numberOfLayers(TString Year="2018")
Definition: DMRtrends.cc:254
void RealisticHitToClusterAssociator::insertHitEnergy ( float  energy,
unsigned int  hitIndex 
)
inline
void RealisticHitToClusterAssociator::insertHitPosition ( float  x,
float  y,
float  z,
unsigned int  hitIndex 
)
inline

Definition at line 54 of file RealisticHitToClusterAssociator.h.

References realisticHits_, x, and y.

Referenced by RealisticSimClusterMapper::buildClusters().

54  {
55  realisticHits_[hitIndex].hitPosition_ = {{x, y, z}};
56  }
void RealisticHitToClusterAssociator::insertLayerId ( unsigned int  layerId,
unsigned int  hitIndex 
)
inline

Definition at line 58 of file RealisticHitToClusterAssociator.h.

References realisticHits_.

Referenced by RealisticSimClusterMapper::buildClusters().

58 { realisticHits_[hitIndex].layerId_ = layerId; }
void RealisticHitToClusterAssociator::insertSimClusterIdAndFraction ( unsigned int  scIdx,
float  fraction,
unsigned int  hitIndex,
float  associatedEnergy 
)
inline

Definition at line 62 of file RealisticHitToClusterAssociator.h.

References HLT_2018_cff::fraction, RealisticHitToClusterAssociator::RealisticHit::hitPosition_, RealisticHitToClusterAssociator::RealisticHit::layerId_, realisticHits_, and realisticSimClusters_.

Referenced by RealisticSimClusterMapper::buildClusters().

62  {
63  realisticHits_[hitIndex].hitToCluster_.emplace_back(RealisticHit::HitToCluster{scIdx, fraction, 0.f, 0.f});
64  realisticSimClusters_[scIdx].setMaxEnergyHit(
65  realisticHits_[hitIndex].layerId_, associatedEnergy, realisticHits_[hitIndex].hitPosition_);
66  }
std::vector< RealisticCluster > realisticSimClusters_
const std::vector<RealisticCluster>& RealisticHitToClusterAssociator::realisticClusters ( ) const
inline

Definition at line 283 of file RealisticHitToClusterAssociator.h.

References realisticSimClusters_.

Referenced by RealisticSimClusterMapper::buildClusters().

283 { return realisticSimClusters_; }
std::vector< RealisticCluster > realisticSimClusters_
float RealisticHitToClusterAssociator::XYdistanceFromMaxHit ( unsigned int  hitId,
unsigned int  clusterId 
)
inline

Definition at line 68 of file RealisticHitToClusterAssociator.h.

References RealisticHitToClusterAssociator::RealisticHit::hitPosition_, cmsLHEtoEOSManager::l, funct::pow(), realisticHits_, realisticSimClusters_, and mathSSE::sqrt().

Referenced by computeAssociation().

68  {
69  auto l = realisticHits_[hitId].layerId_;
70  const auto& maxHitPosition = realisticSimClusters_[clusterId].getMaxEnergyPosition(l);
71  float distanceSquared = std::pow((realisticHits_[hitId].hitPosition_[0] - maxHitPosition[0]), 2) +
72  std::pow((realisticHits_[hitId].hitPosition_[1] - maxHitPosition[1]), 2);
73  return std::sqrt(distanceSquared);
74  }
std::vector< RealisticCluster > realisticSimClusters_
T sqrt(T t)
Definition: SSEVec.h:19
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
float RealisticHitToClusterAssociator::XYdistanceFromPointOnSameLayer ( unsigned int  hitId,
const Hit3DPosition point 
)
inline

Definition at line 76 of file RealisticHitToClusterAssociator.h.

References RealisticHitToClusterAssociator::RealisticHit::hitPosition_, funct::pow(), realisticHits_, and mathSSE::sqrt().

Referenced by filterHitsByDistance().

76  {
77  float distanceSquared = std::pow((realisticHits_[hitId].hitPosition_[0] - point[0]), 2) +
78  std::pow((realisticHits_[hitId].hitPosition_[1] - point[1]), 2);
79  return std::sqrt(distanceSquared);
80  }
T sqrt(T t)
Definition: SSEVec.h:19
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5

Member Data Documentation

std::vector<RealisticHit> RealisticHitToClusterAssociator::realisticHits_
private
std::vector<RealisticCluster> RealisticHitToClusterAssociator::realisticSimClusters_
private