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 
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 
150  double depth_;
151 
154 
156  static const constexpr double depthCorA_ = 0.89;
157  static const constexpr double depthCorB_ = 7.3;
158  static const constexpr double depthCorAp_ = 0.89;
159  static const constexpr double depthCorBp_ = 4.0;
160 
161  static const math::XYZPoint dummyVtx_;
162  };
163 
164  std::ostream& operator<<(std::ostream& out, const PFCluster& cluster);
165 
166 } // namespace reco
167 
168 #endif // DataFormats_ParticleFlowReco_PFCluster_h
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:123
mps_fire.i
i
Definition: mps_fire.py:428
reco::PFCluster::setTime
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
reco::PFCluster::PFCluster
PFCluster()
Definition: PFCluster.h:50
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:118
reco::PFCluster::layer
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:56
reco::PFCluster::vx
double vx() const
Definition: PFCluster.h:117
reco::PFCluster::vz
double vz() const
Definition: PFCluster.h:119
reco::PFCluster::setLayer
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:49
reco::PFCluster::recHitFractions
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
reco::PFCluster::time
float time() const
Definition: PFCluster.h:77
reco::PFCluster::charge
double charge() const
dummy charge
Definition: PFCluster.h:107
reco::CaloCluster::z
double z() const
z coordinate of cluster centroid
Definition: CaloCluster.h:178
reco::PFCluster::timeError
float timeError() const
Definition: PFCluster.h:79
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
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:153
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:74
DivergingColor.frac
float frac
Definition: DivergingColor.py:175
reco::PFCluster::depthCorAp_
static const constexpr double depthCorAp_
Definition: PFCluster.h:158
reco::PFCluster::operator=
PFCluster & operator=(const PFCluster &)
Definition: PFCluster.cc:63
reco::CaloCluster
Definition: CaloCluster.h:31
reco::PFCluster::depth
double depth() const
cluster depth
Definition: PFCluster.h:82
CaloClusterFwd.h
reco::PFCluster::resetHitsAndFractions
void resetHitsAndFractions()
reset only hits and fractions
Definition: PFCluster.cc:28
reco::PFCluster::depth_
double depth_
Definition: PFCluster.h:150
reco::PFCluster::positionREP
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:92
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:112
reco::PFCluster::rechits_
std::vector< reco::PFRecHitFraction > rechits_
vector of rechit fractions (transient)
Definition: PFCluster.h:143
PFLayer::Layer
Layer
layer definition
Definition: PFLayer.h:29
reco::PFCluster::depthCorB_
static const constexpr double depthCorB_
Definition: PFCluster.h:157
reco::PFCluster::depthCorA_
static const constexpr double depthCorA_
depth corrections
Definition: PFCluster.h:156
reco::PFCluster::reset
void reset()
resets clusters parameters
Definition: PFCluster.cc:17
reco::PFCluster::theta
double theta() const
angle
Definition: PFCluster.h:113
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
reco::PFCluster::dummyVtx_
static const math::XYZPoint dummyVtx_
Definition: PFCluster.h:161
reco::PFCluster::getDepthCorrection
static double getDepthCorrection(double energy, bool isBelowPS=false, bool isHadron=false)
Definition: PFCluster.cc:39
reco::PFCluster::calculatePositionREP
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:95
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
reco::PFCluster::depthCorBp_
static const constexpr double depthCorBp_
Definition: PFCluster.h:159
PFRecHitFraction.h
reco::PFCluster::addRecHitFraction
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
reco::PFCluster::EEtoPSAssociation
std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
Definition: PFCluster.h:44
eostools.move
def move(src, dest)
Definition: eostools.py:511
reco::PFCluster::pt
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:110
reco::PFCluster::time_
float time_
Michalis: add timing and depth information.
Definition: PFCluster.h:149
reco::PFCluster::vertex
math::XYZPoint const & vertex() const
dummy vertex access
Definition: PFCluster.h:116
reco::CaloCluster::energy_
double energy_
cluster energy
Definition: CaloCluster.h:223
reco::PFCluster
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
reco::PFCluster::setTimeError
void setTimeError(float timeError)
Definition: PFCluster.h:88
Point3D.h
reco::PFCluster::posrep_
REPPoint posrep_
cluster position: rho, eta, phi (transient)
Definition: PFCluster.h:146
NONE
Definition: TkAlStyle.cc:47
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
reco::PFCluster::setDepth
void setDepth(double depth)
Definition: PFCluster.h:89
reco::PFCluster::REPPoint
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
reco::CaloCluster::x
double x() const
x coordinate of cluster centroid
Definition: CaloCluster.h:172
reco::PFCluster::timeError_
float timeError_
Definition: PFCluster.h:149
CaloCluster.h