CMS 3D CMS Logo

PFlow2DClusterizerWithTime.cc
Go to the documentation of this file.
5 
6 #include "Math/GenVector/VectorUtil.h"
7 
8 #include "vdt/vdtMath.h"
9 
10 #include "TMath.h"
11 
12 #include <iterator>
13 
14 #ifdef PFLOW_DEBUG
15 #define LOGVERB(x) edm::LogVerbatim(x)
16 #define LOGWARN(x) edm::LogWarning(x)
17 #define LOGERR(x) edm::LogError(x)
18 #define LOGDRESSED(x) edm::LogInfo(x)
19 #else
20 #define LOGVERB(x) LogTrace(x)
21 #define LOGWARN(x) edm::LogWarning(x)
22 #define LOGERR(x) edm::LogError(x)
23 #define LOGDRESSED(x) LogDebug(x)
24 #endif
25 
29  _maxIterations(conf.getParameter<unsigned>("maxIterations")),
30  _stoppingTolerance(conf.getParameter<double>("stoppingTolerance")),
31  _showerSigma2(std::pow(conf.getParameter<double>("showerSigma"),2.0)),
32  _timeSigma_eb(std::pow(conf.getParameter<double>("timeSigmaEB"),2.0)),
33  _timeSigma_ee(std::pow(conf.getParameter<double>("timeSigmaEE"),2.0)),
34  _excludeOtherSeeds(conf.getParameter<bool>("excludeOtherSeeds")),
35  _minFracTot(conf.getParameter<double>("minFracTot")),
36  _maxNSigmaTime(std::pow(conf.getParameter<double>("maxNSigmaTime"),2.0)),
37  _minChi2Prob(conf.getParameter<double>("minChi2Prob")),
38  _clusterTimeResFromSeed(conf.getParameter<bool>("clusterTimeResFromSeed")),
39 
40  _layerMap({ {"PS2",(int)PFLayer::PS2},
41  {"PS1",(int)PFLayer::PS1},
42  {"ECAL_ENDCAP",(int)PFLayer::ECAL_ENDCAP},
43  {"ECAL_BARREL",(int)PFLayer::ECAL_BARREL},
44  {"NONE",(int)PFLayer::NONE},
45  {"HCAL_BARREL1",(int)PFLayer::HCAL_BARREL1},
46  {"HCAL_BARREL2_RING0",(int)PFLayer::HCAL_BARREL2},
47  {"HCAL_BARREL2_RING1",100*(int)PFLayer::HCAL_BARREL2},
48  {"HCAL_ENDCAP",(int)PFLayer::HCAL_ENDCAP},
49  {"HF_EM",(int)PFLayer::HF_EM},
50  {"HF_HAD",(int)PFLayer::HF_HAD} }) {
51 
52 
53  const std::vector<edm::ParameterSet>& thresholds =
54  conf.getParameterSetVector("recHitEnergyNorms");
55  for( const auto& pset : thresholds ) {
56  const std::string& det = pset.getParameter<std::string>("detector");
57  const double& rhE_norm = pset.getParameter<double>("recHitEnergyNorm");
58  auto entry = _layerMap.find(det);
59  if( entry == _layerMap.end() ) {
60  throw cms::Exception("InvalidDetectorLayer")
61  << "Detector layer : " << det << " is not in the list of recognized"
62  << " detector layers!";
63  }
64  _recHitEnergyNorms.emplace(_layerMap.find(det)->second,rhE_norm);
65  }
66 
67  if( conf.exists("allCellsPositionCalc") ) {
68  const edm::ParameterSet& acConf =
69  conf.getParameterSet("allCellsPositionCalc");
70  const std::string& algoac =
71  acConf.getParameter<std::string>("algoName");
72  _allCellsPosCalc = std::unique_ptr<PFCPositionCalculatorBase>{PFCPositionCalculatorFactory::get()->create(algoac, acConf)};
73  }
74  // if necessary a third pos calc for convergence testing
75  if( conf.exists("positionCalcForConvergence") ) {
76  const edm::ParameterSet& convConf =
77  conf.getParameterSet("positionCalcForConvergence");
78  const std::string& algoconv =
79  convConf.getParameter<std::string>("algoName");
80  _convergencePosCalc = std::unique_ptr<PFCPositionCalculatorBase>{PFCPositionCalculatorFactory::get()->create(algoconv, convConf)};
81  }
82  if( conf.exists("timeResolutionCalcBarrel") ) {
83  const edm::ParameterSet& timeResConf =
84  conf.getParameterSet("timeResolutionCalcBarrel");
85  _timeResolutionCalcBarrel = std::make_unique<CaloRecHitResolutionProvider>(
86  timeResConf);
87  }
88  if( conf.exists("timeResolutionCalcEndcap") ) {
89  const edm::ParameterSet& timeResConf =
90  conf.getParameterSet("timeResolutionCalcEndcap");
91  _timeResolutionCalcEndcap = std::make_unique<CaloRecHitResolutionProvider>(
92  timeResConf);
93  }
94 }
95 
98  const std::vector<bool>& seedable,
100  reco::PFClusterCollection clustersInTopo;
101  for( const auto& topocluster : input ) {
102  clustersInTopo.clear();
103  seedPFClustersFromTopo(topocluster,seedable,clustersInTopo);
104  const unsigned tolScal =
105  std::pow(std::max(1.0,clustersInTopo.size()-1.0),2.0);
106  growPFClusters(topocluster,seedable,tolScal,0,tolScal,clustersInTopo);
107  // step added by Josh Bendavid, removes low-fraction clusters
108  // did not impact position resolution with fraction cut of 1e-7
109  // decreases the size of each pf cluster considerably
110  prunePFClusters(clustersInTopo);
111  // recalculate the positions of the pruned clusters
112  if( _convergencePosCalc ) {
113  // if defined, use the special position calculation for convergence tests
114  _convergencePosCalc->calculateAndSetPositions(clustersInTopo);
115  } else {
116  if( clustersInTopo.size() == 1 && _allCellsPosCalc ) {
117  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back());
118  } else {
119  _positionCalc->calculateAndSetPositions(clustersInTopo);
120  }
121  }
122  for( auto& clusterout : clustersInTopo ) {
123  output.insert(output.end(),std::move(clusterout));
124  }
125  }
126 }
127 
130  const std::vector<bool>& seedable,
131  reco::PFClusterCollection& initialPFClusters) const {
132  const auto& recHitFractions = topo.recHitFractions();
133  for( const auto& rhf : recHitFractions ) {
134  if( !seedable[rhf.recHitRef().key()] ) continue;
135  initialPFClusters.push_back(reco::PFCluster());
136  reco::PFCluster& current = initialPFClusters.back();
137  current.addRecHitFraction(rhf);
138  current.setSeed(rhf.recHitRef()->detId());
139  if( _convergencePosCalc ) {
140  _convergencePosCalc->calculateAndSetPosition(current);
141  } else {
142  _positionCalc->calculateAndSetPosition(current);
143  }
144  }
145 }
146 
149  const std::vector<bool>& seedable,
150  const unsigned toleranceScaling,
151  const unsigned iter,
152  double diff,
154  if( iter >= _maxIterations ) {
155  LOGDRESSED("PFlow2DClusterizerWithTime:growAndStabilizePFClusters")
156  <<"reached " << _maxIterations << " iterations, terminated position "
157  << "fit with diff = " << diff;
158  }
159  if( iter >= _maxIterations ||
160  diff <= _stoppingTolerance*toleranceScaling) return;
161  // reset the rechits in this cluster, keeping the previous position
162  std::vector<reco::PFCluster::REPPoint> clus_prev_pos;
163  // also calculate and keep the previous time resolution
164  std::vector<double> clus_prev_timeres2;
165 
166  for( auto& cluster : clusters) {
167  const reco::PFCluster::REPPoint& repp = cluster.positionREP();
168  clus_prev_pos.emplace_back(repp.rho(),repp.eta(),repp.phi());
169  if( _convergencePosCalc ) {
170  if( clusters.size() == 1 && _allCellsPosCalc ) {
171  _allCellsPosCalc->calculateAndSetPosition(cluster);
172  } else {
173  _positionCalc->calculateAndSetPosition(cluster);
174  }
175  }
176  double resCluster2;
178  clusterTimeResolutionFromSeed(cluster, resCluster2);
179  else
180  clusterTimeResolution(cluster, resCluster2);
181  clus_prev_timeres2.push_back(resCluster2);
182  cluster.resetHitsAndFractions();
183  }
184  // loop over topo cluster and grow current PFCluster hypothesis
185  std::vector<double> dist2, frac;
186 
187  // Store chi2 values and nhits to calculate chi2 prob. inline
188  std::vector<double> clus_chi2;
189  std::vector<size_t> clus_chi2_nhits;
190 
191  double fractot = 0;
192  for( const reco::PFRecHitFraction& rhf : topo.recHitFractions() ) {
193  const reco::PFRecHitRef& refhit = rhf.recHitRef();
194  int cell_layer = (int)refhit->layer();
195  if( cell_layer == PFLayer::HCAL_BARREL2 &&
196  std::abs(refhit->positionREP().eta()) > 0.34 ) {
197  cell_layer *= 100;
198  }
199  const double recHitEnergyNorm =
200  _recHitEnergyNorms.find(cell_layer)->second;
201  math::XYZPoint topocellpos_xyz(refhit->position());
202  dist2.clear(); frac.clear(); fractot = 0;
203 
204  // add rechits to clusters, calculating fraction based on distance
205  // need position in vector to get cluster time resolution
206  for (size_t iCluster = 0; iCluster < clusters.size(); ++iCluster) {
207  reco::PFCluster& cluster = clusters[iCluster];
208  const math::XYZPoint& clusterpos_xyz = cluster.position();
209  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
210  double d2 = deltav.Mag2()/_showerSigma2;
211 
212  double d2time = dist2Time(cluster, refhit, cell_layer,
213  clus_prev_timeres2[iCluster]);
214  d2 += d2time;
215 
216  if (_minChi2Prob > 0. && !passChi2Prob(iCluster, d2time, clus_chi2,
217  clus_chi2_nhits))
218  d2 = 999.;
219 
220  dist2.emplace_back( d2);
221 
222  if( d2 > 100 ) {
223  LOGDRESSED("PFlow2DClusterizerWithTime:growAndStabilizePFClusters")
224  << "Warning! :: pfcluster-topocell distance is too large! d= "
225  << d2;
226  }
227  // fraction assignment logic
228  double fraction;
229  if( refhit->detId() == cluster.seed() && _excludeOtherSeeds ) {
230  fraction = 1.0;
231  } else if ( seedable[refhit.key()] && _excludeOtherSeeds ) {
232  fraction = 0.0;
233  } else {
234  fraction = cluster.energy()/recHitEnergyNorm * vdt::fast_expf( -0.5*d2 );
235  }
236  fractot += fraction;
237  frac.emplace_back(fraction);
238  }
239  for( unsigned i = 0; i < clusters.size(); ++i ) {
240  if( fractot > _minFracTot ||
241  ( refhit->detId() == clusters[i].seed() && fractot > 0.0 ) ) {
242  frac[i]/=fractot;
243  } else {
244  continue;
245  }
246  // if the fraction has been set to 0, the cell
247  // is now added to the cluster - careful ! (PJ, 19/07/08)
248  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
249  // (PJ, 15/09/08 <- similar to what existed before the
250  // previous bug fix, but keeps the close seeds inside,
251  // even if their fraction was set to zero.)
252  // Also add a protection to keep the seed in the cluster
253  // when the latter gets far from the former. These cases
254  // (about 1% of the clusters) need to be studied, as
255  // they create fake photons, in general.
256  // (PJ, 16/09/08)
257  if( dist2[i] < 100.0 || frac[i] > 0.9999 ) {
258  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit,frac[i]));
259  }
260  }
261  }
262  // recalculate positions and calculate convergence parameter
263  double diff2 = 0.0;
264  for( unsigned i = 0; i < clusters.size(); ++i ) {
265  if( _convergencePosCalc ) {
266  _convergencePosCalc->calculateAndSetPosition(clusters[i]);
267  } else {
268  if( clusters.size() == 1 && _allCellsPosCalc ) {
269  _allCellsPosCalc->calculateAndSetPosition(clusters[i]);
270  } else {
271  _positionCalc->calculateAndSetPosition(clusters[i]);
272  }
273  }
274  const double delta2 =
275  reco::deltaR2(clusters[i].positionREP(),clus_prev_pos[i]);
276  if( delta2 > diff2 ) diff2 = delta2;
277  }
278  diff = std::sqrt(diff2);
279  dist2.clear(); frac.clear(); clus_prev_pos.clear();// avoid badness
280  clus_chi2.clear(); clus_chi2_nhits.clear(); clus_prev_timeres2.clear();
281  growPFClusters(topo,seedable,toleranceScaling,iter+1,diff,clusters);
282 }
283 
286  for( auto& cluster : clusters ) {
287  cluster.pruneUsing( [&](const reco::PFRecHitFraction& rhf)
288  {return rhf.fraction() > _minFractionToKeep;}
289  );
290  }
291 }
292 
294  cluster, double& clusterRes2) const
295 {
296  clusterRes2 = 10000.;
297  for (auto& rhf : cluster.recHitFractions())
298  {
299  const reco::PFRecHit& rh = *(rhf.recHitRef());
300  if (rh.detId() == cluster.seed())
301  {
302  cluster.setTime(rh.time());
303  bool isBarrel = (rh.layer() == PFLayer::HCAL_BARREL1 ||
304  rh.layer() == PFLayer::HCAL_BARREL2 ||
305  rh.layer() == PFLayer::ECAL_BARREL);
306  if (isBarrel)
307  clusterRes2 = _timeResolutionCalcBarrel->timeResolution2(rh.energy());
308  else
309  clusterRes2 = _timeResolutionCalcEndcap->timeResolution2(rh.energy());
310  cluster.setTimeError(std::sqrt(float(clusterRes2)));
311  }
312  }
313 }
314 
316  double& clusterRes2) const
317 {
318  double sumTimeSigma2 = 0.;
319  double sumSigma2 = 0.;
320 
321  for (auto& rhf : cluster.recHitFractions())
322  {
323  const reco::PFRecHit& rh = *(rhf.recHitRef());
324  const double rhf_f = rhf.fraction();
325 
326  if (rhf_f == 0.)
327  continue;
328 
329  bool isBarrel = (rh.layer() == PFLayer::HCAL_BARREL1 ||
330  rh.layer() == PFLayer::HCAL_BARREL2 ||
331  rh.layer() == PFLayer::ECAL_BARREL);
332  double res2 = 10000.;
333  if (isBarrel)
334  res2 = _timeResolutionCalcBarrel->timeResolution2(rh.energy());
335  else
336  res2 = _timeResolutionCalcEndcap->timeResolution2(rh.energy());
337 
338  sumTimeSigma2 += rhf_f * rh.time()/res2;
339  sumSigma2 += rhf_f/res2;
340  }
341  if (sumSigma2 > 0.) {
342  clusterRes2 = 1./sumSigma2;
343  cluster.setTime(sumTimeSigma2/sumSigma2);
344  cluster.setTimeError(std::sqrt(float(clusterRes2)));
345  } else {
346  clusterRes2 = 999999.;
347  }
348 }
349 
351  const reco::PFRecHitRef& refhit, int cell_layer, double prev_timeres2) const
352 {
353  const double deltaT = cluster.time()-refhit->time();
354  const double t2 = deltaT*deltaT;
355  double res2 = 100.;
356 
357  if (cell_layer == PFLayer::HCAL_BARREL1 ||
358  cell_layer == PFLayer::HCAL_BARREL2 ||
359  cell_layer == PFLayer::ECAL_BARREL) {
361  const double resCluster2 = prev_timeres2;
362  res2 = resCluster2 + _timeResolutionCalcBarrel->timeResolution2(
363  refhit->energy());
364  }
365  else {
366  return t2/_timeSigma_eb;
367  }
368  }
369  else if (cell_layer == PFLayer::HCAL_ENDCAP ||
370  cell_layer == PFLayer::HF_EM ||
371  cell_layer == PFLayer::HF_HAD ||
372  cell_layer == PFLayer::ECAL_ENDCAP) {
374  const double resCluster2 = prev_timeres2;
375  res2 = resCluster2 + _timeResolutionCalcEndcap->timeResolution2(
376  refhit->energy());
377  }
378  else {
379  return t2/_timeSigma_ee;
380  }
381  }
382 
383  double distTime2 = t2/res2;
384  if (distTime2 > _maxNSigmaTime)
385  return 999.; // reject hit
386 
387  return distTime2;
388 }
389 
390 bool PFlow2DClusterizerWithTime::passChi2Prob(size_t iCluster, double dist2,
391  std::vector<double>& clus_chi2, std::vector<size_t>& clus_chi2_nhits) const
392 {
393  if (iCluster >= clus_chi2.size()) { // first hit
394  clus_chi2.push_back(dist2);
395  clus_chi2_nhits.push_back(1);
396  return true;
397  }
398  else {
399  double chi2 = clus_chi2[iCluster];
400  size_t nhitsCluster = clus_chi2_nhits[iCluster];
401  chi2 += dist2;
402  if (TMath::Prob(chi2, nhitsCluster) >= _minChi2Prob) {
403  clus_chi2[iCluster] = chi2;
404  clus_chi2_nhits[iCluster] = nhitsCluster + 1;
405  return true;
406  }
407  }
408  return false;
409 }
T getParameter(std::string const &) const
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:131
void prunePFClusters(reco::PFClusterCollection &) const
PFlow2DClusterizerWithTime(const edm::ParameterSet &conf)
float time() const
Definition: PFCluster.h:85
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
float time() const
timing for cleaned hits
Definition: PFRecHit.h:102
unsigned detId() const
rechit detId
Definition: PFRecHit.h:93
bool passChi2Prob(size_t iCluster, double dist2, std::vector< double > &clus_chi2, std::vector< size_t > &clus_chi2_nhits) const
std::unique_ptr< PosCalc > _positionCalc
bool isBarrel(GeomDetEnumerators::SubDetector m)
void setTime(float time, float timeError=0)
Definition: PFCluster.h:92
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
#define LOGDRESSED(x)
key_type key() const
Accessor for product key.
Definition: Ref.h:263
double fraction() const
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
static std::string const input
Definition: EdmProvDump.cc:48
void setSeed(const DetId &id)
Definition: CaloCluster.h:123
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:96
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
void clusterTimeResolution(reco::PFCluster &cluster, double &res) const
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float energy() const
rechit energy
Definition: PFRecHit.h:99
double dist2Time(const reco::PFCluster &, const reco::PFRecHitRef &, int cell_layer, double prev_timeres2) const
double energy() const
cluster energy
Definition: PFCluster.h:82
void setTimeError(float timeError)
Definition: PFCluster.h:93
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:207
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
ParameterSet const & getParameterSet(std::string const &) const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::unordered_map< int, double > _recHitEnergyNorms
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:91
float fast_expf(float x)
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:54
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
void clusterTimeResolutionFromSeed(reco::PFCluster &cluster, double &res) const
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus) override
const std::unordered_map< std::string, int > _layerMap
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:511
T get(const Candidate &c)
Definition: component.h:55