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 class PFClusterAlgo;
23 
24 namespace reco {
25 
46  class PFCluster : public CaloCluster {
47  public:
48  typedef std::vector<std::pair<CaloClusterPtr::key_type, edm::Ptr<PFCluster> > > EEtoPSAssociation;
49  // Next typedef uses double in ROOT 6 rather than Double32_t due to a bug in ROOT 5,
50  // which otherwise would make ROOT5 files unreadable in ROOT6. This does not increase
51  // the size on disk, because due to the bug, double was actually stored on disk in ROOT 5.
52  typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double> > REPPoint;
53 
55 
57  PFCluster(PFLayer::Layer layer, double energy, double x, double y, double z);
58 
60  void reset();
61 
63  void resetHitsAndFractions();
64 
67 
69  const std::vector<reco::PFRecHitFraction>& recHitFractions() const { return rechits_; }
70 
73 
75  PFLayer::Layer layer() const;
76 
78  double energy() const { return energy_; }
79 
81  float time() const { return time_; }
83  float timeError() const { return timeError_; }
84 
86  double depth() const { return depth_; }
87 
88  void setTime(float time, float timeError = 0) {
89  time_ = time;
91  }
93  void setDepth(double depth) { depth_ = depth; }
94 
96  const REPPoint& positionREP() const { return posrep_; }
97 
99  void calculatePositionREP() { posrep_.SetCoordinates(position_.Rho(), position_.Eta(), position_.Phi()); }
100 
102  static double getDepthCorrection(double energy, bool isBelowPS = false, bool isHadron = false);
103 
105  void setColor(int color) { color_ = color; }
106 
108  int color() const { return color_; }
109 
110  PFCluster& operator=(const PFCluster&);
111 
113  static void setDepthCorParameters(int mode, double a, double b, double ap, double bp) {
115  depthCorA_ = a;
116  depthCorB_ = b;
117  depthCorAp_ = ap;
118  depthCorBp_ = bp;
119  }
120 
124 
126  double charge() const { return 0; }
127 
129  double pt() const { return (energy() * sin(position_.theta())); }
130 
132  double theta() const { return position_.theta(); }
133 
135  math::XYZPoint const& vertex() const { return dummyVtx_; }
136  double vx() const { return vertex().x(); }
137  double vy() const { return vertex().y(); }
138  double vz() const { return vertex().z(); }
139 
140 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
141  template <typename pruner>
142  void pruneUsing(pruner prune) {
143  // remove_if+erase algo applied to both vectors...
144  auto iter = std::find_if_not(rechits_.begin(), rechits_.end(), prune);
145  if (iter == rechits_.end())
146  return;
147  auto first = iter - rechits_.begin();
148  for (auto i = first; ++i < int(rechits_.size());) {
149  if (prune(rechits_[i])) {
152  ++first;
153  }
154  }
155  rechits_.erase(rechits_.begin() + first, rechits_.end());
157  }
158 #endif
159 
160  private:
162  std::vector<reco::PFRecHitFraction> rechits_;
163 
166 
169  double depth_;
170 
173 
174 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
175  static std::atomic<int> depthCorMode_;
177 
179  static std::atomic<double> depthCorA_;
180 
182  static std::atomic<double> depthCorB_;
183 
185  static std::atomic<double> depthCorAp_;
186 
188  static std::atomic<double> depthCorBp_;
189 #else
190  static int depthCorMode_;
192 
194  static double depthCorA_;
195 
197  static double depthCorB_;
198 
200  static double depthCorAp_;
201 
203  static double depthCorBp_;
204 #endif
205 
206  static const math::XYZPoint dummyVtx_;
207 
209  int color_;
210  };
211 
212  std::ostream& operator<<(std::ostream& out, const PFCluster& cluster);
213 
214 } // namespace reco
215 
216 #endif
reco::CaloCluster::y
double y() const
y coordinate of cluster centroid
Definition: CaloCluster.h:175
reco::PFCluster::pruneUsing
void pruneUsing(pruner prune)
Definition: PFCluster.h:142
mps_fire.i
i
Definition: mps_fire.py:355
reco::PFCluster::setTime
void setTime(float time, float timeError=0)
Definition: PFCluster.h:88
reco::PFCluster::setDepthCorParameters
static void setDepthCorParameters(int mode, double a, double b, double ap, double bp)
Definition: PFCluster.h:113
reco::PFCluster::PFCluster
PFCluster()
Definition: PFCluster.h:54
reco::PFRecHitFraction
Fraction of a PFRecHit (rechits can be shared between several PFCluster's)
Definition: PFRecHitFraction.h:18
reco::PFCluster::vy
double vy() const
Definition: PFCluster.h:137
reco::PFCluster::layer
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:100
reco::PFCluster::vx
double vx() const
Definition: PFCluster.h:136
cropTnPTrees.frac
frac
Definition: cropTnPTrees.py:18
reco::PFCluster::depthCorAp_
static std::atomic< double > depthCorAp_
Definition: PFCluster.h:185
reco::PFCluster::vz
double vz() const
Definition: PFCluster.h:138
ALCARECOPromptCalibProdSiPixelAli0T_cff.mode
mode
Definition: ALCARECOPromptCalibProdSiPixelAli0T_cff.py:96
reco::PFCluster::setLayer
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:93
reco::PFCluster::recHitFractions
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:69
reco::PFCluster::time
float time() const
Definition: PFCluster.h:81
reco::PFCluster::charge
double charge() const
dummy charge
Definition: PFCluster.h:126
reco::CaloCluster::z
double z() const
z coordinate of cluster centroid
Definition: CaloCluster.h:178
reco::PFCluster::depthCorB_
static std::atomic< double > depthCorB_
Definition: PFCluster.h:182
reco::PFCluster::timeError
float timeError() const
Definition: PFCluster.h:83
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
reco::PFCluster::setColor
void setColor(int color)
set cluster color (for the PFRootEventManager display)
Definition: PFCluster.h:105
dqmdumpme.first
first
Definition: dqmdumpme.py:55
reco::CaloCluster::hitsAndFractions_
std::vector< std::pair< DetId, float > > hitsAndFractions_
Definition: CaloCluster.h:234
PFRecHit.h
reco::PFCluster::layer_
PFLayer::Layer layer_
transient layer
Definition: PFCluster.h:172
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
PFLayer.h
reco::CaloCluster::position_
math::XYZPoint position_
cluster centroid position
Definition: CaloCluster.h:228
reco::PFCluster::energy
double energy() const
cluster energy
Definition: PFCluster.h:78
reco::PFCluster::operator=
PFCluster & operator=(const PFCluster &)
Definition: PFCluster.cc:107
reco::CaloCluster
Definition: CaloCluster.h:31
reco::PFCluster::depth
double depth() const
cluster depth
Definition: PFCluster.h:86
CaloClusterFwd.h
reco::PFCluster::resetHitsAndFractions
void resetHitsAndFractions()
reset only hits and fractions
Definition: PFCluster.cc:72
reco::PFCluster::depth_
double depth_
Definition: PFCluster.h:169
reco::PFCluster::positionREP
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:96
b
double b
Definition: hdecay.h:118
reco::PFCluster::rechits_
std::vector< reco::PFRecHitFraction > rechits_
vector of rechit fractions (transient)
Definition: PFCluster.h:162
PFLayer::Layer
Layer
layer definition
Definition: PFLayer.h:29
reco::PFCluster::reset
void reset()
resets clusters parameters
Definition: PFCluster.cc:61
reco::PFCluster::theta
double theta() const
angle
Definition: PFCluster.h:132
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
a
double a
Definition: hdecay.h:119
reco::PFCluster::dummyVtx_
static const math::XYZPoint dummyVtx_
Definition: PFCluster.h:206
reco::PFCluster::depthCorMode_
static std::atomic< int > depthCorMode_
Definition: PFCluster.h:176
reco::PFCluster::getDepthCorrection
static double getDepthCorrection(double energy, bool isBelowPS=false, bool isHadron=false)
Definition: PFCluster.cc:83
reco::PFCluster::calculatePositionREP
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:99
createfilelist.int
int
Definition: createfilelist.py:10
reco::CaloCluster::particleFlow
Definition: CaloCluster.h:39
reco::operator<<
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
PFLayer
layer definition for PFRecHit and PFCluster
Definition: PFLayer.h:20
PFRecHitFraction.h
reco::PFCluster::addRecHitFraction
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:77
reco::PFCluster::EEtoPSAssociation
std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
Definition: PFCluster.h:48
eostools.move
def move(src, dest)
Definition: eostools.py:511
reco::PFCluster::pt
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:129
reco::PFCluster::time_
float time_
Michalis :Add timing and depth information.
Definition: PFCluster.h:168
reco::PFCluster::vertex
math::XYZPoint const & vertex() const
dummy vertex access
Definition: PFCluster.h:135
reco::CaloCluster::energy_
double energy_
cluster energy
Definition: CaloCluster.h:223
reco::PFCluster
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:46
reco::PFCluster::setTimeError
void setTimeError(float timeError)
Definition: PFCluster.h:92
Point3D.h
reco::PFCluster::posrep_
REPPoint posrep_
cluster position: rho, eta, phi (transient)
Definition: PFCluster.h:165
NONE
Definition: TkAlStyle.cc:47
reco::PFCluster::color_
int color_
color (transient)
Definition: PFCluster.h:209
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
reco::PFCluster::setDepth
void setDepth(double depth)
Definition: PFCluster.h:93
reco::PFCluster::depthCorA_
static std::atomic< double > depthCorA_
Definition: PFCluster.h:179
reco::PFCluster::REPPoint
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:52
reco::PFCluster::color
int color() const
Definition: PFCluster.h:108
reco::PFCluster::depthCorBp_
static std::atomic< double > depthCorBp_
Definition: PFCluster.h:188
reco::CaloCluster::x
double x() const
x coordinate of cluster centroid
Definition: CaloCluster.h:172
reco::PFCluster::timeError_
float timeError_
Definition: PFCluster.h:168
CaloCluster.h