ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::BPhysBGammaFinder Class Reference

#include <BPhysBGammaFinder.h>

Inheritance diagram for DerivationFramework::BPhysBGammaFinder:
Collaboration diagram for DerivationFramework::BPhysBGammaFinder:

Public Member Functions

 BPhysBGammaFinder (const std::string &t, const std::string &n, const IInterface *p)
StatusCode initialize () override
StatusCode finalize () override
virtual StatusCode addBranches (const EventContext &ctx) const override
TVector3 trackMomentum (const xAOD::Vertex &vxCandidate, int trkIndex) const

Private Attributes

SG::ReadHandleKeyArray< xAOD::VertexContainerm_BVertexCollectionsToCheck
SG::ReadDecorHandleKeyArray< xAOD::VertexContainerm_passFlagsToCheck
ToolHandle< Trk::V0Toolsm_v0Tools
ToolHandle< Trk::IVertexFitterm_vertexFitter
ToolHandle< InDet::VertexPointEstimatorm_vertexEstimator
SG::ReadHandleKey< xAOD::TrackParticleContainerm_inputTrackParticleContainerName
SG::ReadHandleKey< xAOD::TrackParticleContainerm_inputLowPtTrackContainerName
SG::WriteHandleKey< xAOD::VertexContainerm_conversionContainerName
float m_maxDeltaQ
float m_Chi2Cut
float m_maxGammaMass

Detailed Description

Definition at line 46 of file BPhysBGammaFinder.h.

Constructor & Destructor Documentation

◆ BPhysBGammaFinder()

DerivationFramework::BPhysBGammaFinder::BPhysBGammaFinder ( const std::string & t,
const std::string & n,
const IInterface * p )

Definition at line 23 of file BPhysBGammaFinder.cxx.

24 : base_class(t,n,p),
25 m_v0Tools("Trk::V0Tools"),
26 m_vertexFitter("Trk::TrkVKalVrtFitter", this),
27 m_vertexEstimator("InDet::VertexPointEstimator", this),
28 m_inputTrackParticleContainerName("InDetTrackParticles"),
29 m_inputLowPtTrackContainerName("LowPtRoITrackParticles"),
30 m_conversionContainerName("BPhysConversionCandidates"),
31 m_maxDeltaQ(700.0),
32 m_Chi2Cut(20.0),
33 m_maxGammaMass(100.0) {
34
35
36 // Declare user-defined properties
37 declareProperty("BVertexContainers", m_BVertexCollectionsToCheck);
38 declareProperty("PassFlagsToCheck", m_passFlagsToCheck);
39 declareProperty("V0Tools", m_v0Tools);
40 declareProperty("VertexFitterTool", m_vertexFitter);
41 declareProperty("VertexEstimator", m_vertexEstimator);
42 declareProperty("InputTrackParticleContainerName", m_inputTrackParticleContainerName);
43 declareProperty("InputLowPtTrackContainerName", m_inputLowPtTrackContainerName);
44 declareProperty("ConversionContainerName", m_conversionContainerName);
45 declareProperty("MaxDeltaQ", m_maxDeltaQ = 700.0); // Maximum mass difference between di-muon+conversion and di-muon
46 declareProperty("Chi2Cut", m_Chi2Cut = 20.0);
47 declareProperty("MaxGammaMass", m_maxGammaMass = 100.0);
48}
SG::WriteHandleKey< xAOD::VertexContainer > m_conversionContainerName
SG::ReadHandleKeyArray< xAOD::VertexContainer > m_BVertexCollectionsToCheck
SG::ReadDecorHandleKeyArray< xAOD::VertexContainer > m_passFlagsToCheck
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_inputLowPtTrackContainerName
ToolHandle< Trk::IVertexFitter > m_vertexFitter
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_inputTrackParticleContainerName
ToolHandle< InDet::VertexPointEstimator > m_vertexEstimator

Member Function Documentation

◆ addBranches()

StatusCode DerivationFramework::BPhysBGammaFinder::addBranches ( const EventContext & ctx) const
overridevirtual

Definition at line 71 of file BPhysBGammaFinder.cxx.

