CMS 3D CMS Logo

PFCluster.h
Go to the documentation of this file.
1 #ifndef DataFormats_ParticleFlowReco_PFCluster_h
2 #define DataFormats_ParticleFlowReco_PFCluster_h
3 
6 
7 #include "Math/GenVector/PositionVector3D.h"
9 #include "Rtypes.h"
10 
14 
15 #include <iostream>
16 #include <vector>
17 #include <algorithm>
18 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
19 #include <atomic>
20 #endif
21 
22 
23 class PFClusterAlgo;
24 
25 namespace reco {
26 
47  class PFCluster : public CaloCluster {
48  public:
49 
50  typedef std::vector<std::pair<CaloClusterPtr::key_type,edm::Ptr<PFCluster> > > EEtoPSAssociation;
51  // Next typedef uses double in ROOT 6 rather than Double32_t due to a bug in ROOT 5,
52  // which otherwise would make ROOT5 files unreadable in ROOT6. This does not increase
53  // the size on disk, because due to the bug, double was actually stored on disk in ROOT 5.
54  typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double> > REPPoint;
55 
57 
60  double x, double y, double z );
61 
63  void reset();
64 
66  void resetHitsAndFractions();
67 
70 
72  const std::vector< reco::PFRecHitFraction >& recHitFractions() const
73  { return rechits_; }
74 
76  void setLayer( PFLayer::Layer layer);
77 
79  PFLayer::Layer layer() const;
80 
82  double energy() const {return energy_;}
83 
85  float time() const {return time_;}
87  float timeError() const { return timeError_; }
88 
90  double depth() const {return depth_;}
91 
92  void setTime(float time, float timeError=0) {time_ = time; timeError_ = timeError; }
94  void setDepth(double depth) {depth_ = depth;}
95 
97  const REPPoint& positionREP() const {return posrep_;}
98 
101  posrep_.SetCoordinates( position_.Rho(),
102  position_.Eta(),
103  position_.Phi() );
104  }
105 
107  static double getDepthCorrection(double energy, bool isBelowPS = false,
108  bool isHadron = false);
109 
111  void setColor(int color) {color_ = color;}
112 
114  int color() const {return color_;}
115 
116 
117  PFCluster& operator=(const PFCluster&);
118 
120  static void setDepthCorParameters(int mode,
121  double a, double b,
122  double ap, double bp ) {
124  depthCorA_ = a;
125  depthCorB_ = b;
126  depthCorAp_ = ap;
127  depthCorBp_ = bp;
128  }
129 
130 
134 
136  double charge() const { return 0;}
137 
139  double pt() const {
140  return (energy() * sin(position_.theta()));
141  }
142 
144  double theta() const {
145  return position_.theta();
146  }
147 
149  math::XYZPoint const & vertex() const {
150  return dummyVtx_;
151  }
152  double vx() const { return vertex().x(); }
153  double vy() const { return vertex().y(); }
154  double vz() const { return vertex().z(); }
155 
156 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
157  template<typename pruner>
158  void pruneUsing(pruner prune) {
159  // remove_if+erase algo applied to both vectors...
160  auto iter = std::find_if_not(rechits_.begin(),rechits_.end(),prune);
161  if (iter==rechits_.end()) return;
162  auto first = iter-rechits_.begin();
163  for (auto i=first; ++i<int(rechits_.size());) {
164  if (prune(rechits_[i])) {
167  ++first;
168  }
169  }
170  rechits_.erase(rechits_.begin()+first,rechits_.end());
172  }
173 #endif
174 
175  private:
176 
178  std::vector< reco::PFRecHitFraction > rechits_;
179 
181  REPPoint posrep_;
182 
185  double depth_;
186 
189 
190 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
191  static std::atomic<int> depthCorMode_;
193 
195  static std::atomic<double> depthCorA_;
196 
198  static std::atomic<double> depthCorB_ ;
199 
201  static std::atomic<double> depthCorAp_;
202 
204  static std::atomic<double> depthCorBp_;
205 #else
206  static int depthCorMode_;
208 
210  static double depthCorA_;
211 
213  static double depthCorB_ ;
214 
216  static double depthCorAp_;
217 
219  static double depthCorBp_;
220 #endif
221 
222  static const math::XYZPoint dummyVtx_;
223 
225  int color_;
226  };
227 
228  std::ostream& operator<<(std::ostream& out,
229  const PFCluster& cluster);
230 
231 
232 }
233 
234 #endif
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:128
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:120
math::XYZPoint position_
cluster centroid position
Definition: CaloCluster.h:217
double energy_
cluster energy
Definition: CaloCluster.h:212
static std::atomic< double > depthCorAp_
Definition: PFCluster.h:201
float time() const
Definition: PFCluster.h:85
void setColor(int color)
set cluster color (for the PFRootEventManager display)
Definition: PFCluster.h:111
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
float timeError_
Definition: PFCluster.h:184
friend std::ostream & operator<<(std::ostream &out, const CaloCluster &cluster)
print me
static void setDepthCorParameters(int mode, double a, double b, double ap, double bp)
Definition: PFCluster.h:120
std::vector< reco::PFRecHitFraction > rechits_
vector of rechit fractions (transient)
Definition: PFCluster.h:178
PFLayer::Layer layer_
transient layer
Definition: PFCluster.h:188
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void pruneUsing(pruner prune)
Definition: PFCluster.h:158
void setTime(float time, float timeError=0)
Definition: PFCluster.h:92
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:139
PFCluster & operator=(const PFCluster &)
Definition: PFCluster.cc:136
std::vector< std::pair< DetId, float > > hitsAndFractions_
Definition: CaloCluster.h:223
static std::atomic< double > depthCorB_
Definition: PFCluster.h:198
double z() const
z coordinate of cluster centroid
Definition: CaloCluster.h:165
double vx() const
Definition: PFCluster.h:152
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:100
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:97
std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
Definition: PFCluster.h:50
static const math::XYZPoint dummyVtx_
Definition: PFCluster.h:222
layer definition for PFRecHit and PFCluster
Definition: PFLayer.h:21
double x() const
x coordinate of cluster centroid
Definition: CaloCluster.h:159
double vy() const
Definition: PFCluster.h:153
float timeError() const
Definition: PFCluster.h:87
void reset()
resets clusters parameters
Definition: PFCluster.cc:79
double energy() const
cluster energy
Definition: PFCluster.h:82
void resetHitsAndFractions()
reset only hits and fractions
Definition: PFCluster.cc:92
Layer
layer definition
Definition: PFLayer.h:31
static double getDepthCorrection(double energy, bool isBelowPS=false, bool isHadron=false)
Definition: PFCluster.cc:108
void setTimeError(float timeError)
Definition: PFCluster.h:93
REPPoint posrep_
cluster position: rho, eta, phi (transient)
Definition: PFCluster.h:181
double vz() const
Definition: PFCluster.h:154
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
static std::atomic< int > depthCorMode_
Definition: PFCluster.h:192
double b
Definition: hdecay.h:120
double theta() const
angle
Definition: PFCluster.h:144
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:99
int color() const
Definition: PFCluster.h:114
static std::atomic< double > depthCorBp_
Definition: PFCluster.h:204
fixed size matrix
double a
Definition: hdecay.h:121
math::XYZPoint const & vertex() const
dummy vertex access
Definition: PFCluster.h:149
int color_
color (transient)
Definition: PFCluster.h:225
float time_
Michalis :Add timing and depth information.
Definition: PFCluster.h:184
void setDepth(double depth)
Definition: PFCluster.h:94
static std::atomic< double > depthCorA_
Definition: PFCluster.h:195
double charge() const
dummy charge
Definition: PFCluster.h:136
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:54
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
double y() const
y coordinate of cluster centroid
Definition: CaloCluster.h:162
def move(src, dest)
Definition: eostools.py:510
double depth() const
cluster depth
Definition: PFCluster.h:90