1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
|
#ifndef DataFormats_ParticleFlowReco_PFCluster_h
#define DataFormats_ParticleFlowReco_PFCluster_h
#include <iostream>
#include <vector>
#include <algorithm>
#include <Rtypes.h>
#include "DataFormats/CaloRecHit/interface/CaloCluster.h"
#include "DataFormats/CaloRecHit/interface/CaloClusterFwd.h"
#include "DataFormats/Math/interface/Point3D.h"
#include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
#include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
#include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
#include "Math/GenVector/PositionVector3D.h"
class PFClusterAlgo;
namespace reco {
/**\class PFCluster
\brief Particle flow cluster, see clustering algorithm in PFClusterAlgo
A particle flow cluster is defined by its energy and position, which are
calculated from a vector of PFRecHitFraction. This calculation is
performed in PFClusterAlgo.
\todo Clean up this class to a common base (talk to Paolo Meridiani)
the extra internal stuff (like the vector of PFRecHitFraction's)
could be moved to a PFClusterExtra.
\todo Now that PFRecHitFraction's hold a reference to the PFRecHit's,
put back the calculation of energy and position to PFCluster.
\todo Add an operator+=
\author Colin Bernet
\date July 2006
*/
class PFCluster : public CaloCluster {
public:
typedef std::vector<std::pair<CaloClusterPtr::key_type, edm::Ptr<PFCluster>>> EEtoPSAssociation;
// Next typedef uses double in ROOT 6 rather than Double32_t due to a bug in ROOT 5,
// which otherwise would make ROOT5 files unreadable in ROOT6. This does not increase
// the size on disk, because due to the bug, double was actually stored on disk in ROOT 5.
typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double>> REPPoint;
PFCluster() : CaloCluster(CaloCluster::particleFlow), time_(-99.0), layer_(PFLayer::NONE) {}
/// constructor
PFCluster(PFLayer::Layer layer, double energy, double x, double y, double z);
/// resets clusters parameters
void reset();
/// reset only hits and fractions
void resetHitsAndFractions();
/// add a given fraction of the rechit
void addRecHitFraction(const reco::PFRecHitFraction& frac);
/// vector of rechit fractions
const std::vector<reco::PFRecHitFraction>& recHitFractions() const { return rechits_; }
/// set layer
void setLayer(PFLayer::Layer layer);
/// cluster layer, see PFLayer.h in this directory
PFLayer::Layer layer() const;
/// cluster energy
double energy() const { return energy_; }
/// \return cluster time
float time() const { return time_; }
/// \return the timing uncertainty
float timeError() const { return timeError_; }
/// cluster depth
double depth() const { return depth_; }
void setTime(float time, float timeError = 0) {
time_ = time;
timeError_ = timeError;
}
void setTimeError(float timeError) { timeError_ = timeError; }
void setDepth(double depth) { depth_ = depth; }
/// cluster position: rho, eta, phi
const REPPoint& positionREP() const { return posrep_; }
/// computes posrep_ once and for all
void calculatePositionREP() { posrep_.SetCoordinates(position_.Rho(), position_.Eta(), position_.Phi()); }
/// \todo move to PFClusterTools
static double getDepthCorrection(double energy, bool isBelowPS = false, bool isHadron = false);
/// some classes to make this fit into a template footprint
/// for RecoPFClusterRefCandidate so we can make jets and MET
/// out of PFClusters.
/// dummy charge
double charge() const { return 0; }
/// transverse momentum, massless approximation
double pt() const { return (energy() * sin(position_.theta())); }
/// angle
double theta() const { return position_.theta(); }
/// dummy vertex access
math::XYZPoint const& vertex() const { return dummyVtx_; }
double vx() const { return vertex().x(); }
double vy() const { return vertex().y(); }
double vz() const { return vertex().z(); }
#if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
template <typename pruner>
void pruneUsing(pruner prune) {
// remove_if+erase algo applied to both vectors...
auto iter = std::find_if_not(rechits_.begin(), rechits_.end(), prune);
if (iter == rechits_.end())
return;
auto first = iter - rechits_.begin();
for (auto i = first; ++i < int(rechits_.size());) {
if (prune(rechits_[i])) {
rechits_[first] = std::move(rechits_[i]);
hitsAndFractions_[first] = std::move(hitsAndFractions_[i]);
++first;
}
}
rechits_.erase(rechits_.begin() + first, rechits_.end());
hitsAndFractions_.erase(hitsAndFractions_.begin() + first, hitsAndFractions_.end());
}
#endif
private:
/// vector of rechit fractions (transient)
std::vector<reco::PFRecHitFraction> rechits_;
/// cluster position: rho, eta, phi (transient)
REPPoint posrep_;
/// Michalis: add timing and depth information
float time_;
float timeError_{0};
double depth_{0};
/// transient layer
PFLayer::Layer layer_;
/// depth corrections
static const constexpr double depthCorA_ = 0.89;
static const constexpr double depthCorB_ = 7.3;
static const constexpr double depthCorAp_ = 0.89;
static const constexpr double depthCorBp_ = 4.0;
static const math::XYZPoint dummyVtx_;
};
std::ostream& operator<<(std::ostream& out, const PFCluster& cluster);
} // namespace reco
#endif // DataFormats_ParticleFlowReco_PFCluster_h
|