CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFCluster.h
Go to the documentation of this file.
1 #ifndef DataFormats_ParticleFlowReco_PFCluster_h
2 #define DataFormats_ParticleFlowReco_PFCluster_h
3 
4 #include <iostream>
5 #include <vector>
6 #include <algorithm>
7 
8 #include <Rtypes.h>
9 
16 #include "Math/GenVector/PositionVector3D.h"
17 
18 class PFClusterAlgo;
19 
20 namespace reco {
21 
42  class PFCluster : public CaloCluster {
43  public:
44  typedef std::vector<std::pair<CaloClusterPtr::key_type, edm::Ptr<PFCluster>>> EEtoPSAssociation;
45  // Next typedef uses double in ROOT 6 rather than Double32_t due to a bug in ROOT 5,
46  // which otherwise would make ROOT5 files unreadable in ROOT6. This does not increase
47  // the size on disk, because due to the bug, double was actually stored on disk in ROOT 5.
48  typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double>> REPPoint;
49 
51 
53  PFCluster(PFLayer::Layer layer, double energy, double x, double y, double z);
54 
56  void reset();
57 
59  void resetHitsAndFractions();
60 
63 
65  const std::vector<reco::PFRecHitFraction>& recHitFractions() const { return rechits_; }
66 
69 
71  PFLayer::Layer layer() const;
72 
74  double energy() const { return energy_; }
75 
77  float time() const { return time_; }
79  float timeError() const { return timeError_; }
80 
82  double depth() const { return depth_; }
83 
84  void setTime(float time, float timeError = 0) {
85  time_ = time;
87  }
89  void setDepth(double depth) { depth_ = depth; }
90 
92  const REPPoint& positionREP() const { return posrep_; }
93 
95  void calculatePositionREP() { posrep_.SetCoordinates(position_.Rho(), position_.Eta(), position_.Phi()); }
96 
98  static double getDepthCorrection(double energy, bool isBelowPS = false, bool isHadron = false);
99 
100  PFCluster& operator=(const PFCluster&);
101 
105 
107  double charge() const { return 0; }
108 
110  double pt() const { return (energy() * sin(position_.theta())); }
111 
113  double theta() const { return position_.theta(); }
114 
116  math::XYZPoint const& vertex() const { return dummyVtx_; }
117  double vx() const { return vertex().x(); }
118  double vy() const { return vertex().y(); }
119  double vz() const { return vertex().z(); }
120 
121 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
122  template <typename pruner>
123  void pruneUsing(pruner prune) {
124  // remove_if+erase algo applied to both vectors...
125  auto iter = std::find_if_not(rechits_.begin(), rechits_.end(), prune);
126  if (iter == rechits_.end())
127  return;
128  auto first = iter - rechits_.begin();
129  for (auto i = first; ++i < int(rechits_.size());) {
130  if (prune(rechits_[i])) {
133  ++first;
134  }
135  }
136  rechits_.erase(rechits_.begin() + first, rechits_.end());
138  }
139 #endif
140 
141  private:
143  std::vector<reco::PFRecHitFraction> rechits_;
144 
147 
149  float time_;
150  float timeError_{0};
151  double depth_{0};
152 
155 
157  static const constexpr double depthCorA_ = 0.89;
158  static const constexpr double depthCorB_ = 7.3;
159  static const constexpr double depthCorAp_ = 0.89;
160  static const constexpr double depthCorBp_ = 4.0;
161 
162  static const math::XYZPoint dummyVtx_;
163  };
164 
165  std::ostream& operator<<(std::ostream& out, const PFCluster& cluster);
166 
167 } // namespace reco
168 
169 #endif // DataFormats_ParticleFlowReco_PFCluster_h
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:56
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:49
math::XYZPoint position_
cluster centroid position
Definition: CaloCluster.h:228
double energy_
cluster energy
Definition: CaloCluster.h:223
float time() const
Definition: PFCluster.h:77
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
float timeError_
Definition: PFCluster.h:150
PFLayer::Layer layer_
transient layer
Definition: PFCluster.h:154
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void pruneUsing(pruner prune)
Definition: PFCluster.h:123
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:110
PFCluster & operator=(const PFCluster &)
Definition: PFCluster.cc:63
double z() const
z coordinate of cluster centroid
Definition: CaloCluster.h:178
std::vector< reco::PFRecHitFraction > rechits_
vector of rechit fractions (transient)
Definition: PFCluster.h:143
double vx() const
Definition: PFCluster.h:117
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:95
static const constexpr double depthCorAp_
Definition: PFCluster.h:159
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:92
def move
Definition: eostools.py:511
static const math::XYZPoint dummyVtx_
Definition: PFCluster.h:162
layer definition for PFRecHit and PFCluster
Definition: PFLayer.h:20
double x() const
x coordinate of cluster centroid
Definition: CaloCluster.h:172
double vy() const
Definition: PFCluster.h:118
float timeError() const
Definition: PFCluster.h:79
void reset()
resets clusters parameters
Definition: PFCluster.cc:17
double energy() const
cluster energy
Definition: PFCluster.h:74
std::vector< std::pair< DetId, float > > hitsAndFractions_
Definition: CaloCluster.h:234
void resetHitsAndFractions()
reset only hits and fractions
Definition: PFCluster.cc:28
Layer
layer definition
Definition: PFLayer.h:29
static const constexpr double depthCorB_
Definition: PFCluster.h:158
static double getDepthCorrection(double energy, bool isBelowPS=false, bool isHadron=false)
Definition: PFCluster.cc:39
void setTimeError(float timeError)
Definition: PFCluster.h:88
static const constexpr double depthCorA_
depth corrections
Definition: PFCluster.h:157
REPPoint posrep_
cluster position: rho, eta, phi (transient)
Definition: PFCluster.h:146
std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
Definition: PFCluster.h:44
double vz() const
Definition: PFCluster.h:119
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double theta() const
angle
Definition: PFCluster.h:113
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
static const constexpr double depthCorBp_
Definition: PFCluster.h:160
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
math::XYZPoint const & vertex() const
dummy vertex access
Definition: PFCluster.h:116
float time_
Michalis: add timing and depth information.
Definition: PFCluster.h:149
void setDepth(double depth)
Definition: PFCluster.h:89
double charge() const
dummy charge
Definition: PFCluster.h:107
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
double y() const
y coordinate of cluster centroid
Definition: CaloCluster.h:175
double depth() const
cluster depth
Definition: PFCluster.h:82