71 {
72
73 std::vector<const xAOD::Vertex*> BVertices;
74 std::vector<const xAOD::TrackParticle*> BVertexTracks;
75
76
77 // Output conversion container
78 std::unique_ptr<xAOD::VertexContainer> conversionContainer(new xAOD::VertexContainer());
79 std::unique_ptr<xAOD::VertexAuxContainer> conversionAuxContainer(new xAOD::VertexAuxContainer());
80 conversionContainer->setStore(conversionAuxContainer.get());
81
82 // Retrieve track particles from StoreGate
83 SG::ReadHandle<xAOD::TrackParticleContainer> inputTrackParticles{m_inputTrackParticleContainerName, ctx};
84 ATH_MSG_DEBUG( "Track particle container size " << inputTrackParticles->size() );
85
86 std::vector<const xAOD::TrackParticle*> trackPair(2);
87
88 static const SG::Decorator< std::vector< VertexLink > > BGammaLinks( "BGammaLinks" );
89
90 auto flaghandles = m_passFlagsToCheck.makeHandles(ctx);
91 // Retrieve vertex containers
92 for (SG::ReadHandle<xAOD::VertexContainer>& BVtxContainer : m_BVertexCollectionsToCheck.makeHandles(ctx)) {
93 if (!BVtxContainer.isValid()) { // replaces the CHECK
94 ATH_MSG_ERROR("Failed to retrieve VertexContainer " << BVtxContainer.key());
95 return StatusCode::FAILURE; // or ATH_CHECK if you prefer the macro style
96 }
97
98 ATH_MSG_DEBUG( "Vertex Container (" << BVtxContainer.key() << ") contains " << BVtxContainer->size() << " vertices" );
99
100
101 for (const xAOD::Vertex* vertex : *BVtxContainer) {
102 auto &vect = BGammaLinks(*vertex) = std::vector< VertexLink >();
103
104 bool passedHypothesis = false;
105 BVertexTracks.clear();
106
107 for (const auto &flag : flaghandles) {
108 bool pass = flag(*vertex);
109 if (pass) passedHypothesis = true;
110 }
111
112 if (!passedHypothesis) continue;
113 xAOD::BPhysHypoHelper Bc("Bc", vertex);
114
115 // link to Bc+ vertex
116 std::vector<const xAOD::Vertex*> precedingVertices(1, vertex);
117
118 // Collect up B-vertex tracks
119 for (size_t i = 0; i < vertex->nTrackParticles(); ++i) BVertexTracks.push_back(vertex->trackParticle(i));
120
121 // Track Selection
122 // Track1 Loop
123 for (xAOD::TrackParticleContainer::const_iterator tpIt1 = inputTrackParticles->begin(); tpIt1 != inputTrackParticles->end(); ++tpIt1) {
124 trackPair[0] = *tpIt1;
125
126 auto itr1 = std::find(BVertexTracks.begin(), BVertexTracks.end(), trackPair[0]);
127 if (itr1 != BVertexTracks.end()) continue;
128
129 const Trk::Perigee& trackPerigee1 = trackPair[0]->perigeeParameters();
130
131 // Track2 Loop
132 for (xAOD::TrackParticleContainer::const_iterator tpIt2 = tpIt1 + 1; tpIt2 != inputTrackParticles->end(); ++tpIt2) {
133 trackPair[1] = *tpIt2;
134 if (trackPair[0] == trackPair[1]) continue;
135
136 auto itr2 = std::find(BVertexTracks.begin(), BVertexTracks.end(), trackPair[1]);
137 if (itr2 != BVertexTracks.end()) continue;
138
139 const Trk::Perigee& trackPerigee2 = trackPair[1]->perigeeParameters();
140
141 // Track pair selection
142 TLorentzVector e1, e2, gamma_m, BcStar;
143 e1.SetPtEtaPhiM(trackPair[0]->pt(), trackPair[0]->eta(), trackPair[0]->phi(), Trk::electron);
144 e2.SetPtEtaPhiM(trackPair[1]->pt(), trackPair[1]->eta(), trackPair[1]->phi(), Trk::electron);
145
146 gamma_m = e1 + e2;
147 if (gamma_m.M() > m_maxGammaMass) continue;
148
149 TLorentzVector mu1 = Bc.refTrk(0, Trk::muon);
150 TLorentzVector mu2 = Bc.refTrk(1, Trk::muon);
151 TLorentzVector mu3 = Bc.refTrk(2, Trk::muon);
152
153 BcStar = mu1 + mu2 + mu3 + e1 + e2;
154 double Q = BcStar.M() - Bc.mass() - 2 * Trk::electron;
155 if (Q > m_maxDeltaQ) continue;
156
157 // Estimate starting point + cuts on compatiblity of tracks
158 int sflag = 0;
159 int errorcode = 0;
160 Amg::Vector3D startingPoint = m_vertexEstimator->getCirclesIntersectionPoint(&trackPerigee1, &trackPerigee2, sflag, errorcode);
161 if (errorcode != 0) startingPoint = Amg::Vector3D::Zero(3);
162
163 std::vector<float> RefTrackPx, RefTrackPy, RefTrackPz, RefTrackE;
164 std::vector<float> OrigTrackPx, OrigTrackPy, OrigTrackPz, OrigTrackE;
165
166
167 // Do the vertex fit
168 auto convVertexCandidate = m_vertexFitter->fit(ctx, trackPair, startingPoint);
169
170 // Check for successful fit
171 if (convVertexCandidate) {
172 if (convVertexCandidate->chiSquared() / convVertexCandidate->numberDoF() > m_Chi2Cut) continue;
173
174 xAOD::BPhysHelper Photon(convVertexCandidate.get());
175 // set link to the parent Bc+ vertex
176 Photon.setPrecedingVertices(precedingVertices, BVtxContainer.cptr());
177
178 // Parameters at vertex
179 convVertexCandidate->clearTracks();
180 ElementLink<xAOD::TrackParticleContainer> newLink1;
181 newLink1.setElement(*tpIt1);
182 newLink1.setStorableObject(*inputTrackParticles);
183 ElementLink<xAOD::TrackParticleContainer> newLink2;
184 newLink2.setElement(*tpIt2);
185 newLink2.setStorableObject(*inputTrackParticles);
186 convVertexCandidate->addTrackAtVertex(newLink1);
187 convVertexCandidate->addTrackAtVertex(newLink2);
188
189 std::vector<Amg::Vector3D> positionList;
190
191 //Get photon momentum 3-vector
192 Amg::Vector3D momentum = m_v0Tools->V0Momentum(convVertexCandidate.get());
193
194 TLorentzVector photon, electron1, electron2, ph;
195 electron1.SetVectM( trackMomentum( *convVertexCandidate, 0 ), Trk::electron );
196 electron2.SetVectM( trackMomentum( *convVertexCandidate, 1 ), Trk::electron );
197 photon = electron1 + electron2;
198 ph.SetXYZM(momentum.x(), momentum.y(), momentum.z(), 0.);
199
200 // Use to keep track of which dimuon(s) gave a chi_c/b candidate
201 static const SG::Accessor<std::vector<float> > RefTrackPxAcc("RefTrackPx");
202 static const SG::Accessor<std::vector<float> > RefTrackPyAcc("RefTrackPy");
203 static const SG::Accessor<std::vector<float> > RefTrackPzAcc("RefTrackPz");
204 std::vector<float> B_Px = RefTrackPxAcc(*vertex);
205 std::vector<float> B_Py = RefTrackPyAcc(*vertex);
206 std::vector<float> B_Pz = RefTrackPzAcc(*vertex);
207
208 TLorentzVector muon1, muon2, muon3;
209 muon1.SetXYZM(B_Px.at(0), B_Py.at(0), B_Pz.at(0), Trk::muon);
210 muon2.SetXYZM(B_Px.at(1), B_Py.at(1), B_Pz.at(1), Trk::muon);
211 muon3.SetXYZM(B_Px.at(2), B_Py.at(2), B_Pz.at(2), Trk::muon);
212
213 TLorentzVector B_m = muon1 + muon2 + muon3;
214
215 const double deltaQ = (B_m + photon).M() - Bc.mass() - 2 * Trk::electron;
216 const double mass = photon.M();
217
218 RefTrackPx.push_back(trackMomentum(*convVertexCandidate, 0).Px());
219 RefTrackPx.push_back(trackMomentum(*convVertexCandidate, 1).Px());
220
221 RefTrackPy.push_back(trackMomentum(*convVertexCandidate, 0).Py());
222 RefTrackPy.push_back(trackMomentum(*convVertexCandidate, 1).Py());
223
224 RefTrackPz.push_back(trackMomentum(*convVertexCandidate, 0).Pz());
225 RefTrackPz.push_back(trackMomentum(*convVertexCandidate, 1).Pz());
226
227 for (size_t i = 0; i < B_Px.size(); i++) {
228 RefTrackPx.push_back(B_Px.at(i));
229 RefTrackPy.push_back(B_Py.at(i));
230 RefTrackPz.push_back(B_Pz.at(i));
231 }
232
233 RefTrackE.push_back(electron1.E());
234 RefTrackE.push_back(electron2.E());
235 RefTrackE.push_back(muon1.E());
236 RefTrackE.push_back(muon2.E());
237 RefTrackE.push_back(muon3.E());
238
239 OrigTrackPx.push_back(e1.Px());
240 OrigTrackPx.push_back(e2.Px());
241
242 OrigTrackPy.push_back(e1.Py());
243 OrigTrackPy.push_back(e2.Py());
244
245 OrigTrackPz.push_back(e1.Pz());
246 OrigTrackPz.push_back(e2.Pz());
247
248 OrigTrackE.push_back(e1.E());
249 OrigTrackE.push_back(e2.E());
250
251
252 ATH_MSG_DEBUG( "pt = " << photon.Pt() << " ph " << ph.Pt() << " mass " << photon.M() << " px size " << RefTrackPx.size() );
253 ATH_MSG_DEBUG( "Candidate DeltaM = " << (B_m + photon).M() << " MeV DiMuon " << " ( Mass = " << B_m.M() << " MeV )");
254
255 // Decorate selected conversions
256 ATH_MSG_DEBUG( "Decorating conversion vertices" );
257
258 static const SG::Accessor<float> pxAcc("px");
259 static const SG::Accessor<float> pyAcc("py");
260 static const SG::Accessor<float> pzAcc("pz");
261 pxAcc(*convVertexCandidate) = momentum.x();
262 pyAcc(*convVertexCandidate) = momentum.y();
263 pzAcc(*convVertexCandidate) = momentum.z();
264
265 static const SG::Accessor<float> deltaQAcc("deltaQ");
266 static const SG::Accessor<float> gamma_massAcc("gamma_mass");
267 static const SG::Accessor< std::vector<float> > RefTrackEAcc("RefTrackE");
268 deltaQAcc(*convVertexCandidate) = deltaQ;
269 gamma_massAcc(*convVertexCandidate) = mass;
270 RefTrackPxAcc(*convVertexCandidate) = RefTrackPx;
271 RefTrackPyAcc(*convVertexCandidate) = RefTrackPy;
272 RefTrackPzAcc(*convVertexCandidate) = RefTrackPz;
273 RefTrackEAcc(*convVertexCandidate) = RefTrackE;
274
275 static const SG::Accessor< std::vector<float> > OrigTrackPxAcc("OrigTrackPx");
276 static const SG::Accessor< std::vector<float> > OrigTrackPyAcc("OrigTrackPy");
277 static const SG::Accessor< std::vector<float> > OrigTrackPzAcc("OrigTrackPz");
278 static const SG::Accessor< std::vector<float> > OrigTrackEAcc("OrigTrackE");
279 OrigTrackPxAcc(*convVertexCandidate) = OrigTrackPx;
280 OrigTrackPyAcc(*convVertexCandidate) = OrigTrackPy;
281 OrigTrackPzAcc(*convVertexCandidate) = OrigTrackPz;
282 OrigTrackEAcc(*convVertexCandidate) = OrigTrackE;
283
284 static const SG::Accessor<Char_t> passed_GammaAcc("passed_Gamma");
285 passed_GammaAcc(*convVertexCandidate) = true; // Used in event skimming
286
287
288 // add cross-link to the original Bc+ vertex
289 VertexLink BGammaLink;
290 BGammaLink.setElement(convVertexCandidate.get());
291 BGammaLink.setStorableObject(*conversionContainer);
292 conversionContainer->push_back( std::move(convVertexCandidate) );
293 vect.push_back(std::move(BGammaLink));
294 }
295 else {
296 ATH_MSG_DEBUG( "Vertex Fit Failed" );
297 }
298
299 } // end of Track2 Loop
300 } // end of Track1 Loop
301
302 } // end of Bc loop
303
304 } // end of vertex container loop
305
306 SG::WriteHandle<xAOD::VertexContainer> wh(m_conversionContainerName, ctx);
307 ATH_CHECK( wh.record(std::move(conversionContainer), std::move(conversionAuxContainer)) );
308
309 return StatusCode::SUCCESS;
310}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
TVector3 trackMomentum(const xAOD::Vertex &vxCandidate, int trkIndex) const
Eigen::Matrix< double, 3, 1 > Vector3D
ElementLink< xAOD::VertexContainer > VertexLink
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:573
SG::Decorator< T, ALLOC > Decorator
Helper class to provide type-safe access to aux data, specialized for JaggedVecElt.
Definition AuxElement.h:576
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
double e2(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 2nd sampling
double e1(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 1st sampling
bool flag
Definition master.py:29
str wh
Definition parseDir.py:45
VertexAuxContainer_v1 VertexAuxContainer
Definition of the current jet auxiliary container.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
Photon_v1 Photon
Definition of the current "egamma version".

◆ finalize()

StatusCode DerivationFramework::BPhysBGammaFinder::finalize ( )
override

Definition at line 66 of file BPhysBGammaFinder.cxx.

66 {
67 return StatusCode::SUCCESS;
68}

◆ initialize()

StatusCode DerivationFramework::BPhysBGammaFinder::initialize ( )
override

Definition at line 51 of file BPhysBGammaFinder.cxx.

51 {
52
53 ATH_MSG_DEBUG("in initialize()");
54
55 ATH_CHECK( m_v0Tools.retrieve() );
56 ATH_CHECK( m_vertexFitter.retrieve() );
57 ATH_CHECK( m_vertexEstimator.retrieve() );
62 return StatusCode::SUCCESS;
63}

◆ trackMomentum()

TVector3 DerivationFramework::BPhysBGammaFinder::trackMomentum ( const xAOD::Vertex & vxCandidate,
int trkIndex ) const

Definition at line 314 of file BPhysBGammaFinder.cxx.

314 {
315
316 double px = 0.;
317 double py = 0.;
318 double pz = 0.;
319 const Trk::TrackParameters* aPerigee = vxCandidate.vxTrackAtVertex()[trkIndex].perigeeAtVertex();
320 px = aPerigee->momentum()[Trk::px];
321 py = aPerigee->momentum()[Trk::py];
322 pz = aPerigee->momentum()[Trk::pz];
323
324 return TVector3(px,py,pz);
325}
const Amg::Vector3D & momentum() const
Access method for the momentum.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
@ pz
global momentum (cartesian)
Definition ParamDefs.h:61
@ px
Definition ParamDefs.h:59
@ py
Definition ParamDefs.h:60
ParametersBase< TrackParametersDim, Charged > TrackParameters

Member Data Documentation

◆ m_BVertexCollectionsToCheck

SG::ReadHandleKeyArray<xAOD::VertexContainer> DerivationFramework::BPhysBGammaFinder::m_BVertexCollectionsToCheck
private

Definition at line 60 of file BPhysBGammaFinder.h.

◆ m_Chi2Cut

float DerivationFramework::BPhysBGammaFinder::m_Chi2Cut
private

Definition at line 72 of file BPhysBGammaFinder.h.

◆ m_conversionContainerName

SG::WriteHandleKey<xAOD::VertexContainer> DerivationFramework::BPhysBGammaFinder::m_conversionContainerName
private

Definition at line 69 of file BPhysBGammaFinder.h.

◆ m_inputLowPtTrackContainerName

SG::ReadHandleKey<xAOD::TrackParticleContainer> DerivationFramework::BPhysBGammaFinder::m_inputLowPtTrackContainerName
private

Definition at line 68 of file BPhysBGammaFinder.h.

◆ m_inputTrackParticleContainerName

SG::ReadHandleKey<xAOD::TrackParticleContainer> DerivationFramework::BPhysBGammaFinder::m_inputTrackParticleContainerName
private

Definition at line 67 of file BPhysBGammaFinder.h.

◆ m_maxDeltaQ

float DerivationFramework::BPhysBGammaFinder::m_maxDeltaQ
private

Definition at line 71 of file BPhysBGammaFinder.h.

◆ m_maxGammaMass

float DerivationFramework::BPhysBGammaFinder::m_maxGammaMass
private

Definition at line 73 of file BPhysBGammaFinder.h.

◆ m_passFlagsToCheck

SG::ReadDecorHandleKeyArray<xAOD::VertexContainer> DerivationFramework::BPhysBGammaFinder::m_passFlagsToCheck
private

Definition at line 61 of file BPhysBGammaFinder.h.

◆ m_v0Tools

ToolHandle<Trk::V0Tools> DerivationFramework::BPhysBGammaFinder::m_v0Tools
private

Definition at line 63 of file BPhysBGammaFinder.h.

◆ m_vertexEstimator

ToolHandle<InDet::VertexPointEstimator> DerivationFramework::BPhysBGammaFinder::m_vertexEstimator
private

Definition at line 65 of file BPhysBGammaFinder.h.

◆ m_vertexFitter

ToolHandle<Trk::IVertexFitter> DerivationFramework::BPhysBGammaFinder::m_vertexFitter
private

Definition at line 64 of file BPhysBGammaFinder.h.


The documentation for this class was generated from the following files: