CMS 3D CMS Logo

Basic2DGenericPFlowPositionCalc.cc
Go to the documentation of this file.
9 
10 #include "vdt/vdtMath.h"
11 
12 #include <boost/iterator/function_output_iterator.hpp>
13 
14 #include <cmath>
15 #include <iterator>
16 #include <memory>
17 #include <tuple>
18 
20 public:
23  _posCalcNCrystals(conf.getParameter<int>("posCalcNCrystals")),
24  _minAllowedNorm(conf.getParameter<double>("minAllowedNormalization")) {
25  std::vector<int> detectorEnum;
26  std::vector<int> depths;
27  std::vector<double> logWeightDenom;
28  std::vector<float> logWeightDenomInv;
29 
30  if (conf.exists("logWeightDenominatorByDetector")) {
31  const std::vector<edm::ParameterSet>& logWeightDenominatorByDetectorPSet =
32  conf.getParameterSetVector("logWeightDenominatorByDetector");
33 
34  for (const auto& pset : logWeightDenominatorByDetectorPSet) {
35  if (!pset.exists("detector")) {
36  throw cms::Exception("logWeightDenominatorByDetectorPSet") << "logWeightDenominator : detector not specified";
37  }
38 
39  const std::string& det = pset.getParameter<std::string>("detector");
40 
41  if (det == std::string("HCAL_BARREL1") || det == std::string("HCAL_ENDCAP")) {
42  std::vector<int> depthsT = pset.getParameter<std::vector<int> >("depths");
43  std::vector<double> logWeightDenomT = pset.getParameter<std::vector<double> >("logWeightDenominator");
44  if (logWeightDenomT.size() != depthsT.size()) {
45  throw cms::Exception("logWeightDenominator") << "logWeightDenominator mismatch with the numbers of depths";
46  }
47  for (unsigned int i = 0; i < depthsT.size(); ++i) {
48  if (det == std::string("HCAL_BARREL1"))
49  detectorEnum.push_back(1);
50  if (det == std::string("HCAL_ENDCAP"))
51  detectorEnum.push_back(2);
52  depths.push_back(depthsT[i]);
53  logWeightDenom.push_back(logWeightDenomT[i]);
54  }
55  }
56  }
57  } else {
58  detectorEnum.push_back(0);
59  depths.push_back(0);
60  logWeightDenom.push_back(conf.getParameter<double>("logWeightDenominator"));
61  }
62 
63  for (unsigned int i = 0; i < depths.size(); ++i) {
64  logWeightDenomInv.push_back(1. / logWeightDenom[i]);
65  }
66 
67  // _logWeightDenom = std::make_pair(depths,logWeightDenomInv);
68  _logWeightDenom = std::make_tuple(detectorEnum, depths, logWeightDenomInv);
69 
70  _timeResolutionCalcBarrel.reset(nullptr);
71  if (conf.exists("timeResolutionCalcBarrel")) {
72  const edm::ParameterSet& timeResConf = conf.getParameterSet("timeResolutionCalcBarrel");
73  _timeResolutionCalcBarrel = std::make_unique<CaloRecHitResolutionProvider>(timeResConf);
74  }
75  _timeResolutionCalcEndcap.reset(nullptr);
76  if (conf.exists("timeResolutionCalcEndcap")) {
77  const edm::ParameterSet& timeResConf = conf.getParameterSet("timeResolutionCalcEndcap");
78  _timeResolutionCalcEndcap = std::make_unique<CaloRecHitResolutionProvider>(timeResConf);
79  }
80 
81  switch (_posCalcNCrystals) {
82  case 5:
83  case 9:
84  case -1:
85  break;
86  default:
87  edm::LogError("Basic2DGenericPFlowPositionCalc") << "posCalcNCrystals not valid";
88  assert(0); // bug
89  }
90  }
91 
94 
97 
98 private:
99  const int _posCalcNCrystals;
100  std::tuple<std::vector<int>, std::vector<int>, std::vector<float> > _logWeightDenom;
101  const float _minAllowedNorm;
102 
103  std::unique_ptr<CaloRecHitResolutionProvider> _timeResolutionCalcBarrel;
104  std::unique_ptr<CaloRecHitResolutionProvider> _timeResolutionCalcEndcap;
105 
107 };
108 
110 
111 namespace {
112  inline bool isBarrel(int cell_layer) {
113  return (cell_layer == PFLayer::HCAL_BARREL1 || cell_layer == PFLayer::HCAL_BARREL2 ||
114  cell_layer == PFLayer::ECAL_BARREL);
115  }
116 } // namespace
117 
120 }
121 
123  for (reco::PFCluster& cluster : clusters) {
125  }
126 }
127 
129  if (!cluster.seed()) {
130  throw cms::Exception("ClusterWithNoSeed") << " Found a cluster with no seed: " << cluster;
131  }
132  double cl_energy = 0;
133  double cl_time = 0;
134  double cl_timeweight = 0.0;
135  double max_e = 0.0;
136  PFLayer::Layer max_e_layer = PFLayer::NONE;
137  // find the seed and max layer and also calculate time
138  //Michalis : Even if we dont use timing in clustering here we should fill
139  //the time information for the cluster. This should use the timing resolution(1/E)
140  //so the weight should be fraction*E^2
141  //calculate a simplistic depth now. The log weighted will be done
142  //in different stage
143 
144  auto const recHitCollection =
145  &(*cluster.recHitFractions()[0].recHitRef()) - cluster.recHitFractions()[0].recHitRef().key();
146  auto nhits = cluster.recHitFractions().size();
147  struct LHit {
148  reco::PFRecHit const* hit;
149  float energy;
150  float fraction;
151  };
152  declareDynArray(LHit, nhits, hits);
153  for (auto i = 0U; i < nhits; ++i) {
154  auto const& hf = cluster.recHitFractions()[i];
155  auto k = hf.recHitRef().key();
156  auto p = recHitCollection + k;
157  hits[i] = {p, (*p).energy(), float(hf.fraction())};
158  }
159 
161  LHit mySeed = {};
162  for (auto const& rhf : hits) {
163  const reco::PFRecHit& refhit = *rhf.hit;
164  if (refhit.detId() == cluster.seed())
165  mySeed = rhf;
166  const auto rh_fraction = rhf.fraction;
167  const auto rh_rawenergy = rhf.energy;
168  const auto rh_energy = rh_rawenergy * rh_fraction;
169 #ifdef PF_DEBUG
170  if UNLIKELY (edm::isNotFinite(rh_energy)) {
171  throw cms::Exception("PFClusterAlgo") << "rechit " << refhit.detId() << " has a NaN energy... "
172  << "The input of the particle flow clustering seems to be corrupted.";
173  }
174 #endif
175  cl_energy += rh_energy;
176  // If time resolution is given, calculated weighted average
177  if (resGiven) {
178  double res2 = 1.e-4;
179  int cell_layer = (int)refhit.layer();
180  res2 = isBarrel(cell_layer) ? 1. / _timeResolutionCalcBarrel->timeResolution2(rh_rawenergy)
181  : 1. / _timeResolutionCalcEndcap->timeResolution2(rh_rawenergy);
182  cl_time += rh_fraction * refhit.time() * res2;
183  cl_timeweight += rh_fraction * res2;
184  } else { // assume resolution = 1/E**2
185  const double rh_rawenergy2 = rh_rawenergy * rh_rawenergy;
186  cl_timeweight += rh_rawenergy2 * rh_fraction;
187  cl_time += rh_rawenergy2 * rh_fraction * refhit.time();
188  }
189 
190  if (rh_energy > max_e) {
191  max_e = rh_energy;
192  max_e_layer = refhit.layer();
193  }
194  }
195 
196  cluster.setEnergy(cl_energy);
197  cluster.setTime(cl_time / cl_timeweight);
198  if (resGiven) {
199  cluster.setTimeError(std::sqrt(1.0f / float(cl_timeweight)));
200  }
201  cluster.setLayer(max_e_layer);
202 
203  // calculate the position
204  bool single_depth = true;
205  int ref_depth = -1;
206  double depth = 0.0;
207  double position_norm = 0.0;
208  double x(0.0), y(0.0), z(0.0);
209  if (nullptr != mySeed.hit) {
210  auto seedNeighbours = mySeed.hit->neighbours();
211  switch (_posCalcNCrystals) {
212  case 5:
213  seedNeighbours = mySeed.hit->neighbours4();
214  break;
215  case 9:
216  seedNeighbours = mySeed.hit->neighbours8();
217  break;
218  default:
219  break;
220  }
221 
222  auto compute = [&](LHit const& rhf) {
223  const reco::PFRecHit& refhit = *rhf.hit;
224 
225  int cell_layer = (int)refhit.layer();
226  float threshold = 0;
227 
228  for (unsigned int j = 0; j < (std::get<2>(_logWeightDenom)).size(); ++j) {
229  // barrel is detecor type1
230  int detectorEnum = std::get<0>(_logWeightDenom)[j];
231  int depth = std::get<1>(_logWeightDenom)[j];
232 
233  if ((cell_layer == PFLayer::HCAL_BARREL1 && detectorEnum == 1 && refhit.depth() == depth) ||
234  (cell_layer == PFLayer::HCAL_ENDCAP && detectorEnum == 2 && refhit.depth() == depth) || detectorEnum == 0)
235  threshold = std::get<2>(_logWeightDenom)[j];
236  }
237 
238  if (ref_depth < 0)
239  ref_depth = refhit.depth(); // Initialize reference depth
240  else if (refhit.depth() != ref_depth) {
241  // Found a rechit with a different depth
242  single_depth = false;
243  }
244  const auto rh_energy = rhf.energy * rhf.fraction;
245  const auto norm =
246  (rhf.fraction < _minFractionInCalc ? 0.0f : std::max(0.0f, vdt::fast_logf(rh_energy * threshold)));
247  const auto rhpos_xyz = refhit.position() * norm;
248  x += rhpos_xyz.x();
249  y += rhpos_xyz.y();
250  z += rhpos_xyz.z();
251  depth += refhit.depth() * norm;
252  position_norm += norm;
253  };
254 
255  if (_posCalcNCrystals != -1) // sorted to make neighbour search faster (maybe)
256  std::sort(hits.begin(), hits.end(), [](LHit const& a, LHit const& b) { return a.hit < b.hit; });
257 
258  if (_posCalcNCrystals == -1)
259  for (auto const& rhf : hits)
260  compute(rhf);
261  else { // only seed and its neighbours
262  compute(mySeed);
263  // search seedNeighbours to find energy fraction in cluster (sic)
264  unInitDynArray(reco::PFRecHit const*, seedNeighbours.size(), nei);
265  for (auto k : seedNeighbours) {
266  nei.push_back(recHitCollection + k);
267  }
268  std::sort(nei.begin(), nei.end());
269  struct LHitLess {
270  auto operator()(LHit const& a, reco::PFRecHit const* b) const { return a.hit < b; }
271  auto operator()(reco::PFRecHit const* b, LHit const& a) const { return b < a.hit; }
272  };
274  hits.begin(), hits.end(), nei.begin(), nei.end(), boost::make_function_output_iterator(compute), LHitLess());
275  }
276  } else {
277  throw cms::Exception("Basic2DGenerticPFlowPositionCalc")
278  << "Cluster seed hit is null, something is wrong with PFlow RecHit!";
279  }
280 
281  if (position_norm < _minAllowedNorm) {
282  edm::LogError("WeirdClusterNormalization") << "PFCluster too far from seeding cell: set position to (0,0,0).";
283  cluster.setPosition(math::XYZPoint(0, 0, 0));
284  cluster.calculatePositionREP();
285  } else {
286  const double norm_inverse = 1.0 / position_norm;
287  x *= norm_inverse;
288  y *= norm_inverse;
289  z *= norm_inverse;
290  if (single_depth)
291  depth = ref_depth;
292  else
293  depth *= norm_inverse;
294  cluster.setPosition(math::XYZPoint(x, y, z));
295  cluster.setDepth(depth);
296  cluster.calculatePositionREP();
297  }
298 }
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:49
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void setPosition(const math::XYZPoint &p)
Definition: CaloCluster.h:140
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
Basic2DGenericPFlowPositionCalc(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
bool exists(std::string const &parameterName) const
checks if a parameter exists
ParameterSet const & getParameterSet(std::string const &) const
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:96
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
void setEnergy(double energy)
Definition: CaloCluster.h:136
float time() const
timing for cleaned hits
Definition: PFRecHit.h:102
Log< level::Error, false > LogError
PositionType const & position() const
rechit cell centre x, y, z
Definition: PFRecHit.h:117
assert(be >=bs)
void calculateAndSetPositionActual(reco::PFCluster &) const
unsigned detId() const
rechit detId
Definition: PFRecHit.h:93
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:95
T sqrt(T t)
Definition: SSEVec.h:19
#define unInitDynArray(T, n, x)
Definition: DynArray.h:88
double f[11][100]
Layer
layer definition
Definition: PFLayer.h:29
void setTimeError(float timeError)
Definition: PFCluster.h:88
void calculateAndSetPositions(reco::PFClusterCollection &) override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
def compute(min, max)
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
double b
Definition: hdecay.h:118
void calculateAndSetPosition(reco::PFCluster &) override
std::tuple< std::vector< int >, std::vector< int >, std::vector< float > > _logWeightDenom
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
Basic2DGenericPFlowPositionCalc & operator=(const Basic2DGenericPFlowPositionCalc &)=delete
VParameterSet const & getParameterSetVector(std::string const &name) const
double a
Definition: hdecay.h:119
void setDepth(double depth)
Definition: PFCluster.h:89
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
int depth() const
depth for segemntation
Definition: PFRecHit.h:105
#define UNLIKELY(x)
Definition: Likely.h:21
float fast_logf(float x)
std::vector< std::string > set_intersection(std::vector< std::string > const &v1, std::vector< std::string > const &v2)