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
|
#ifndef DataFormats_Candidate_VertexCompositePtrCandidate_H
#define DataFormats_Candidate_VertexCompositePtrCandidate_H
/** \class reco::VertexCompositePtrCandidate
*
* A composite Candidate with error
* matrix and other vertex fix information.
*
* \author Luca Lista, INFN
*
*
*/
#include "DataFormats/Candidate/interface/VertexCompositePtrCandidateFwd.h"
#include "DataFormats/Candidate/interface/CompositePtrCandidate.h"
namespace reco {
class VertexCompositePtrCandidate : public CompositePtrCandidate {
public:
enum { dimension4D = 4 };
/// covariance error matrix (3x3)
typedef math::Error<dimension4D>::type CovarianceMatrix4D;
/// matix size
enum { size4D = dimension4D * (dimension4D + 1) / 2 };
VertexCompositePtrCandidate() : CompositePtrCandidate(), chi2_(0), ndof_(0), time_(0) {}
/// constructor from values
VertexCompositePtrCandidate(
Charge q, const LorentzVector &p4, const Point &vtx, int pdgId = 0, int status = 0, bool integerCharge = true)
: CompositePtrCandidate(q, p4, vtx, pdgId, status, integerCharge), chi2_(0), ndof_(0), time_(0) {}
VertexCompositePtrCandidate(Charge q,
const LorentzVector &p4,
const Point &vtx,
double time,
int pdgId = 0,
int status = 0,
bool integerCharge = true)
: CompositePtrCandidate(q, p4, vtx, pdgId, status, integerCharge), chi2_(0), ndof_(0), time_(time) {}
/// constructor from values
VertexCompositePtrCandidate(Charge q,
const LorentzVector &p4,
const Point &vtx,
const CovarianceMatrix &err,
double chi2,
double ndof,
int pdgId = 0,
int status = 0,
bool integerCharge = true);
VertexCompositePtrCandidate(Charge q,
const LorentzVector &p4,
const Point &vtx,
double time,
const CovarianceMatrix4D &err,
double chi2,
double ndof,
int pdgId = 0,
int status = 0,
bool integerCharge = true);
/// constructor from values
explicit VertexCompositePtrCandidate(const Candidate &p) : CompositePtrCandidate(p), chi2_(0), ndof_(0), time_(0) {}
/// constructor from values
explicit VertexCompositePtrCandidate(const CompositePtrCandidate &p)
: CompositePtrCandidate(p), chi2_(0), ndof_(0), time_(0) {}
/// destructor
~VertexCompositePtrCandidate() override;
/// returns a clone of the candidate
VertexCompositePtrCandidate *clone() const override;
/// chi-squares
double vertexChi2() const override { return chi2_; }
/** Number of degrees of freedom
* Meant to be Double32_t for soft-assignment fitters:
* tracks may contribute to the vertex with fractional weights.
* The ndof is then = to the sum of the track weights.
* see e.g. CMS NOTE-2006/032, CMS NOTE-2004/002
*/
double vertexNdof() const override { return ndof_; }
/// chi-squared divided by n.d.o.f.
double vertexNormalizedChi2() const override { return chi2_ / ndof_; }
/// (i, j)-th element of error matrix, i, j = 0, ... 3
double vertexCovariance(int i, int j) const override { return covariance_[idx(i, j)]; }
using reco::LeafCandidate::vertexCovariance; // avoid hiding the
/// return SMatrix 4D
CovarianceMatrix4D vertexCovariance4D() const {
CovarianceMatrix4D m;
fillVertexCovariance(m);
return m;
}
/// fill SMatrix
void fillVertexCovariance(CovarianceMatrix &v) const override;
/// 4D version
void fillVertexCovariance(CovarianceMatrix4D &v) const;
/// set chi2 and ndof
void setChi2AndNdof(double chi2, double ndof) {
chi2_ = chi2;
ndof_ = ndof;
}
/// set covariance matrix
void setCovariance(const CovarianceMatrix &m);
/// set covariance matrix
void setCovariance(const CovarianceMatrix4D &m);
// set time
void setTime(double time) { time_ = time; }
/// the following functions are implemented to have a more consistent interface with the one of reco::Vertex
typedef math::Error<dimension>::type Error;
typedef math::Error<dimension4D>::type Error4D;
const Point &position() const { return vertex(); }
double t() const { return time_; }
double tError() const { return std::sqrt(vertexCovariance(3, 3)); }
Error error() const {
Error m;
fillVertexCovariance(m);
return m;
}
/// return SMatrix 4D
Error4D error4D() const {
Error4D m;
fillVertexCovariance(m);
return m;
}
private:
/// chi-sqared
Double32_t chi2_;
/// number of degrees of freedom
Double32_t ndof_;
/// covariance matrix (4x4) as vector
Double32_t covariance_[size4D];
/// vertex time
Double32_t time_;
/// position index
index idx(index i, index j) const {
int a = (i <= j ? i : j), b = (i <= j ? j : i);
return b * (b + 1) / 2 + a;
}
};
} // namespace reco
#endif
|