test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Basic2DGenericPFlowPositionCalc.cc
Go to the documentation of this file.
2 
5 
6 #include <cmath>
8 #include<iterator>
9 #include <boost/function_output_iterator.hpp>
10 
11 #include "vdt/vdtMath.h"
12 
13 
14 namespace {
15  inline
16  bool isBarrel(int cell_layer){
17  return (cell_layer == PFLayer::HCAL_BARREL1 ||
18  cell_layer == PFLayer::HCAL_BARREL2 ||
19  cell_layer == PFLayer::ECAL_BARREL);
20  }
21 }
22 
26 }
27 
30  for( reco::PFCluster& cluster : clusters ) {
32  }
33 }
34 
37  if( !cluster.seed() ) {
38  throw cms::Exception("ClusterWithNoSeed")
39  << " Found a cluster with no seed: " << cluster;
40  }
41  double cl_energy = 0;
42  double cl_time = 0;
43  double cl_timeweight=0.0;
44  double max_e = 0.0;
45  PFLayer::Layer max_e_layer = PFLayer::NONE;
46  // find the seed and max layer and also calculate time
47  //Michalis : Even if we dont use timing in clustering here we should fill
48  //the time information for the cluster. This should use the timing resolution(1/E)
49  //so the weight should be fraction*E^2
50  //calculate a simplistic depth now. The log weighted will be done
51  //in different stage
52 
53 
54  auto const recHitCollection = &(*cluster.recHitFractions()[0].recHitRef()) - cluster.recHitFractions()[0].recHitRef().key();
55  auto nhits = cluster.recHitFractions().size();
56  struct LHit{ reco::PFRecHit const * hit; float energy; float fraction;};
57  declareDynArray(LHit,nhits,hits);
58  for(auto i=0U; i<nhits; ++i) {
59  auto const & hf = cluster.recHitFractions()[i];
60  auto k = hf.recHitRef().key();
61  auto p = recHitCollection+k;
62  hits[i]= {p,(*p).energy(), float(hf.fraction())};
63  }
64 
65 
66  bool resGiven = bool(_timeResolutionCalcBarrel) & bool(_timeResolutionCalcEndcap);
67  LHit mySeed={nullptr};
68  for( auto const & rhf : hits ) {
69  const reco::PFRecHit & refhit = *rhf.hit;
70  if( refhit.detId() == cluster.seed() ) mySeed = rhf;
71  const auto rh_fraction = rhf.fraction;
72  const auto rh_rawenergy = rhf.energy;
73  const auto rh_energy = rh_rawenergy * rh_fraction;
74 #ifdef PF_DEBUG
75  if unlikely( edm::isNotFinite(rh_energy) ) {
76  throw cms::Exception("PFClusterAlgo")
77  <<"rechit " << refhit.detId() << " has a NaN energy... "
78  << "The input of the particle flow clustering seems to be corrupted.";
79  }
80 #endif
81  cl_energy += rh_energy;
82  // If time resolution is given, calculated weighted average
83  if ( resGiven ) {
84  double res2 = 1.e-4;
85  int cell_layer = (int)refhit.layer();
86  res2 = isBarrel(cell_layer) ? 1./_timeResolutionCalcBarrel->timeResolution2(rh_rawenergy) :
87  1./_timeResolutionCalcEndcap->timeResolution2(rh_rawenergy);
88  cl_time += rh_fraction*refhit.time()*res2;
89  cl_timeweight += rh_fraction*res2;
90  }
91  else { // assume resolution = 1/E**2
92  const double rh_rawenergy2 = rh_rawenergy*rh_rawenergy;
93  cl_timeweight+=rh_rawenergy2*rh_fraction;
94  cl_time += rh_rawenergy2*rh_fraction*refhit.time();
95  }
96 
97  if( rh_energy > max_e ) {
98  max_e = rh_energy;
99  max_e_layer = refhit.layer();
100  }
101  }
102 
103  cluster.setEnergy(cl_energy);
104  cluster.setTime(cl_time/cl_timeweight);
105  cluster.setLayer(max_e_layer);
106 
107  // calculate the position
108  double depth = 0.0;
109  double position_norm = 0.0;
110  double x(0.0),y(0.0),z(0.0);
111  if( nullptr != mySeed.hit ) {
112  auto seedNeighbours = mySeed.hit->neighbours();
113  switch( _posCalcNCrystals ) {
114  case 5:
115  seedNeighbours = mySeed.hit->neighbours4();
116  break;
117  case 9:
118  seedNeighbours = mySeed.hit->neighbours8();
119  break;
120  default:
121  break;
122  }
123 
124  auto compute = [&](LHit const& rhf) {
125  const reco::PFRecHit & refhit = *rhf.hit;
126  const auto rh_energy = rhf.energy * rhf.fraction;
127  const auto norm = ( rhf.fraction < _minFractionInCalc ?
128  0.0f :
129  std::max(0.0f,vdt::fast_logf(rh_energy*_logWeightDenom)) );
130  const auto rhpos_xyz = refhit.position()*norm;
131  x += rhpos_xyz.x();
132  y += rhpos_xyz.y();
133  z += rhpos_xyz.z();
134  depth += refhit.depth()*norm;
135  position_norm += norm;
136  };
137 
138  if(_posCalcNCrystals != -1) // sorted to make neighbour search faster (maybe)
139  std::sort(hits.begin(),hits.end(),[](LHit const& a, LHit const& b) { return a.hit<b.hit;});
140 
141 
142  if(_posCalcNCrystals == -1)
143  for( auto const & rhf : hits ) compute(rhf);
144  else { // only seed and its neighbours
145  compute(mySeed);
146  // search seedNeighbours to find energy fraction in cluster (sic)
147  unInitDynArray(reco::PFRecHit const *,seedNeighbours.size(),nei);
148  for(auto k : seedNeighbours){
149  nei.push_back(recHitCollection+k);
150  }
151  std::sort(nei.begin(),nei.end());
152  struct LHitLess {
153  auto operator()(LHit const &a, reco::PFRecHit const * b) const {return a.hit<b;}
154  auto operator()(reco::PFRecHit const * b, LHit const &a) const {return b<a.hit;}
155  };
156  std::set_intersection(hits.begin(),hits.end(),nei.begin(),nei.end(),
157  boost::make_function_output_iterator(compute),
158  LHitLess()
159  );
160  }
161  } else {
162  throw cms::Exception("Basic2DGenerticPFlowPositionCalc")
163  << "Cluster seed hit is null, something is wrong with PFlow RecHit!";
164  }
165 
166  if( position_norm < _minAllowedNorm ) {
167  edm::LogError("WeirdClusterNormalization")
168  << "PFCluster too far from seeding cell: set position to (0,0,0).";
169  cluster.setPosition(math::XYZPoint(0,0,0));
170  cluster.calculatePositionREP();
171  } else {
172  const double norm_inverse = 1.0/position_norm;
173  x *= norm_inverse;
174  y *= norm_inverse;
175  z *= norm_inverse;
176  depth *= norm_inverse;
177  cluster.setPosition(math::XYZPoint(x,y,z));
178  cluster.setDepth(depth);
179  cluster.calculatePositionREP();
180  }
181 }
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:120
int i
Definition: DBlmapReader.cc:9
void calculateAndSetPositionActual(reco::PFCluster &) const
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
float time() const
timing for cleaned hits
Definition: PFRecHit.h:118
unsigned detId() const
rechit detId
Definition: PFRecHit.h:108
void setPosition(const math::XYZPoint &p)
Definition: CaloCluster.h:112
bool isBarrel(GeomDetEnumerators::SubDetector m)
void setTime(float time, float timeError=0)
Definition: PFCluster.h:92
void setEnergy(double energy)
Definition: CaloCluster.h:108
PositionType const & position() const
rechit cell centre x, y, z
Definition: PFRecHit.h:129
int depth() const
depth for segemntation
Definition: PFRecHit.h:121
void calculateAndSetPositions(reco::PFClusterCollection &)
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:111
#define unlikely(x)
bool isNotFinite(T x)
Definition: isFinite.h:10
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:100
#define unInitDynArray(T, n, x)
Definition: DynArray.h:58
double f[11][100]
float energy() const
rechit energy
Definition: PFRecHit.h:114
Layer
layer definition
Definition: PFLayer.h:31
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:202
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
double b
Definition: hdecay.h:120
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
double a
Definition: hdecay.h:121
void setDepth(double depth)
Definition: PFCluster.h:94
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
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
float fast_logf(float x)