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

#include <JpsiPlusDpstCascade.h>

Inheritance diagram for DerivationFramework::JpsiPlusDpstCascade:
Collaboration diagram for DerivationFramework::JpsiPlusDpstCascade:

Public Member Functions

 JpsiPlusDpstCascade (const std::string &t, const std::string &n, const IInterface *p)
 ~JpsiPlusDpstCascade ()
virtual StatusCode initialize () override
StatusCode performSearch (std::vector< Trk::VxCascadeInfo * > *cascadeinfoContainer, const EventContext &ctx) const
virtual StatusCode addBranches (const EventContext &ctx) const override

Private Attributes

std::string m_vertexContainerKey
std::string m_vertexD0ContainerKey
std::vector< std::string > m_cascadeOutputsKeys
std::string m_VxPrimaryCandidateName
 Name of primary vertex container // FIXME Use Handles.
double m_jpsiMassLower
double m_jpsiMassUpper
double m_jpsipiMassLower
double m_jpsipiMassUpper
double m_D0MassLower
double m_D0MassUpper
double m_DstMassLower
double m_DstMassUpper
double m_MassLower
double m_MassUpper
double m_vtx0MassHypo
double m_vtx1MassHypo
double m_vtx0Daug1MassHypo
double m_vtx0Daug2MassHypo
double m_vtx0Daug3MassHypo
double m_vtx1Daug1MassHypo
double m_vtx1Daug2MassHypo
double m_mass_jpsi
int m_Dx_pid
bool m_constrD0
bool m_constrJpsi
double m_chi2cut
SG::ReadHandleKey< xAOD::EventInfom_eventInfo_key {this, "EventInfo", "EventInfo", "Input event information"}
ToolHandle< Trk::TrkVKalVrtFitterm_iVertexFitter
ToolHandle< Analysis::PrimaryVertexRefitterm_pvRefitter
ToolHandle< Trk::V0Toolsm_V0Tools
ToolHandle< DerivationFramework::CascadeToolsm_CascadeTools
ServiceHandle< IPartPropSvc > m_partPropSvc {this, "PartPropSvc", "PartPropSvc"}
bool m_refitPV
std::string m_refPVContainerName
std::string m_hypoName
 name of the mass hypothesis.
int m_PV_max
int m_DoVertexType
size_t m_PV_minNTracks

Detailed Description

Definition at line 36 of file JpsiPlusDpstCascade.h.

Constructor & Destructor Documentation

◆ JpsiPlusDpstCascade()

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

Definition at line 397 of file JpsiPlusDpstCascade.cxx.

397 : base_class(t,n,p),
400 m_cascadeOutputsKeys{ "JpsiPlusDpstCascadeVtx1", "JpsiPlusDpstCascadeVtx2" },
401 m_VxPrimaryCandidateName("PrimaryVertices"),
402 m_jpsiMassLower(0.0),
403 m_jpsiMassUpper(10000.0),
405 m_jpsipiMassUpper(10000.0),
406 m_D0MassLower(0.0),
407 m_D0MassUpper(10000.0),
408 m_DstMassLower(0.0),
409 m_DstMassUpper(10000.0),
410 m_MassLower(0.0),
411 m_MassUpper(20000.0),
412 m_vtx0MassHypo(-1),
413 m_vtx1MassHypo(-1),
419 m_mass_jpsi(-1),
420 m_Dx_pid(421),
421 m_constrD0(true),
422 m_constrJpsi(true),
423 m_chi2cut(-1.0),
424 m_iVertexFitter("Trk::TrkVKalVrtFitter"),
425 m_pvRefitter("Analysis::PrimaryVertexRefitter"),
426 m_V0Tools("Trk::V0Tools"),
427 m_CascadeTools("DerivationFramework::CascadeTools")
428 {
429 declareProperty("JpsipiVertices", m_vertexContainerKey);
430 declareProperty("D0Vertices", m_vertexD0ContainerKey);
431 declareProperty("VxPrimaryCandidateName", m_VxPrimaryCandidateName);
432 declareProperty("RefPVContainerName", m_refPVContainerName = "RefittedPrimaryVertices");
433 declareProperty("JpsiMassLowerCut", m_jpsiMassLower);
434 declareProperty("JpsiMassUpperCut", m_jpsiMassUpper);
435 declareProperty("JpsipiMassLowerCut", m_jpsipiMassLower);
436 declareProperty("JpsipiMassUpperCut", m_jpsipiMassUpper);
437 declareProperty("D0MassLowerCut", m_D0MassLower);
438 declareProperty("D0MassUpperCut", m_D0MassUpper);
439 declareProperty("DstMassLowerCut", m_DstMassLower);
440 declareProperty("DstMassUpperCut", m_DstMassUpper);
441 declareProperty("MassLowerCut", m_MassLower);
442 declareProperty("MassUpperCut", m_MassUpper);
443 declareProperty("HypothesisName", m_hypoName = "Bc");
444 declareProperty("Vtx0MassHypo", m_vtx0MassHypo);
445 declareProperty("Vtx1MassHypo", m_vtx1MassHypo);
446 declareProperty("Vtx0Daug1MassHypo", m_vtx0Daug1MassHypo);
447 declareProperty("Vtx0Daug2MassHypo", m_vtx0Daug2MassHypo);
448 declareProperty("Vtx0Daug3MassHypo", m_vtx0Daug3MassHypo);
449 declareProperty("Vtx1Daug1MassHypo", m_vtx1Daug1MassHypo);
450 declareProperty("Vtx1Daug2MassHypo", m_vtx1Daug2MassHypo);
451 declareProperty("JpsiMass", m_mass_jpsi);
452 declareProperty("DxHypothesis", m_Dx_pid);
453 declareProperty("ApplyD0MassConstraint", m_constrD0);
454 declareProperty("ApplyJpsiMassConstraint", m_constrJpsi);
455 declareProperty("Chi2Cut", m_chi2cut);
456 declareProperty("RefitPV", m_refitPV = true);
457 declareProperty("MaxnPV", m_PV_max = 999);
458 declareProperty("MinNTracksInPV", m_PV_minNTracks = 0);
459 declareProperty("DoVertexType", m_DoVertexType = 7);
460 declareProperty("TrkVertexFitterTool", m_iVertexFitter);
461 declareProperty("PVRefitter", m_pvRefitter);
462 declareProperty("V0Tools", m_V0Tools);
463 declareProperty("CascadeTools", m_CascadeTools);
464 declareProperty("CascadeVertexCollections", m_cascadeOutputsKeys);
465 }
std::string m_VxPrimaryCandidateName
Name of primary vertex container // FIXME Use Handles.
ToolHandle< Trk::TrkVKalVrtFitter > m_iVertexFitter
std::vector< std::string > m_cascadeOutputsKeys
std::string m_hypoName
name of the mass hypothesis.
ToolHandle< Analysis::PrimaryVertexRefitter > m_pvRefitter
ToolHandle< DerivationFramework::CascadeTools > m_CascadeTools

◆ ~JpsiPlusDpstCascade()

DerivationFramework::JpsiPlusDpstCascade::~JpsiPlusDpstCascade ( )

Definition at line 467 of file JpsiPlusDpstCascade.cxx.

467{ }

Member Function Documentation

◆ addBranches()

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

Definition at line 60 of file JpsiPlusDpstCascade.cxx.

61 {
62 std::vector<Trk::VxCascadeInfo*> cascadeinfoContainer;
63 constexpr int topoN = 2;
64 std::array<xAOD::VertexContainer*, topoN> Vtxwritehandles;
65 std::array<xAOD::VertexAuxContainer*, topoN> Vtxwritehandlesaux;
66 if(m_cascadeOutputsKeys.size() !=topoN) { ATH_MSG_FATAL("Incorrect number of VtxContainers"); return StatusCode::FAILURE; }
67
68 for(int i =0; i<topoN;i++){
69 Vtxwritehandles[i] = new xAOD::VertexContainer();
70 Vtxwritehandlesaux[i] = new xAOD::VertexAuxContainer();
71 Vtxwritehandles[i]->setStore(Vtxwritehandlesaux[i]);
72 ATH_CHECK(evtStore()->record(Vtxwritehandles[i] , m_cascadeOutputsKeys[i] )); // FIXME Use Handles
73 ATH_CHECK(evtStore()->record(Vtxwritehandlesaux[i], m_cascadeOutputsKeys[i] + "Aux.")); // FIXME Use Handles
74 }
75
76 //----------------------------------------------------
77 // retrieve primary vertices
78 //----------------------------------------------------
79 const xAOD::Vertex * primaryVertex{};
80 const xAOD::VertexContainer *pvContainer{};
81 ATH_CHECK(evtStore()->retrieve(pvContainer, m_VxPrimaryCandidateName)); // FIXME Use Handles
82 ATH_MSG_DEBUG("Found " << m_VxPrimaryCandidateName << " in StoreGate!");
83
84 if (pvContainer->size()==0){
85 ATH_MSG_WARNING("You have no primary vertices: " << pvContainer->size());
86 return StatusCode::RECOVERABLE;
87 } else {
88 primaryVertex = (*pvContainer)[0];
89 }
90
91 //----------------------------------------------------
92 // Try to retrieve refitted primary vertices
93 //----------------------------------------------------
94 xAOD::VertexContainer* refPvContainer{};
95 xAOD::VertexAuxContainer* refPvAuxContainer{};
96 if (m_refitPV) {
98 // refitted PV container exists. Get it from the store gate
99 ATH_CHECK(evtStore()->retrieve(refPvContainer , m_refPVContainerName )); // FIXME Use Handles
100 ATH_CHECK(evtStore()->retrieve(refPvAuxContainer, m_refPVContainerName + "Aux.")); // FIXME Use Handles
101 } else {
102 // refitted PV container does not exist. Create a new one.
103 refPvContainer = new xAOD::VertexContainer;
104 refPvAuxContainer = new xAOD::VertexAuxContainer;
105 refPvContainer->setStore(refPvAuxContainer);
106 ATH_CHECK(evtStore()->record(refPvContainer , m_refPVContainerName)); // FIXME Use Handles
107 ATH_CHECK(evtStore()->record(refPvAuxContainer, m_refPVContainerName+"Aux.")); // FIXME Use Handles
108 }
109 }
110
111 ATH_CHECK(performSearch(&cascadeinfoContainer, ctx));
112
113 SG::ReadHandle<xAOD::EventInfo> evt(m_eventInfo_key, ctx);
114 if(not evt.isValid()) ATH_MSG_ERROR("Cannot Retrieve " << evt.key() );
115 BPhysPVCascadeTools helper(&(*m_CascadeTools), evt.cptr());
116 helper.SetMinNTracksInPV(m_PV_minNTracks);
117
118 // Decorators for the main vertex: chi2, ndf, pt and pt error, plus the D0 vertex variables
119 SG::AuxElement::Decorator<VertexLinkVector> CascadeLinksDecor("CascadeVertexLinks");
120 SG::AuxElement::Decorator<VertexLinkVector> JpsipiLinksDecor("JpsipiVertexLinks");
121 SG::AuxElement::Decorator<VertexLinkVector> D0LinksDecor("D0VertexLinks");
122 SG::AuxElement::Decorator<float> chi2_decor("ChiSquared");
123 SG::AuxElement::Decorator<float> ndof_decor("NumberDoF");
124 SG::AuxElement::Decorator<float> Pt_decor("Pt");
125 SG::AuxElement::Decorator<float> PtErr_decor("PtErr");
126 SG::AuxElement::Decorator<float> Mass_svdecor("D0_mass");
127 SG::AuxElement::Decorator<float> MassErr_svdecor("D0_massErr");
128 SG::AuxElement::Decorator<float> Pt_svdecor("D0_Pt");
129 SG::AuxElement::Decorator<float> PtErr_svdecor("D0_PtErr");
130 SG::AuxElement::Decorator<float> Lxy_svdecor("D0_Lxy");
131 SG::AuxElement::Decorator<float> LxyErr_svdecor("D0_LxyErr");
132 SG::AuxElement::Decorator<float> Tau_svdecor("D0_Tau");
133 SG::AuxElement::Decorator<float> TauErr_svdecor("D0_TauErr");
134
135 SG::AuxElement::Decorator<float> MassMumu_decor("Mumu_mass");
136 SG::AuxElement::Decorator<float> MassKpi_svdecor("Kpi_mass");
137 SG::AuxElement::Decorator<float> MassJpsi_decor("Jpsi_mass");
138 SG::AuxElement::Decorator<float> MassPiD0_decor("PiD0_mass");
139
140 ATH_MSG_DEBUG("cascadeinfoContainer size " << cascadeinfoContainer.size());
141
142 // Get Jpsi+pi container and identify the input Jpsi+pi
143 const xAOD::VertexContainer *jpsipiContainer{};
144 ATH_CHECK(evtStore()->retrieve(jpsipiContainer , m_vertexContainerKey )); // FIXME Use Handles
145 const xAOD::VertexContainer *d0Container{};
146 ATH_CHECK(evtStore()->retrieve(d0Container , m_vertexD0ContainerKey )); // FIXME Use Handles
147
148 for (Trk::VxCascadeInfo* x : cascadeinfoContainer) {
149 if(x==nullptr) {
150 ATH_MSG_ERROR("cascadeinfoContainer is null");
151 //x is dereferenced if we pass this
152 return StatusCode::FAILURE;
153 }
154
155 // the cascade fitter returns:
156 // std::vector<xAOD::Vertex*>, each xAOD::Vertex contains the refitted track parameters (perigee at the vertex position)
157 // vertices[iv] the links to the original TPs and a covariance of size 3+5*NTRK; the chi2 of the total fit
158 // is split between the cascade vertices as per track contribution
159 // std::vector< std::vector<TLorentzVector> >, each std::vector<TLorentzVector> contains the refitted momenta (TLorentzVector)
160 // momenta[iv][...] of all tracks in the corresponding vertex, including any pseudotracks (from cascade vertices)
161 // originating in this vertex; the masses are as assigned in the cascade fit
162 // std::vector<Amg::MatrixX>, the corresponding covariance matrices in momentum space
163 // covariance[iv]
164 // int nDoF, double Chi2
165 //
166 // the invariant mass, pt, lifetime etc. errors should be calculated using the covariance matrices in momentum space as these
167 // take into account the full track-track and track-vertex correlations
168 //
169 // in the case of Jpsi+V0: vertices[0] is the V0 vertex, vertices[1] is the B/Lambda_b(bar) vertex, containing the 2 Jpsi tracks.
170 // The covariance terms between the two vertices are not stored. In momentum space momenta[0] contains the 2 V0 tracks,
171 // their momenta add up to the momentum of the 3rd track in momenta[1], the first two being the Jpsi tracks
172
173 const std::vector<xAOD::Vertex*> &cascadeVertices = x->vertices();
174 if(cascadeVertices.size()!=topoN)
175 ATH_MSG_ERROR("Incorrect number of vertices");
176 if(cascadeVertices[0] == nullptr || cascadeVertices[1] == nullptr) ATH_MSG_ERROR("Error null vertex");
177 // Keep vertices (bear in mind that they come in reverse order!)
178 for(int i =0;i<topoN;i++) Vtxwritehandles[i]->push_back(cascadeVertices[i]);
179
180 x->setSVOwnership(false); // Prevent Container from deleting vertices
181 const auto mainVertex = cascadeVertices[1]; // this is the B_c+/- vertex
182 const std::vector< std::vector<TLorentzVector> > &moms = x->getParticleMoms();
183
184 // Set links to cascade vertices
185 std::vector<const xAOD::Vertex*> verticestoLink;
186 verticestoLink.push_back(cascadeVertices[0]);
187 if(Vtxwritehandles[1] == nullptr) ATH_MSG_ERROR("Vtxwritehandles[1] is null");
188 if(!BPhysPVCascadeTools::LinkVertices(CascadeLinksDecor, verticestoLink, Vtxwritehandles[0], cascadeVertices[1]))
189 ATH_MSG_ERROR("Error decorating with cascade vertices");
190
191 // Identify the input Jpsi+pi
192 const xAOD::Vertex* jpsipiVertex = BPhysPVCascadeTools::FindVertex<3>(jpsipiContainer, cascadeVertices[1]);
193 ATH_MSG_DEBUG("1 pt Jpsi+pi tracks " << cascadeVertices[1]->trackParticle(0)->pt() << ", " << cascadeVertices[1]->trackParticle(1)->pt() << ", " << cascadeVertices[1]->trackParticle(2)->pt());
194 if (jpsipiVertex) ATH_MSG_DEBUG("2 pt Jpsi+pi tracks " << jpsipiVertex->trackParticle(0)->pt() << ", " << jpsipiVertex->trackParticle(1)->pt() << ", " << jpsipiVertex->trackParticle(2)->pt());
195
196 // Identify the input D0
197 const xAOD::Vertex* d0Vertex = BPhysPVCascadeTools::FindVertex<2>(d0Container, cascadeVertices[0]);;
198 ATH_MSG_DEBUG("1 pt D0 tracks " << cascadeVertices[0]->trackParticle(0)->pt() << ", " << cascadeVertices[0]->trackParticle(1)->pt());
199 if (d0Vertex) ATH_MSG_DEBUG("2 pt D0 tracks " << d0Vertex->trackParticle(0)->pt() << ", " << d0Vertex->trackParticle(1)->pt());
200
201 // Set links to input vertices
202 std::vector<const xAOD::Vertex*> jpsipiVerticestoLink;
203 if (jpsipiVertex) jpsipiVerticestoLink.push_back(jpsipiVertex);
204 else ATH_MSG_WARNING("Could not find linking Jpsi+pi");
205 if(!BPhysPVCascadeTools::LinkVertices(JpsipiLinksDecor, jpsipiVerticestoLink, jpsipiContainer, cascadeVertices[1]))
206 ATH_MSG_ERROR("Error decorating with Jpsi+pi vertices");
207
208 std::vector<const xAOD::Vertex*> d0VerticestoLink;
209 if (d0Vertex) d0VerticestoLink.push_back(d0Vertex);
210 else ATH_MSG_WARNING("Could not find linking D0");
211 if(!BPhysPVCascadeTools::LinkVertices(D0LinksDecor, d0VerticestoLink, d0Container, cascadeVertices[1]))
212 ATH_MSG_ERROR("Error decorating with D0 vertices");
213
214 bool tagD0(true);
215 if (jpsipiVertex){
216 if(abs(m_Dx_pid)==421 && (jpsipiVertex->trackParticle(2)->charge()==-1)) tagD0 = false;
217 }
218
219 double mass_b = m_vtx0MassHypo;
220 double mass_d0 = m_vtx1MassHypo;
221 std::vector<double> massesJpsipi;
222 massesJpsipi.push_back(m_vtx0Daug1MassHypo);
223 massesJpsipi.push_back(m_vtx0Daug2MassHypo);
224 massesJpsipi.push_back(m_vtx0Daug3MassHypo);
225 std::vector<double> massesD0;
226 if(tagD0){
227 massesD0.push_back(m_vtx1Daug1MassHypo);
228 massesD0.push_back(m_vtx1Daug2MassHypo);
229 }else{ // Change the oreder of masses for D*-->D0bar pi-, D0bar->K+pi-
230 massesD0.push_back(m_vtx1Daug2MassHypo);
231 massesD0.push_back(m_vtx1Daug1MassHypo);
232 }
233 std::vector<double> Masses;
234 Masses.push_back(m_vtx0Daug1MassHypo);
235 Masses.push_back(m_vtx0Daug2MassHypo);
236 Masses.push_back(m_vtx0Daug3MassHypo);
237 Masses.push_back(m_vtx1MassHypo);
238
239 // loop over candidates -- Don't apply PV_minNTracks requirement here
240 // because it may result in exclusion of the high-pt PV.
241 // get good PVs
242
243 xAOD::BPhysHypoHelper vtx(m_hypoName, mainVertex);
244
245 // Get refitted track momenta from all vertices, charged tracks only
247
248 // Decorate main vertex
249 //
250 // 1.a) mass, mass error
251 BPHYS_CHECK( vtx.setMass(m_CascadeTools->invariantMass(moms[1])) );
252 BPHYS_CHECK( vtx.setMassErr(m_CascadeTools->invariantMassError(moms[1],x->getCovariance()[1])) );
253 // 1.b) pt and pT error (the default pt of mainVertex is != the pt of the full cascade fit!)
254 Pt_decor(*mainVertex) = m_CascadeTools->pT(moms[1]);
255 PtErr_decor(*mainVertex) = m_CascadeTools->pTError(moms[1],x->getCovariance()[1]);
256 // 1.c) chi2 and ndof (the default chi2 of mainVertex is != the chi2 of the full cascade fit!)
257 chi2_decor(*mainVertex) = x->fitChi2();
258 ndof_decor(*mainVertex) = x->nDoF();
259
260 float massMumu = 0.;
261 if (jpsipiVertex) {
262 TLorentzVector p4_mu1, p4_mu2;
263 p4_mu1.SetPtEtaPhiM(jpsipiVertex->trackParticle(0)->pt(),
264 jpsipiVertex->trackParticle(0)->eta(),
265 jpsipiVertex->trackParticle(0)->phi(), m_vtx0Daug1MassHypo);
266 p4_mu2.SetPtEtaPhiM(jpsipiVertex->trackParticle(1)->pt(),
267 jpsipiVertex->trackParticle(1)->eta(),
268 jpsipiVertex->trackParticle(1)->phi(), m_vtx0Daug2MassHypo);
269 massMumu = (p4_mu1 + p4_mu2).M();
270 }
271 MassMumu_decor(*mainVertex) = massMumu;
272
273 float massKpi = 0.;
274 if (d0Vertex) {
275 TLorentzVector p4_ka, p4_pi;
276 if(tagD0){
277 p4_pi.SetPtEtaPhiM(d0Vertex->trackParticle(0)->pt(),
278 d0Vertex->trackParticle(0)->eta(),
279 d0Vertex->trackParticle(0)->phi(), m_vtx1Daug1MassHypo);
280 p4_ka.SetPtEtaPhiM(d0Vertex->trackParticle(1)->pt(),
281 d0Vertex->trackParticle(1)->eta(),
282 d0Vertex->trackParticle(1)->phi(), m_vtx1Daug2MassHypo);
283 }else{ // Change the oreder of masses for D*-->D0bar pi-, D0bar->K+pi-
284 p4_pi.SetPtEtaPhiM(d0Vertex->trackParticle(1)->pt(),
285 d0Vertex->trackParticle(1)->eta(),
286 d0Vertex->trackParticle(1)->phi(), m_vtx1Daug1MassHypo);
287 p4_ka.SetPtEtaPhiM(d0Vertex->trackParticle(0)->pt(),
288 d0Vertex->trackParticle(0)->eta(),
289 d0Vertex->trackParticle(0)->phi(), m_vtx1Daug2MassHypo);
290 }
291 massKpi = (p4_ka + p4_pi).M();
292 }
293 MassKpi_svdecor(*mainVertex) = massKpi;
294 MassJpsi_decor(*mainVertex) = (moms[1][0] + moms[1][1]).M();
295 MassPiD0_decor(*mainVertex) = (moms[1][2] + moms[1][3]).M();
296
297
298 ATH_CHECK(helper.FillCandwithRefittedVertices(m_refitPV, pvContainer,
299 refPvContainer, &(*m_pvRefitter), m_PV_max, m_DoVertexType, x, 1, mass_b, vtx));
300
301 // 4) decorate the main vertex with D0 vertex mass, pt, lifetime and lxy values (plus errors)
302 // D0 points to the main vertex, so lifetime and lxy are w.r.t the main vertex
303 Mass_svdecor(*mainVertex) = m_CascadeTools->invariantMass(moms[0]);
304 MassErr_svdecor(*mainVertex) = m_CascadeTools->invariantMassError(moms[0],x->getCovariance()[0]);
305 Pt_svdecor(*mainVertex) = m_CascadeTools->pT(moms[0]);
306 PtErr_svdecor(*mainVertex) = m_CascadeTools->pTError(moms[0],x->getCovariance()[0]);
307 Lxy_svdecor(*mainVertex) = m_CascadeTools->lxy(moms[0],cascadeVertices[0],cascadeVertices[1]);
308 LxyErr_svdecor(*mainVertex) = m_CascadeTools->lxyError(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1]);
309 Tau_svdecor(*mainVertex) = m_CascadeTools->tau(moms[0],cascadeVertices[0],cascadeVertices[1]);
310 TauErr_svdecor(*mainVertex) = m_CascadeTools->tauError(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1]);
311
312 // Some checks in DEBUG mode
313 ATH_MSG_DEBUG("chi2 " << x->fitChi2()
314 << " chi2_1 " << m_V0Tools->chisq(cascadeVertices[0])
315 << " chi2_2 " << m_V0Tools->chisq(cascadeVertices[1])
316 << " vprob " << m_CascadeTools->vertexProbability(x->nDoF(),x->fitChi2()));
317 ATH_MSG_DEBUG("ndf " << x->nDoF() << " ndf_1 " << m_V0Tools->ndof(cascadeVertices[0]) << " ndf_2 " << m_V0Tools->ndof(cascadeVertices[1]));
318 ATH_MSG_DEBUG("V0Tools mass_d0 " << m_V0Tools->invariantMass(cascadeVertices[0],massesD0)
319 << " error " << m_V0Tools->invariantMassError(cascadeVertices[0],massesD0)
320 << " mass_J " << m_V0Tools->invariantMass(cascadeVertices[1],massesJpsipi)
321 << " error " << m_V0Tools->invariantMassError(cascadeVertices[1],massesJpsipi));
322 // masses and errors, using track masses assigned in the fit
323 double Mass_B = m_CascadeTools->invariantMass(moms[1]);
324 double Mass_D0 = m_CascadeTools->invariantMass(moms[0]);
325 double Mass_B_err = m_CascadeTools->invariantMassError(moms[1],x->getCovariance()[1]);
326 double Mass_D0_err = m_CascadeTools->invariantMassError(moms[0],x->getCovariance()[0]);
327 ATH_MSG_DEBUG("Mass_B " << Mass_B << " Mass_D0 " << Mass_D0);
328 ATH_MSG_DEBUG("Mass_B_err " << Mass_B_err << " Mass_D0_err " << Mass_D0_err);
329 double mprob_B = m_CascadeTools->massProbability(mass_b,Mass_B,Mass_B_err);
330 double mprob_D0 = m_CascadeTools->massProbability(mass_d0,Mass_D0,Mass_D0_err);
331 ATH_MSG_DEBUG("mprob_B " << mprob_B << " mprob_D0 " << mprob_D0);
332 // masses and errors, assigning user defined track masses
333 ATH_MSG_DEBUG("Mass_b " << m_CascadeTools->invariantMass(moms[1],Masses)
334 << " Mass_d0 " << m_CascadeTools->invariantMass(moms[0],massesD0));
335 ATH_MSG_DEBUG("Mass_b_err " << m_CascadeTools->invariantMassError(moms[1],x->getCovariance()[1],Masses)
336 << " Mass_d0_err " << m_CascadeTools->invariantMassError(moms[0],x->getCovariance()[0],massesD0));
337 ATH_MSG_DEBUG("pt_b " << m_CascadeTools->pT(moms[1])
338 << " pt_d " << m_CascadeTools->pT(moms[0])
339 << " pt_d0 " << m_V0Tools->pT(cascadeVertices[0]));
340 ATH_MSG_DEBUG("ptErr_b " << m_CascadeTools->pTError(moms[1],x->getCovariance()[1])
341 << " ptErr_d " << m_CascadeTools->pTError(moms[0],x->getCovariance()[0])
342 << " ptErr_d0 " << m_V0Tools->pTError(cascadeVertices[0]));
343 ATH_MSG_DEBUG("lxy_B " << m_V0Tools->lxy(cascadeVertices[1],primaryVertex) << " lxy_D " << m_V0Tools->lxy(cascadeVertices[0],cascadeVertices[1]));
344 ATH_MSG_DEBUG("lxy_b " << m_CascadeTools->lxy(moms[1],cascadeVertices[1],primaryVertex) << " lxy_d " << m_CascadeTools->lxy(moms[0],cascadeVertices[0],cascadeVertices[1]));
345 ATH_MSG_DEBUG("lxyErr_b " << m_CascadeTools->lxyError(moms[1],x->getCovariance()[1],cascadeVertices[1],primaryVertex)
346 << " lxyErr_d " << m_CascadeTools->lxyError(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1])
347 << " lxyErr_d0 " << m_V0Tools->lxyError(cascadeVertices[0],cascadeVertices[1]));
348 ATH_MSG_DEBUG("tau_B " << m_CascadeTools->tau(moms[1],cascadeVertices[1],primaryVertex,mass_b)
349 << " tau_d0 " << m_V0Tools->tau(cascadeVertices[0],cascadeVertices[1],massesD0));
350 ATH_MSG_DEBUG("tau_b " << m_CascadeTools->tau(moms[1],cascadeVertices[1],primaryVertex)
351 << " tau_d " << m_CascadeTools->tau(moms[0],cascadeVertices[0],cascadeVertices[1])
352 << " tau_D " << m_CascadeTools->tau(moms[0],cascadeVertices[0],cascadeVertices[1],mass_d0));
353 ATH_MSG_DEBUG("tauErr_b " << m_CascadeTools->tauError(moms[1],x->getCovariance()[1],cascadeVertices[1],primaryVertex)
354 << " tauErr_d " << m_CascadeTools->tauError(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1])
355 << " tauErr_d0 " << m_V0Tools->tauError(cascadeVertices[0],cascadeVertices[1],massesD0));
356 ATH_MSG_DEBUG("TauErr_b " << m_CascadeTools->tauError(moms[1],x->getCovariance()[1],cascadeVertices[1],primaryVertex,mass_b)
357 << " TauErr_d " << m_CascadeTools->tauError(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1],mass_d0)
358 << " TauErr_d0 " << m_V0Tools->tauError(cascadeVertices[0],cascadeVertices[1],massesD0,mass_d0));
359
360 ATH_MSG_DEBUG("CascadeTools main vert wrt PV " << " CascadeTools SV " << " V0Tools SV");
361 ATH_MSG_DEBUG("a0z " << m_CascadeTools->a0z(moms[1],cascadeVertices[1],primaryVertex)
362 << ", " << m_CascadeTools->a0z(moms[0],cascadeVertices[0],cascadeVertices[1])
363 << ", " << m_V0Tools->a0z(cascadeVertices[0],cascadeVertices[1]));
364 ATH_MSG_DEBUG("a0zErr " << m_CascadeTools->a0zError(moms[1],x->getCovariance()[1],cascadeVertices[1],primaryVertex)
365 << ", " << m_CascadeTools->a0zError(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1])
366 << ", " << m_V0Tools->a0zError(cascadeVertices[0],cascadeVertices[1]));
367 ATH_MSG_DEBUG("a0xy " << m_CascadeTools->a0xy(moms[1],cascadeVertices[1],primaryVertex)
368 << ", " << m_CascadeTools->a0xy(moms[0],cascadeVertices[0],cascadeVertices[1])
369 << ", " << m_V0Tools->a0xy(cascadeVertices[0],cascadeVertices[1]));
370 ATH_MSG_DEBUG("a0xyErr " << m_CascadeTools->a0xyError(moms[1],x->getCovariance()[1],cascadeVertices[1],primaryVertex)
371 << ", " << m_CascadeTools->a0xyError(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1])
372 << ", " << m_V0Tools->a0xyError(cascadeVertices[0],cascadeVertices[1]));
373 ATH_MSG_DEBUG("a0 " << m_CascadeTools->a0(moms[1],cascadeVertices[1],primaryVertex)
374 << ", " << m_CascadeTools->a0(moms[0],cascadeVertices[0],cascadeVertices[1])
375 << ", " << m_V0Tools->a0(cascadeVertices[0],cascadeVertices[1]));
376 ATH_MSG_DEBUG("a0Err " << m_CascadeTools->a0Error(moms[1],x->getCovariance()[1],cascadeVertices[1],primaryVertex)
377 << ", " << m_CascadeTools->a0Error(moms[0],x->getCovariance()[0],cascadeVertices[0],cascadeVertices[1])
378 << ", " << m_V0Tools->a0Error(cascadeVertices[0],cascadeVertices[1]));
379 ATH_MSG_DEBUG("x0 " << m_V0Tools->vtx(cascadeVertices[0]).x() << " y0 " << m_V0Tools->vtx(cascadeVertices[0]).y() << " z0 " << m_V0Tools->vtx(cascadeVertices[0]).z());
380 ATH_MSG_DEBUG("x1 " << m_V0Tools->vtx(cascadeVertices[1]).x() << " y1 " << m_V0Tools->vtx(cascadeVertices[1]).y() << " z1 " << m_V0Tools->vtx(cascadeVertices[1]).z());
381 ATH_MSG_DEBUG("X0 " << primaryVertex->x() << " Y0 " << primaryVertex->y() << " Z0 " << primaryVertex->z());
382 ATH_MSG_DEBUG("rxy0 " << m_V0Tools->rxy(cascadeVertices[0]) << " rxyErr0 " << m_V0Tools->rxyError(cascadeVertices[0]));
383 ATH_MSG_DEBUG("rxy1 " << m_V0Tools->rxy(cascadeVertices[1]) << " rxyErr1 " << m_V0Tools->rxyError(cascadeVertices[1]));
384 ATH_MSG_DEBUG("Rxy0 wrt PV " << m_V0Tools->rxy(cascadeVertices[0],primaryVertex) << " RxyErr0 wrt PV " << m_V0Tools->rxyError(cascadeVertices[0],primaryVertex));
385 ATH_MSG_DEBUG("Rxy1 wrt PV " << m_V0Tools->rxy(cascadeVertices[1],primaryVertex) << " RxyErr1 wrt PV " << m_V0Tools->rxyError(cascadeVertices[1],primaryVertex));
386 ATH_MSG_DEBUG("number of covariance matrices " << (x->getCovariance()).size());
387 } // loop over cascadeinfoContainer
388
389 // Deleting cascadeinfo since this won't be stored.
390 // Vertices have been kept in m_cascadeOutputs and should be owned by their container
391 for (auto x : cascadeinfoContainer) delete x;
392
393 return StatusCode::SUCCESS;
394 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define BPHYS_CHECK(EXP)
Useful CHECK macro.
#define x
size_type size() const noexcept
Returns the number of elements in the collection.
static bool LinkVertices(SG::AuxElement::Decorator< VertexLinkVector > &decor, const std::vector< const xAOD::Vertex * > &vertices, const xAOD::VertexContainer *vertexContainer, const xAOD::Vertex *vert)
static const xAOD::Vertex * FindVertex(const xAOD::VertexContainer *c, const xAOD::Vertex *v)
static void SetVectorInfo(xAOD::BPhysHelper &, const Trk::VxCascadeInfo *)
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfo_key
StatusCode performSearch(std::vector< Trk::VxCascadeInfo * > *cascadeinfoContainer, const EventContext &ctx) const
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
float charge() const
Returns the charge.
float z() const
Returns the z position.
const TrackParticle * trackParticle(size_t i) const
Get the pointer to a given track that was used in vertex reco.
float y() const
Returns the y position.
float x() const
Returns the x position.
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
retrieve(aClass, aKey=None)
Definition PyKernel.py:110
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.

