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
|
/*
* See header file for a description of this class.
*
* \author Paolo Ronchese INFN Padova
*
*/
//-----------------------
// This Class' Header --
//-----------------------
#include "HeavyFlavorAnalysis/SpecificDecay/interface/BPHDecayToV0DiffMassBuilder.h"
//-------------------------------
// Collaborating Class Headers --
//-------------------------------
#include "HeavyFlavorAnalysis/SpecificDecay/interface/BPHDecayGenericBuilderBase.h"
#include "HeavyFlavorAnalysis/SpecificDecay/interface/BPHDecayToTkpTknSymChargeBuilder.h"
#include "HeavyFlavorAnalysis/RecoDecay/interface/BPHRecoBuilder.h"
#include "HeavyFlavorAnalysis/RecoDecay/interface/BPHPlusMinusCandidate.h"
#include "DataFormats/Candidate/interface/Candidate.h"
//---------------
// C++ Headers --
//---------------
#include <cmath>
using namespace std;
//-------------------
// Initializations --
//-------------------
//----------------
// Constructors --
//----------------
BPHDecayToV0DiffMassBuilder::BPHDecayToV0DiffMassBuilder(const BPHEventSetupWrapper& es,
const string& daug1Name,
double daug1Mass,
double daug1Sigma,
const string& daug2Name,
double daug2Mass,
double daug2Sigma,
const BPHRecoBuilder::BPHGenericCollection* daug1Collection,
const BPHRecoBuilder::BPHGenericCollection* daug2Collection,
double expectedMass)
: BPHDecayGenericBuilderBase(es),
BPHDecayToV0Builder(es, daug1Name, daug2Name, daug1Collection, daug2Collection),
BPHDecayToTkpTknSymChargeBuilder(es,
daug1Name,
daug1Mass,
daug1Sigma,
daug2Name,
daug2Mass,
daug2Sigma,
daug1Collection,
daug2Collection,
expectedMass),
p1Mass(daug1Mass),
p2Mass(daug2Mass),
p1Sigma(daug1Sigma),
p2Sigma(daug2Sigma),
expMass(expectedMass) {}
BPHDecayToV0DiffMassBuilder::BPHDecayToV0DiffMassBuilder(const BPHEventSetupWrapper& es,
const string& daug1Name,
double daug1Mass,
double daug1Sigma,
const string& daug2Name,
double daug2Mass,
double daug2Sigma,
const vector<reco::VertexCompositeCandidate>* v0Collection,
double expectedMass,
const string& searchList)
: BPHDecayGenericBuilderBase(es),
BPHDecayToV0Builder(es, daug1Name, daug2Name, v0Collection, searchList),
BPHDecayToTkpTknSymChargeBuilder(
es, daug1Name, daug1Mass, daug1Sigma, daug2Name, daug2Mass, daug2Sigma, nullptr, nullptr, expectedMass),
p1Mass(daug1Mass),
p2Mass(daug2Mass),
p1Sigma(daug1Sigma),
p2Sigma(daug2Sigma),
expMass(expectedMass) {}
BPHDecayToV0DiffMassBuilder::BPHDecayToV0DiffMassBuilder(const BPHEventSetupWrapper& es,
const string& daug1Name,
double daug1Mass,
double daug1Sigma,
const string& daug2Name,
double daug2Mass,
double daug2Sigma,
const vector<reco::VertexCompositePtrCandidate>* vpCollection,
double expectedMass,
const string& searchList)
: BPHDecayGenericBuilderBase(es),
BPHDecayToV0Builder(es, daug1Name, daug2Name, vpCollection, searchList),
BPHDecayToTkpTknSymChargeBuilder(
es, daug1Name, daug1Mass, daug1Sigma, daug2Name, daug2Mass, daug2Sigma, nullptr, nullptr, expectedMass),
p1Mass(daug1Mass),
p2Mass(daug2Mass),
p1Sigma(daug1Sigma),
p2Sigma(daug2Sigma),
expMass(expectedMass) {}
//--------------
// Operations --
//--------------
void BPHDecayToV0DiffMassBuilder::buildFromBPHGenericCollection() {
BPHDecayToTkpTknSymChargeBuilder::build();
return;
}
BPHPlusMinusCandidatePtr BPHDecayToV0DiffMassBuilder::buildCandidate(const reco::Candidate* c1,
const reco::Candidate* c2,
const void* v0,
v0Type type) {
BPHPlusMinusCandidatePtr candX = BPHPlusMinusCandidateWrap::create(evSetup);
BPHPlusMinusCandidatePtr candY = BPHPlusMinusCandidateWrap::create(evSetup);
BPHPlusMinusCandidate* cptrX = candX.get();
BPHPlusMinusCandidate* cptrY = candY.get();
cptrX->add(p1Name, c1, sList, p1Mass, p1Sigma);
cptrX->add(p2Name, c2, sList, p2Mass, p2Sigma);
cptrY->add(p1Name, c2, sList, p1Mass, p1Sigma);
cptrY->add(p2Name, c1, sList, p2Mass, p2Sigma);
double mv0 = 0.0;
switch (type) {
case VertexCompositeCandidate:
mv0 = static_cast<const reco::VertexCompositeCandidate*>(v0)->mass();
break;
case VertexCompositePtrCandidate:
mv0 = static_cast<const reco::VertexCompositePtrCandidate*>(v0)->mass();
break;
default:
mv0 = expMass;
break;
}
double m1 = 0.0;
double m2 = 0.0;
if (p1Mass > p2Mass) {
m1 = c1->mass();
m2 = c2->mass();
} else {
m1 = c2->mass();
m2 = c1->mass();
}
// check daughter masses in V0 CompositeCandidate
double mcut = (p1Mass + p2Mass) / 2;
if ((m1 > mcut) && (m2 < mcut))
return candX;
if ((m1 < mcut) && (m2 > mcut))
return candY;
// choose combination having the best invariant mass
return (fabs(mv0 - cptrX->mass()) < fabs(mv0 - cptrY->mass()) ? candX : candY);
}
|