◆ initialize()

StatusCode DerivationFramework::JpsiPlusDpstCascade::initialize ( )
overridevirtual

Definition at line 28 of file JpsiPlusDpstCascade.cxx.

28 {
29
30 // retrieving vertex Fitter
31 ATH_CHECK( m_iVertexFitter.retrieve());
32
33 // retrieving the V0 tools
34 ATH_CHECK( m_V0Tools.retrieve());
35
36 // retrieving the Cascade tools
37 ATH_CHECK( m_CascadeTools.retrieve());
38
39 // Get the beam spot service
40 ATH_CHECK(m_eventInfo_key.initialize());
41
42 ATH_CHECK( m_partPropSvc.retrieve() );
43 auto pdt = m_partPropSvc->PDT();
44
45 // retrieve particle masses
49
55
56 return StatusCode::SUCCESS;
57 }
static double getParticleMass(const HepPDT::ParticleDataTable *pdt, int pdg)
ServiceHandle< IPartPropSvc > m_partPropSvc
static const int KPLUS
static const int MUON
static const int BCPLUS
static const int PIPLUS
static const int JPSI
static const int D0

◆ performSearch()

StatusCode DerivationFramework::JpsiPlusDpstCascade::performSearch ( std::vector< Trk::VxCascadeInfo * > * cascadeinfoContainer,
const EventContext & ctx ) const

Definition at line 469 of file JpsiPlusDpstCascade.cxx.

470 {
471 ATH_MSG_DEBUG( "JpsiPlusDpstCascade::performSearch" );
472 assert(cascadeinfoContainer!=nullptr);
473
474 // Get TrackParticle container (for setting links to the original tracks)
475 const xAOD::TrackParticleContainer *trackContainer{};
476 ATH_CHECK(evtStore()->retrieve(trackContainer , "InDetTrackParticles" )); // FIXME Use Handles
477
478 // Get Jpsi+pi container
479 const xAOD::VertexContainer *jpsipiContainer{};
480 ATH_CHECK(evtStore()->retrieve(jpsipiContainer , m_vertexContainerKey )); // FIXME Use Handles
481
482 // Get D0 container
483 const xAOD::VertexContainer *d0Container{};
484 ATH_CHECK(evtStore()->retrieve(d0Container , m_vertexD0ContainerKey )); // FIXME Use Handles
485
486 double mass_d0 = m_vtx1MassHypo;
487 std::vector<const xAOD::TrackParticle*> tracksJpsipi;
488 std::vector<const xAOD::TrackParticle*> tracksJpsi;
489 std::vector<const xAOD::TrackParticle*> tracksD0;
490 std::vector<const xAOD::TrackParticle*> tracksBc;
491 std::vector<double> massesJpsipi;
492 massesJpsipi.push_back(m_vtx0Daug1MassHypo);
493 massesJpsipi.push_back(m_vtx0Daug2MassHypo);
494 massesJpsipi.push_back(m_vtx0Daug3MassHypo);
495 std::vector<double> massesD0;
496 massesD0.push_back(m_vtx1Daug1MassHypo);
497 massesD0.push_back(m_vtx1Daug2MassHypo);
498 std::vector<double> massesD0b; // Change the oreder of masses for D*-->D0bar pi-, D0bar->K+pi-
499 massesD0b.push_back(m_vtx1Daug2MassHypo);
500 massesD0b.push_back(m_vtx1Daug1MassHypo);
501 std::vector<double> Masses;
502 Masses.push_back(m_vtx0Daug1MassHypo);
503 Masses.push_back(m_vtx0Daug2MassHypo);
504 Masses.push_back(m_vtx0Daug3MassHypo);
505 Masses.push_back(m_vtx1MassHypo);
506
507 // Select J/psi pi+ candidates before calling cascade fit
508 std::vector<const xAOD::Vertex*> selectedJpsipiCandidates;
509 for(auto vxcItr=jpsipiContainer->cbegin(); vxcItr!=jpsipiContainer->cend(); ++vxcItr) {
510
511 // Check the passed flag first
512 const xAOD::Vertex* vtx = *vxcItr;
513 SG::AuxElement::Accessor<Char_t> flagAcc1("passed_Jpsipi");
514 if(flagAcc1.isAvailable(*vtx)){
515 if(!flagAcc1(*vtx)) continue;
516 }
517
518 // Check J/psi candidate invariant mass and skip if need be
519 TLorentzVector p4Mup_in, p4Mum_in;
520 p4Mup_in.SetPtEtaPhiM((*vxcItr)->trackParticle(0)->pt(),
521 (*vxcItr)->trackParticle(0)->eta(),
522 (*vxcItr)->trackParticle(0)->phi(), m_vtx0Daug1MassHypo);
523 p4Mum_in.SetPtEtaPhiM((*vxcItr)->trackParticle(1)->pt(),
524 (*vxcItr)->trackParticle(1)->eta(),
525 (*vxcItr)->trackParticle(1)->phi(), m_vtx0Daug2MassHypo);
526 double mass_Jpsi = (p4Mup_in + p4Mum_in).M();
527 ATH_MSG_DEBUG("Jpsi mass " << mass_Jpsi);
528 if (mass_Jpsi < m_jpsiMassLower || mass_Jpsi > m_jpsiMassUpper) {
529 ATH_MSG_DEBUG(" Original Jpsi candidate rejected by the mass cut: mass = "
530 << mass_Jpsi << " != (" << m_jpsiMassLower << ", " << m_jpsiMassUpper << ")" );
531 continue;
532 }
533
534 // Check J/psi pi+ candidate invariant mass and skip if need be
535 double mass_Jpsipi = m_V0Tools->invariantMass(*vxcItr, massesJpsipi);
536 ATH_MSG_DEBUG("Jpsipi mass " << mass_Jpsipi);
537 if (mass_Jpsipi < m_jpsipiMassLower || mass_Jpsipi > m_jpsipiMassUpper) {
538 ATH_MSG_DEBUG(" Original Jpsipi candidate rejected by the mass cut: mass = "
539 << mass_Jpsipi << " != (" << m_jpsipiMassLower << ", " << m_jpsipiMassUpper << ")" );
540 continue;
541 }
542
543 selectedJpsipiCandidates.push_back(*vxcItr);
544 }
545 if(selectedJpsipiCandidates.size()<1) return StatusCode::SUCCESS;
546
547 // Select the D0/D0b candidates before calling cascade fit
548 std::vector<const xAOD::Vertex*> selectedD0Candidates;
549 for(auto vxcItr=d0Container->cbegin(); vxcItr!=d0Container->cend(); ++vxcItr) {
550
551 // Check the passed flag first
552 const xAOD::Vertex* vtx = *vxcItr;
553 SG::AuxElement::Accessor<Char_t> flagAcc1("passed_D0");
554 SG::AuxElement::Accessor<Char_t> flagAcc2("passed_D0b");
555 bool isD0(true);
556 bool isD0b(true);
557 if(flagAcc1.isAvailable(*vtx)){
558 if(!flagAcc1(*vtx)) isD0 = false;
559 }
560 if(flagAcc2.isAvailable(*vtx)){
561 if(!flagAcc2(*vtx)) isD0b = false;
562 }
563 if(!(isD0||isD0b)) continue;
564
565 // Ensure the total charge is correct
566 if ((*vxcItr)->trackParticle(0)->charge() != 1 || (*vxcItr)->trackParticle(1)->charge() != -1) {
567 ATH_MSG_DEBUG(" Original D0/D0-bar candidate rejected by the charge requirement: "
568 << (*vxcItr)->trackParticle(0)->charge() << ", " << (*vxcItr)->trackParticle(1)->charge() );
569 continue;
570 }
571
572 // Check D0/D0bar candidate invariant mass and skip if need be
573 double mass_D0 = m_V0Tools->invariantMass(*vxcItr,massesD0);
574 double mass_D0b = m_V0Tools->invariantMass(*vxcItr,massesD0b);
575 ATH_MSG_DEBUG("D0 mass " << mass_D0 << ", D0b mass "<<mass_D0b);
576 if ((mass_D0 < m_D0MassLower || mass_D0 > m_D0MassUpper) && (mass_D0b < m_D0MassLower || mass_D0b > m_D0MassUpper)) {
577 ATH_MSG_DEBUG(" Original D0 candidate rejected by the mass cut: mass = "
578 << mass_D0 << " != (" << m_D0MassLower << ", " << m_D0MassUpper << ") "
579 << mass_D0b << " != (" << m_D0MassLower << ", " << m_D0MassUpper << ") " );
580 continue;
581 }
582
583 selectedD0Candidates.push_back(*vxcItr);
584 }
585 if(selectedD0Candidates.size()<1) return StatusCode::SUCCESS;
586
587 // Select J/psi D*+ candidates
588 // Iterate over Jpsi+pi vertices
589 for(auto jpsipiItr=selectedJpsipiCandidates.cbegin(); jpsipiItr!=selectedJpsipiCandidates.cend(); ++jpsipiItr) {
590
591 size_t jpsipiTrkNum = (*jpsipiItr)->nTrackParticles();
592 tracksJpsipi.clear();
593 tracksJpsi.clear();
594 for( unsigned int it=0; it<jpsipiTrkNum; it++) tracksJpsipi.push_back((*jpsipiItr)->trackParticle(it));
595 for( unsigned int it=0; it<jpsipiTrkNum-1; it++) tracksJpsi.push_back((*jpsipiItr)->trackParticle(it));
596
597 if (tracksJpsipi.size() != 3 || massesJpsipi.size() != 3 ) {
598 ATH_MSG_INFO("problems with Jpsi+pi input");
599 }
600
601 bool tagD0(true);
602 if(abs(m_Dx_pid)==421 && (*jpsipiItr)->trackParticle(2)->charge()==-1) tagD0 = false;
603
604 TLorentzVector p4_pi1; // Momentum of soft pion
605 p4_pi1.SetPtEtaPhiM((*jpsipiItr)->trackParticle(2)->pt(),
606 (*jpsipiItr)->trackParticle(2)->eta(),
607 (*jpsipiItr)->trackParticle(2)->phi(), m_vtx0Daug3MassHypo);
608
609 // Iterate over D0/D0bar vertices
610 for(auto d0Itr=selectedD0Candidates.cbegin(); d0Itr!=selectedD0Candidates.cend(); ++d0Itr) {
611
612 // Check identical tracks in input
613 if(std::find(tracksJpsipi.cbegin(), tracksJpsipi.cend(), (*d0Itr)->trackParticle(0)) != tracksJpsipi.cend()) continue;
614 if(std::find(tracksJpsipi.cbegin(), tracksJpsipi.cend(), (*d0Itr)->trackParticle(1)) != tracksJpsipi.cend()) continue;
615
616
617 TLorentzVector p4_ka, p4_pi2;
618 if(tagD0){ // for D*+
619 p4_pi2.SetPtEtaPhiM((*d0Itr)->trackParticle(0)->pt(),
620 (*d0Itr)->trackParticle(0)->eta(),
621 (*d0Itr)->trackParticle(0)->phi(), m_vtx1Daug1MassHypo);
622 p4_ka.SetPtEtaPhiM( (*d0Itr)->trackParticle(1)->pt(),
623 (*d0Itr)->trackParticle(1)->eta(),
624 (*d0Itr)->trackParticle(1)->phi(), m_vtx1Daug2MassHypo);
625 }else{ // change the order in the case of D*-
626 p4_pi2.SetPtEtaPhiM((*d0Itr)->trackParticle(1)->pt(),
627 (*d0Itr)->trackParticle(1)->eta(),
628 (*d0Itr)->trackParticle(1)->phi(), m_vtx1Daug1MassHypo);
629 p4_ka.SetPtEtaPhiM( (*d0Itr)->trackParticle(0)->pt(),
630 (*d0Itr)->trackParticle(0)->eta(),
631 (*d0Itr)->trackParticle(0)->phi(), m_vtx1Daug2MassHypo);
632 }
633 // Check D*+/- candidate invariant mass and skip if need be
634 double mass_Dst= (p4_pi1 + p4_ka + p4_pi2).M();
635 ATH_MSG_DEBUG("D*+/- mass " << mass_Dst);
636 if (mass_Dst < m_DstMassLower || mass_Dst > m_DstMassUpper) {
637 ATH_MSG_DEBUG(" Original D*+/- candidate rejected by the mass cut: mass = "
638 << mass_Dst << " != (" << m_DstMassLower << ", " << m_DstMassUpper << ")" );
639 continue;
640 }
641
642 size_t d0TrkNum = (*d0Itr)->nTrackParticles();
643 tracksD0.clear();
644 for( unsigned int it=0; it<d0TrkNum; it++) tracksD0.push_back((*d0Itr)->trackParticle(it));
645 if (tracksD0.size() != 2 || massesD0.size() != 2 ) {
646 ATH_MSG_INFO("problems with D0 input");
647 }
648
649 ATH_MSG_DEBUG("using tracks" << tracksJpsipi[0] << ", " << tracksJpsipi[1] << ", " << tracksJpsipi[2] << ", " << tracksD0[0] << ", " << tracksD0[1]);
650 ATH_MSG_DEBUG("Charge of Jpsi+pi tracks: "<<(*jpsipiItr)->trackParticle(0)->charge()<<", "<<(*jpsipiItr)->trackParticle(1)->charge()<<", "<<(*jpsipiItr)->trackParticle(2)->charge());
651 ATH_MSG_DEBUG("Charge of D0 tracks: "<<(*d0Itr)->trackParticle(0)->charge()<<", "<<(*d0Itr)->trackParticle(1)->charge());
652
653 tracksBc.clear();
654 for( unsigned int it=0; it<jpsipiTrkNum; it++) tracksBc.push_back((*jpsipiItr)->trackParticle(it));
655 for( unsigned int it=0; it<d0TrkNum; it++) tracksBc.push_back((*d0Itr)->trackParticle(it));
656
657
658 // Apply the user's settings to the fitter
659 // Reset
660 std::unique_ptr<Trk::IVKalState> state (m_iVertexFitter->makeState());
661 // Robustness
662 int robustness = 0;
663 m_iVertexFitter->setRobustness(robustness, *state);
664 // Build up the topology
665 // Vertex list
666 std::vector<Trk::VertexID> vrtList;
667 // D0 vertex
668 Trk::VertexID vID;
669 if (m_constrD0) {
670 if(tagD0) vID = m_iVertexFitter->startVertex(tracksD0,massesD0,*state,mass_d0);
671 else vID = m_iVertexFitter->startVertex(tracksD0,massesD0b,*state,mass_d0);
672 } else {
673 if(tagD0) vID = m_iVertexFitter->startVertex(tracksD0,massesD0,*state);
674 else vID = m_iVertexFitter->startVertex(tracksD0,massesD0b,*state);
675 }
676 vrtList.push_back(vID);
677 // B vertex including Jpsi+pi
678 Trk::VertexID vID2 = m_iVertexFitter->nextVertex(tracksJpsipi,massesJpsipi,vrtList,*state);
679 if (m_constrJpsi) {
680 std::vector<Trk::VertexID> cnstV;
681 cnstV.clear();
682 if ( !m_iVertexFitter->addMassConstraint(vID2,tracksJpsi,cnstV,*state,m_mass_jpsi).isSuccess() ) {
683 ATH_MSG_WARNING("addMassConstraint failed");
684 //return StatusCode::FAILURE;
685 }
686 }
687 // Do the work
688 std::unique_ptr<Trk::VxCascadeInfo> result(m_iVertexFitter->fitCascade(*state));
689
690 if (result != nullptr) {
691
692 // reset links to original tracks
693 BPhysPVCascadeTools::PrepareVertexLinks(result.get(), trackContainer);
694 ATH_MSG_DEBUG("storing tracks " << ((result->vertices())[0])->trackParticle(0) << ", "
695 << ((result->vertices())[0])->trackParticle(1) << ", "
696 << ((result->vertices())[1])->trackParticle(0) << ", "
697 << ((result->vertices())[1])->trackParticle(1) << ", "
698 << ((result->vertices())[1])->trackParticle(2));
699 // necessary to prevent memory leak
700 result->setSVOwnership(true);
701
702 // Chi2/DOF cut
703 double bChi2DOF = result->fitChi2()/result->nDoF();
704 ATH_MSG_DEBUG("Candidate chi2/DOF is " << bChi2DOF);
705 bool chi2CutPassed = (m_chi2cut <= 0.0 || bChi2DOF < m_chi2cut);
706
707 const std::vector< std::vector<TLorentzVector> > &moms = result->getParticleMoms();
708 double mass = m_CascadeTools->invariantMass(moms[1]);
709 if(chi2CutPassed) {
710 if (mass >= m_MassLower && mass <= m_MassUpper) {
711 cascadeinfoContainer->push_back(result.release());
712 } else {
713 ATH_MSG_DEBUG("Candidate rejected by the mass cut: mass = "
714 << mass << " != (" << m_MassLower << ", " << m_MassUpper << ")" );
715 }
716 }
717 }
718
719 } //Iterate over D0 vertices
720
721 } //Iterate over Jpsi+pi vertices
722
723 ATH_MSG_DEBUG("cascadeinfoContainer size " << cascadeinfoContainer->size());
724
725 return StatusCode::SUCCESS;
726 }
#define ATH_MSG_INFO(x)
const_iterator cbegin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
const_iterator cend() const noexcept
Return a const_iterator pointing past the end of the collection.
static void PrepareVertexLinks(Trk::VxCascadeInfo *result, const xAOD::TrackParticleContainer *importedTrackCollection)
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".

Member Data Documentation

◆ m_cascadeOutputsKeys

std::vector<std::string> DerivationFramework::JpsiPlusDpstCascade::m_cascadeOutputsKeys
private

Definition at line 48 of file JpsiPlusDpstCascade.h.

◆ m_CascadeTools

ToolHandle< DerivationFramework::CascadeTools > DerivationFramework::JpsiPlusDpstCascade::m_CascadeTools
private

Definition at line 81 of file JpsiPlusDpstCascade.h.

◆ m_chi2cut

double DerivationFramework::JpsiPlusDpstCascade::m_chi2cut
private

Definition at line 75 of file JpsiPlusDpstCascade.h.

◆ m_constrD0

bool DerivationFramework::JpsiPlusDpstCascade::m_constrD0
private

Definition at line 73 of file JpsiPlusDpstCascade.h.

◆ m_constrJpsi

bool DerivationFramework::JpsiPlusDpstCascade::m_constrJpsi
private

Definition at line 74 of file JpsiPlusDpstCascade.h.

◆ m_D0MassLower

double DerivationFramework::JpsiPlusDpstCascade::m_D0MassLower
private

Definition at line 56 of file JpsiPlusDpstCascade.h.

◆ m_D0MassUpper

double DerivationFramework::JpsiPlusDpstCascade::m_D0MassUpper
private

Definition at line 57 of file JpsiPlusDpstCascade.h.

◆ m_DoVertexType

int DerivationFramework::JpsiPlusDpstCascade::m_DoVertexType
private

Definition at line 90 of file JpsiPlusDpstCascade.h.

◆ m_DstMassLower

double DerivationFramework::JpsiPlusDpstCascade::m_DstMassLower
private

Definition at line 58 of file JpsiPlusDpstCascade.h.

◆ m_DstMassUpper

double DerivationFramework::JpsiPlusDpstCascade::m_DstMassUpper
private

Definition at line 59 of file JpsiPlusDpstCascade.h.

◆ m_Dx_pid

int DerivationFramework::JpsiPlusDpstCascade::m_Dx_pid
private

Definition at line 72 of file JpsiPlusDpstCascade.h.

◆ m_eventInfo_key

SG::ReadHandleKey<xAOD::EventInfo> DerivationFramework::JpsiPlusDpstCascade::m_eventInfo_key {this, "EventInfo", "EventInfo", "Input event information"}
private

Definition at line 77 of file JpsiPlusDpstCascade.h.

77{this, "EventInfo", "EventInfo", "Input event information"};

◆ m_hypoName

std::string DerivationFramework::JpsiPlusDpstCascade::m_hypoName
private

name of the mass hypothesis.

E.g. Jpsi, Upsi, etc. Will be used as a prefix for decorations

Definition at line 86 of file JpsiPlusDpstCascade.h.

◆ m_iVertexFitter

ToolHandle< Trk::TrkVKalVrtFitter > DerivationFramework::JpsiPlusDpstCascade::m_iVertexFitter
private

Definition at line 78 of file JpsiPlusDpstCascade.h.

◆ m_jpsiMassLower

double DerivationFramework::JpsiPlusDpstCascade::m_jpsiMassLower
private

Definition at line 52 of file JpsiPlusDpstCascade.h.

◆ m_jpsiMassUpper

double DerivationFramework::JpsiPlusDpstCascade::m_jpsiMassUpper
private

Definition at line 53 of file JpsiPlusDpstCascade.h.

◆ m_jpsipiMassLower

double DerivationFramework::JpsiPlusDpstCascade::m_jpsipiMassLower
private

Definition at line 54 of file JpsiPlusDpstCascade.h.

◆ m_jpsipiMassUpper

double DerivationFramework::JpsiPlusDpstCascade::m_jpsipiMassUpper
private

Definition at line 55 of file JpsiPlusDpstCascade.h.

◆ m_mass_jpsi

double DerivationFramework::JpsiPlusDpstCascade::m_mass_jpsi
private

Definition at line 71 of file JpsiPlusDpstCascade.h.

◆ m_MassLower

double DerivationFramework::JpsiPlusDpstCascade::m_MassLower
private

Definition at line 60 of file JpsiPlusDpstCascade.h.

◆ m_MassUpper

double DerivationFramework::JpsiPlusDpstCascade::m_MassUpper
private

Definition at line 61 of file JpsiPlusDpstCascade.h.

◆ m_partPropSvc

ServiceHandle<IPartPropSvc> DerivationFramework::JpsiPlusDpstCascade::m_partPropSvc {this, "PartPropSvc", "PartPropSvc"}
private

Definition at line 82 of file JpsiPlusDpstCascade.h.

82{this, "PartPropSvc", "PartPropSvc"};

◆ m_PV_max

int DerivationFramework::JpsiPlusDpstCascade::m_PV_max
private

Definition at line 89 of file JpsiPlusDpstCascade.h.

◆ m_PV_minNTracks

size_t DerivationFramework::JpsiPlusDpstCascade::m_PV_minNTracks
private

Definition at line 91 of file JpsiPlusDpstCascade.h.

◆ m_pvRefitter

ToolHandle< Analysis::PrimaryVertexRefitter > DerivationFramework::JpsiPlusDpstCascade::m_pvRefitter
private

Definition at line 79 of file JpsiPlusDpstCascade.h.

◆ m_refitPV

bool DerivationFramework::JpsiPlusDpstCascade::m_refitPV
private

Definition at line 84 of file JpsiPlusDpstCascade.h.

◆ m_refPVContainerName

std::string DerivationFramework::JpsiPlusDpstCascade::m_refPVContainerName
private

Definition at line 85 of file JpsiPlusDpstCascade.h.

◆ m_V0Tools

ToolHandle< Trk::V0Tools > DerivationFramework::JpsiPlusDpstCascade::m_V0Tools
private

Definition at line 80 of file JpsiPlusDpstCascade.h.

◆ m_vertexContainerKey

std::string DerivationFramework::JpsiPlusDpstCascade::m_vertexContainerKey
private

Definition at line 46 of file JpsiPlusDpstCascade.h.

◆ m_vertexD0ContainerKey

std::string DerivationFramework::JpsiPlusDpstCascade::m_vertexD0ContainerKey
private

Definition at line 47 of file JpsiPlusDpstCascade.h.

◆ m_vtx0Daug1MassHypo

double DerivationFramework::JpsiPlusDpstCascade::m_vtx0Daug1MassHypo
private

Definition at line 64 of file JpsiPlusDpstCascade.h.

◆ m_vtx0Daug2MassHypo

double DerivationFramework::JpsiPlusDpstCascade::m_vtx0Daug2MassHypo
private

Definition at line 65 of file JpsiPlusDpstCascade.h.

◆ m_vtx0Daug3MassHypo

double DerivationFramework::JpsiPlusDpstCascade::m_vtx0Daug3MassHypo
private

Definition at line 66 of file JpsiPlusDpstCascade.h.

◆ m_vtx0MassHypo

double DerivationFramework::JpsiPlusDpstCascade::m_vtx0MassHypo
private

Definition at line 62 of file JpsiPlusDpstCascade.h.

◆ m_vtx1Daug1MassHypo

double DerivationFramework::JpsiPlusDpstCascade::m_vtx1Daug1MassHypo
private

Definition at line 67 of file JpsiPlusDpstCascade.h.

◆ m_vtx1Daug2MassHypo

double DerivationFramework::JpsiPlusDpstCascade::m_vtx1Daug2MassHypo
private

Definition at line 68 of file JpsiPlusDpstCascade.h.

◆ m_vtx1MassHypo

double DerivationFramework::JpsiPlusDpstCascade::m_vtx1MassHypo
private

Definition at line 63 of file JpsiPlusDpstCascade.h.

◆ m_VxPrimaryCandidateName

std::string DerivationFramework::JpsiPlusDpstCascade::m_VxPrimaryCandidateName
private

Name of primary vertex container // FIXME Use Handles.

Definition at line 50 of file JpsiPlusDpstCascade.h.


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