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

#include <JpsiPlusPsiCascade.h>

Inheritance diagram for DerivationFramework::JpsiPlusPsiCascade:
Collaboration diagram for DerivationFramework::JpsiPlusPsiCascade:

Public Member Functions

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

Private Attributes

SG::ReadHandleKey< xAOD::VertexContainerm_vertexContainerKey
SG::ReadHandleKey< xAOD::VertexContainerm_vertexPsiContainerKey
std::vector< std::string > m_vertexJpsiHypoNames
std::vector< std::string > m_vertexPsiHypoNames
SG::WriteHandleKeyArray< xAOD::VertexContainerm_cascadeOutputsKeys
SG::ReadHandleKey< xAOD::VertexContainerm_VxPrimaryCandidateName
 Name of primary vertex container.
SG::ReadHandleKey< xAOD::TrackParticleContainerm_trackContainerName
SG::ReadHandleKey< xAOD::EventInfom_eventInfo_key
double m_jpsiMassLower
double m_jpsiMassUpper
double m_diTrackMassLower
double m_diTrackMassUpper
double m_psiMassLower
double m_psiMassUpper
double m_jpsi2MassLower
double m_jpsi2MassUpper
double m_MassLower
double m_MassUpper
int m_vtx1Daug_num
double m_vtx1Daug1MassHypo
double m_vtx1Daug2MassHypo
double m_vtx1Daug3MassHypo
double m_vtx1Daug4MassHypo
double m_vtx2Daug1MassHypo
double m_vtx2Daug2MassHypo
double m_mass_jpsi
double m_mass_diTrk
double m_mass_psi
double m_mass_jpsi2
bool m_constrPsi
bool m_constrJpsi
bool m_constrDiTrk
bool m_constrJpsi2
double m_chi2cut_Psi
double m_chi2cut_Jpsi
double m_chi2cut
unsigned int m_maxCandidates
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
SG::WriteHandleKey< xAOD::VertexContainerm_refPVContainerName
std::string m_hypoName
int m_PV_max
int m_DoVertexType
size_t m_PV_minNTracks

Detailed Description

Definition at line 33 of file JpsiPlusPsiCascade.h.

Constructor & Destructor Documentation

◆ JpsiPlusPsiCascade()

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

Definition at line 252 of file JpsiPlusPsiCascade.cxx.

252 : base_class(type,name,parent),
255 m_cascadeOutputsKeys({"JpsiPlusPsiCascadeVtx1", "JpsiPlusPsiCascadeVtx2", "JpsiPlusPsiCascadeVtx3"}),
256 m_VxPrimaryCandidateName("PrimaryVertices"),
257 m_trackContainerName("InDetTrackParticles"),
258 m_eventInfo_key("EventInfo"),
259 m_jpsiMassLower(0.0),
260 m_jpsiMassUpper(20000.0),
261 m_diTrackMassLower(-1.0),
262 m_diTrackMassUpper(-1.0),
263 m_psiMassLower(0.0),
264 m_psiMassUpper(25000.0),
265 m_jpsi2MassLower(0.0),
266 m_jpsi2MassUpper(20000.0),
267 m_MassLower(0.0),
268 m_MassUpper(31000.0),
276 m_mass_jpsi(-1),
277 m_mass_diTrk(-1),
278 m_mass_psi(-1),
279 m_mass_jpsi2(-1),
280 m_constrPsi(false),
281 m_constrJpsi(false),
282 m_constrDiTrk(false),
283 m_constrJpsi2(false),
284 m_chi2cut_Psi(-1.0),
285 m_chi2cut_Jpsi(-1.0),
286 m_chi2cut(-1.0),
288 m_iVertexFitter("Trk::TrkVKalVrtFitter"),
289 m_pvRefitter("Analysis::PrimaryVertexRefitter", this),
290 m_V0Tools("Trk::V0Tools"),
291 m_CascadeTools("DerivationFramework::CascadeTools")
292 {
293 declareProperty("JpsiVertices", m_vertexContainerKey);
294 declareProperty("PsiVertices", m_vertexPsiContainerKey);
295 declareProperty("JpsiVtxHypoNames", m_vertexJpsiHypoNames);
296 declareProperty("PsiVtxHypoNames", m_vertexPsiHypoNames);
297 declareProperty("VxPrimaryCandidateName", m_VxPrimaryCandidateName);
298 declareProperty("TrackContainerName", m_trackContainerName);
299 declareProperty("RefPVContainerName", m_refPVContainerName = "RefittedPrimaryVertices");
300 declareProperty("JpsiMassLowerCut", m_jpsiMassLower);
301 declareProperty("JpsiMassUpperCut", m_jpsiMassUpper);
302 declareProperty("DiTrackMassLower", m_diTrackMassLower); // only effective when m_vtx1Daug_num=4
303 declareProperty("DiTrackMassUpper", m_diTrackMassUpper); // only effective when m_vtx1Daug_num=4
304 declareProperty("PsiMassLowerCut", m_psiMassLower);
305 declareProperty("PsiMassUpperCut", m_psiMassUpper);
306 declareProperty("Jpsi2MassLowerCut", m_jpsi2MassLower);
307 declareProperty("Jpsi2MassUpperCut", m_jpsi2MassUpper);
308 declareProperty("MassLowerCut", m_MassLower);
309 declareProperty("MassUpperCut", m_MassUpper);
310 declareProperty("HypothesisName", m_hypoName = "TQ");
311 declareProperty("NumberOfPsiDaughters", m_vtx1Daug_num); // 3 or 4 only
312 declareProperty("Vtx1Daug1MassHypo", m_vtx1Daug1MassHypo);
313 declareProperty("Vtx1Daug2MassHypo", m_vtx1Daug2MassHypo);
314 declareProperty("Vtx1Daug3MassHypo", m_vtx1Daug3MassHypo);
315 declareProperty("Vtx1Daug4MassHypo", m_vtx1Daug4MassHypo);
316 declareProperty("Vtx2Daug1MassHypo", m_vtx2Daug1MassHypo);
317 declareProperty("Vtx2Daug2MassHypo", m_vtx2Daug2MassHypo);
318 declareProperty("JpsiMass", m_mass_jpsi);
319 declareProperty("DiTrackMass", m_mass_diTrk);
320 declareProperty("PsiMass", m_mass_psi);
321 declareProperty("Jpsi2Mass", m_mass_jpsi2);
322 declareProperty("ApplyJpsiMassConstraint", m_constrJpsi);
323 declareProperty("ApplyDiTrackMassConstraint", m_constrDiTrk); // only effective when m_vtx1Daug_num=4
324 declareProperty("ApplyPsiMassConstraint", m_constrPsi);
325 declareProperty("ApplyJpsi2MassConstraint", m_constrJpsi2);
326 declareProperty("Chi2CutPsi", m_chi2cut_Psi);
327 declareProperty("Chi2CutJpsi", m_chi2cut_Jpsi);
328 declareProperty("Chi2Cut", m_chi2cut);
329 declareProperty("MaxCandidates", m_maxCandidates);
330 declareProperty("RefitPV", m_refitPV = true);
331 declareProperty("MaxnPV", m_PV_max = 1000);
332 declareProperty("MinNTracksInPV", m_PV_minNTracks = 0);
333 declareProperty("DoVertexType", m_DoVertexType = 7);
334 declareProperty("TrkVertexFitterTool", m_iVertexFitter);
335 declareProperty("PVRefitter", m_pvRefitter);
336 declareProperty("V0Tools", m_V0Tools);
337 declareProperty("CascadeTools", m_CascadeTools);
338 declareProperty("CascadeVertexCollections", m_cascadeOutputsKeys);
339 }
SG::WriteHandleKeyArray< xAOD::VertexContainer > m_cascadeOutputsKeys
SG::WriteHandleKey< xAOD::VertexContainer > m_refPVContainerName
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexPsiContainerKey
std::vector< std::string > m_vertexJpsiHypoNames
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trackContainerName
SG::ReadHandleKey< xAOD::VertexContainer > m_VxPrimaryCandidateName
Name of primary vertex container.
ToolHandle< Analysis::PrimaryVertexRefitter > m_pvRefitter
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainerKey
ToolHandle< Trk::TrkVKalVrtFitter > m_iVertexFitter
ToolHandle< DerivationFramework::CascadeTools > m_CascadeTools
std::vector< std::string > m_vertexPsiHypoNames
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfo_key

◆ ~JpsiPlusPsiCascade()

virtual DerivationFramework::JpsiPlusPsiCascade::~JpsiPlusPsiCascade ( )
virtualdefault

Member Function Documentation

◆ addBranches()

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

Definition at line 62 of file JpsiPlusPsiCascade.cxx.

62 {
63 if (m_vtx1Daug_num != 3 && m_vtx1Daug_num != 4) {
64 ATH_MSG_FATAL("Incorrect number of Psi daughters (should be 3 or 4)");
65 return StatusCode::FAILURE;
66 }
67
68 constexpr int topoN = 3;
69 if(m_cascadeOutputsKeys.size() != topoN) {
70 ATH_MSG_FATAL("Incorrect number of VtxContainers");
71 return StatusCode::FAILURE;
72 }
73 std::array<SG::WriteHandle<xAOD::VertexContainer>, topoN> VtxWriteHandles; int ikey(0);
74 for(const SG::WriteHandleKey<xAOD::VertexContainer>& key : m_cascadeOutputsKeys) {
75 VtxWriteHandles[ikey] = SG::WriteHandle<xAOD::VertexContainer>(key, ctx);
76 ATH_CHECK( VtxWriteHandles[ikey].record(std::make_unique<xAOD::VertexContainer>(), std::make_unique<xAOD::VertexAuxContainer>()) );
77 ikey++;
78 }
79
80 //----------------------------------------------------
81 // retrieve primary vertices
82 //----------------------------------------------------
83 SG::ReadHandle<xAOD::VertexContainer> pvContainer(m_VxPrimaryCandidateName, ctx);
84 ATH_CHECK( pvContainer.isValid() );
85 if (pvContainer.cptr()->size()==0) {
86 ATH_MSG_WARNING("You have no primary vertices: " << pvContainer.cptr()->size());
87 return StatusCode::RECOVERABLE;
88 }
89
90 //----------------------------------------------------
91 // Record refitted primary vertices
92 //----------------------------------------------------
93 SG::WriteHandle<xAOD::VertexContainer> refPvContainer;
94 if(m_refitPV) {
95 refPvContainer = SG::WriteHandle<xAOD::VertexContainer>(m_refPVContainerName, ctx);
96 ATH_CHECK( refPvContainer.record(std::make_unique<xAOD::VertexContainer>(), std::make_unique<xAOD::VertexAuxContainer>()) );
97 }
98
99 std::vector<Trk::VxCascadeInfo*> cascadeinfoContainer;
100 std::vector<Trk::VxCascadeInfo*> cascadeinfoContainer_noConstr;
101 ATH_CHECK(performSearch(&cascadeinfoContainer,&cascadeinfoContainer_noConstr,ctx));
102
103 SG::ReadHandle<xAOD::EventInfo> evt(m_eventInfo_key, ctx);
104 ATH_CHECK( evt.isValid() );
105 BPhysPVCascadeTools helper(&(*m_CascadeTools), evt.cptr());
106 helper.SetMinNTracksInPV(m_PV_minNTracks);
107
108 // Decorators for the vertices
109 SG::AuxElement::Decorator<VertexLinkVector> CascadeLinksDecor("CascadeVertexLinks");
110 SG::AuxElement::Decorator<VertexLinkVector> PsiLinksDecor("PsiVertexLinks");
111 SG::AuxElement::Decorator<VertexLinkVector> JpsiLinksDecor("JpsiVertexLinks");
112 SG::AuxElement::Decorator<float> chi2_decor("ChiSquared");
113 SG::AuxElement::Decorator<int> ndof_decor("nDoF");
114 SG::AuxElement::Decorator<float> chi2_nc_decor("ChiSquared_nc");
115 SG::AuxElement::Decorator<int> ndof_nc_decor("nDoF_nc");
116 SG::AuxElement::Decorator<float> Pt_decor("Pt");
117 SG::AuxElement::Decorator<float> PtErr_decor("PtErr");
118 SG::AuxElement::Decorator<float> chi2_SV1_decor("ChiSquared_SV1");
119 SG::AuxElement::Decorator<float> chi2_nc_SV1_decor("ChiSquared_nc_SV1");
120 SG::AuxElement::Decorator<float> chi2_V1_decor("ChiSquared_V1");
121 SG::AuxElement::Decorator<int> ndof_V1_decor("nDoF_V1");
122 SG::AuxElement::Decorator<float> lxy_SV1_decor("lxy_SV1");
123 SG::AuxElement::Decorator<float> lxyErr_SV1_decor("lxyErr_SV1");
124 SG::AuxElement::Decorator<float> a0xy_SV1_decor("a0xy_SV1");
125 SG::AuxElement::Decorator<float> a0xyErr_SV1_decor("a0xyErr_SV1");
126 SG::AuxElement::Decorator<float> a0z_SV1_decor("a0z_SV1");
127 SG::AuxElement::Decorator<float> a0zErr_SV1_decor("a0zErr_SV1");
128 SG::AuxElement::Decorator<float> chi2_SV2_decor("ChiSquared_SV2");
129 SG::AuxElement::Decorator<float> chi2_nc_SV2_decor("ChiSquared_nc_SV2");
130 SG::AuxElement::Decorator<float> chi2_V2_decor("ChiSquared_V2");
131 SG::AuxElement::Decorator<int> ndof_V2_decor("nDoF_V2");
132 SG::AuxElement::Decorator<float> lxy_SV2_decor("lxy_SV2");
133 SG::AuxElement::Decorator<float> lxyErr_SV2_decor("lxyErr_SV2");
134 SG::AuxElement::Decorator<float> a0xy_SV2_decor("a0xy_SV2");
135 SG::AuxElement::Decorator<float> a0xyErr_SV2_decor("a0xyErr_SV2");
136 SG::AuxElement::Decorator<float> a0z_SV2_decor("a0z_SV2");
137 SG::AuxElement::Decorator<float> a0zErr_SV2_decor("a0zErr_SV2");
138
139 // Get the containers and identify the input Jpsi and Psi
140 SG::ReadHandle<xAOD::VertexContainer> psiContainer(m_vertexPsiContainerKey, ctx);
141 ATH_CHECK( psiContainer.isValid() );
142 SG::ReadHandle<xAOD::VertexContainer> jpsiContainer(m_vertexContainerKey, ctx);
143 ATH_CHECK( jpsiContainer.isValid() );
144
145 for(size_t ic=0; ic<cascadeinfoContainer.size(); ic++) {
146 Trk::VxCascadeInfo* cascade_info = cascadeinfoContainer[ic];
147 if(cascade_info==nullptr) {
148 ATH_MSG_ERROR("CascadeInfo is null");
149 continue;
150 }
151
152 Trk::VxCascadeInfo* cascade_info_noConstr = cascadeinfoContainer_noConstr[ic];
153
154 const std::vector<xAOD::Vertex*> &cascadeVertices = cascade_info->vertices();
155 if(cascadeVertices.size() != topoN) ATH_MSG_ERROR("Incorrect number of vertices");
156 if(cascadeVertices[0]==nullptr || cascadeVertices[1]==nullptr || cascadeVertices[2]==nullptr) ATH_MSG_ERROR("Error null vertex");
157 // Keep vertices
158 for(int i=0; i<topoN; i++) VtxWriteHandles[i].ptr()->push_back(cascadeVertices[i]);
159
160 cascade_info->setSVOwnership(false); // Prevent Container from deleting vertices
161 const auto mainVertex = cascadeVertices[2]; // this is the mother vertex
162 const std::vector< std::vector<TLorentzVector> > &moms = cascade_info->getParticleMoms();
163
164 // Set links to cascade vertices
165 std::vector<VertexLink> precedingVertexLinks;
166 VertexLink vertexLink1;
167 vertexLink1.setElement(cascadeVertices[0]);
168 vertexLink1.setStorableObject(*VtxWriteHandles[0].ptr());
169 if( vertexLink1.isValid() ) precedingVertexLinks.push_back( vertexLink1 );
170 VertexLink vertexLink2;
171 vertexLink2.setElement(cascadeVertices[1]);
172 vertexLink2.setStorableObject(*VtxWriteHandles[1].ptr());
173 if( vertexLink2.isValid() ) precedingVertexLinks.push_back( vertexLink2 );
174 CascadeLinksDecor(*mainVertex) = precedingVertexLinks;
175
176 // Identify the input Jpsi
177 const xAOD::Vertex* jpsiVertex = BPhysPVCascadeTools::FindVertex<2>(jpsiContainer.cptr(), cascadeVertices[1]);
178 // Identify the input Psi
179 const xAOD::Vertex* psiVertex(0);
180 if(m_vtx1Daug_num==4) psiVertex = BPhysPVCascadeTools::FindVertex<4>(psiContainer.cptr(), cascadeVertices[0]);
181 else psiVertex = BPhysPVCascadeTools::FindVertex<3>(psiContainer.cptr(), cascadeVertices[0]);
182
183 // Set links to input vertices
184 std::vector<const xAOD::Vertex*> jpsiVerticestoLink;
185 if(jpsiVertex) jpsiVerticestoLink.push_back(jpsiVertex);
186 else ATH_MSG_WARNING("Could not find linking Jpsi");
187 if(!BPhysPVCascadeTools::LinkVertices(JpsiLinksDecor, jpsiVerticestoLink, jpsiContainer.cptr(), mainVertex)) ATH_MSG_ERROR("Error decorating with Jpsi vertex");
188
189 std::vector<const xAOD::Vertex*> psiVerticestoLink;
190 if(psiVertex) psiVerticestoLink.push_back(psiVertex);
191 else ATH_MSG_WARNING("Could not find linking Psi");
192 if(!BPhysPVCascadeTools::LinkVertices(PsiLinksDecor, psiVerticestoLink, psiContainer.cptr(), mainVertex)) ATH_MSG_ERROR("Error decorating with Psi vertex");
193
194 xAOD::BPhysHypoHelper vtx(m_hypoName, mainVertex);
195
196 // Get refitted track momenta from all vertices, charged tracks only
197 BPhysPVCascadeTools::SetVectorInfo(vtx, cascade_info);
198 vtx.setPass(true);
199
200 //
201 // Decorate main vertex
202 //
203 // mass, mass error
204 // https://gitlab.cern.ch/atlas/athena/-/blob/21.2/Tracking/TrkVertexFitter/TrkVKalVrtFitter/TrkVKalVrtFitter/VxCascadeInfo.h
205 BPHYS_CHECK( vtx.setMass(m_CascadeTools->invariantMass(moms[2])) );
206 BPHYS_CHECK( vtx.setMassErr(m_CascadeTools->invariantMassError(moms[2],cascade_info->getCovariance()[2])) );
207 // pt and pT error (the default pt of mainVertex is != the pt of the full cascade fit!)
208 Pt_decor(*mainVertex) = m_CascadeTools->pT(moms[2]);
209 PtErr_decor(*mainVertex) = m_CascadeTools->pTError(moms[2],cascade_info->getCovariance()[2]);
210 // chi2 and ndof (the default chi2 of mainVertex is != the chi2 of the full cascade fit!)
211 chi2_decor(*mainVertex) = cascade_info->fitChi2();
212 ndof_decor(*mainVertex) = cascade_info->nDoF();
213 chi2_nc_decor(*mainVertex) = cascade_info_noConstr ? cascade_info_noConstr->fitChi2() : -999999.;
214 ndof_nc_decor(*mainVertex) = cascade_info_noConstr ? cascade_info_noConstr->nDoF() : -1;
215
216 // decorate the Psi vertex
217 chi2_SV1_decor(*cascadeVertices[0]) = m_V0Tools->chisq(cascadeVertices[0]);
218 chi2_nc_SV1_decor(*cascadeVertices[0]) = cascade_info_noConstr ? m_V0Tools->chisq(cascade_info_noConstr->vertices()[0]) : -999999.;
219 chi2_V1_decor(*cascadeVertices[0]) = m_V0Tools->chisq(psiVertex);
220 ndof_V1_decor(*cascadeVertices[0]) = m_V0Tools->ndof(psiVertex);
221 lxy_SV1_decor(*cascadeVertices[0]) = m_CascadeTools->lxy(moms[0],cascadeVertices[0],mainVertex);
222 lxyErr_SV1_decor(*cascadeVertices[0]) = m_CascadeTools->lxyError(moms[0],cascade_info->getCovariance()[0],cascadeVertices[0],mainVertex);
223 a0z_SV1_decor(*cascadeVertices[0]) = m_CascadeTools->a0z(moms[0],cascadeVertices[0],mainVertex);
224 a0zErr_SV1_decor(*cascadeVertices[0]) = m_CascadeTools->a0zError(moms[0],cascade_info->getCovariance()[0],cascadeVertices[0],mainVertex);
225 a0xy_SV1_decor(*cascadeVertices[0]) = m_CascadeTools->a0xy(moms[0],cascadeVertices[0],mainVertex);
226 a0xyErr_SV1_decor(*cascadeVertices[0]) = m_CascadeTools->a0xyError(moms[0],cascade_info->getCovariance()[0],cascadeVertices[0],mainVertex);
227
228 // decorate the Jpsi vertex
229 chi2_SV2_decor(*cascadeVertices[1]) = m_V0Tools->chisq(cascadeVertices[1]);
230 chi2_nc_SV2_decor(*cascadeVertices[1]) = cascade_info_noConstr ? m_V0Tools->chisq(cascade_info_noConstr->vertices()[1]) : -999999.;
231 chi2_V2_decor(*cascadeVertices[1]) = m_V0Tools->chisq(jpsiVertex);
232 ndof_V2_decor(*cascadeVertices[1]) = m_V0Tools->ndof(jpsiVertex);
233 lxy_SV2_decor(*cascadeVertices[1]) = m_CascadeTools->lxy(moms[1],cascadeVertices[1],mainVertex);
234 lxyErr_SV2_decor(*cascadeVertices[1]) = m_CascadeTools->lxyError(moms[1],cascade_info->getCovariance()[1],cascadeVertices[1],mainVertex);
235 a0z_SV2_decor(*cascadeVertices[1]) = m_CascadeTools->a0z(moms[1],cascadeVertices[1],mainVertex);
236 a0zErr_SV2_decor(*cascadeVertices[1]) = m_CascadeTools->a0zError(moms[1],cascade_info->getCovariance()[1],cascadeVertices[1],mainVertex);
237 a0xy_SV2_decor(*cascadeVertices[1]) = m_CascadeTools->a0xy(moms[1],cascadeVertices[1],mainVertex);
238 a0xyErr_SV2_decor(*cascadeVertices[1]) = m_CascadeTools->a0xyError(moms[1],cascade_info->getCovariance()[1],cascadeVertices[1],mainVertex);
239
240 double Mass_Moth = m_CascadeTools->invariantMass(moms[2]); //size=2
241 ATH_CHECK(helper.FillCandwithRefittedVertices(m_refitPV, pvContainer.cptr(), m_refitPV ? refPvContainer.ptr() : 0, &(*m_pvRefitter), m_PV_max, m_DoVertexType, cascade_info, 2, Mass_Moth, vtx));
242 } // loop over cascadeinfoContainer
243
244 // Deleting cascadeinfo since this won't be stored.
245 // Vertices have been kept in m_cascadeOutputs and should be owned by their container
246 for (auto cascade_info : cascadeinfoContainer) delete cascade_info;
247 for (auto cascade_info_noConstr : cascadeinfoContainer_noConstr) delete cascade_info_noConstr;
248
249 return StatusCode::SUCCESS;
250 }
#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 BPHYS_CHECK(EXP)
Useful CHECK macro.
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 *)
StatusCode performSearch(std::vector< Trk::VxCascadeInfo * > *cascadeinfoContainer, std::vector< Trk::VxCascadeInfo * > *cascadeinfoContainer_noConstr, const EventContext &ctx) const
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
pointer_type ptr()
Dereference the pointer.
double fitChi2() const
const std::vector< Amg::MatrixX > & getCovariance() const
const std::vector< std::vector< TLorentzVector > > & getParticleMoms() const
const std::vector< xAOD::Vertex * > & vertices() const
void setSVOwnership(bool Ownership)
ElementLink< xAOD::VertexContainer > VertexLink
void * ptr(T *p)
Definition SGImplSvc.cxx:74
int ic
Definition grepfile.py:33
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ initialize()

StatusCode DerivationFramework::JpsiPlusPsiCascade::initialize ( )
overridevirtual

Definition at line 22 of file JpsiPlusPsiCascade.cxx.

22 {
23 // retrieving vertex Fitter
24 ATH_CHECK( m_iVertexFitter.retrieve() );
25
26 // retrieve PV refitter
27 ATH_CHECK( m_pvRefitter.retrieve() );
28
29 // retrieving the V0 tools
30 ATH_CHECK( m_V0Tools.retrieve() );
31
32 // retrieving the Cascade tools
33 ATH_CHECK( m_CascadeTools.retrieve() );
34
35 ATH_CHECK( m_vertexContainerKey.initialize() );
36 ATH_CHECK( m_vertexPsiContainerKey.initialize() );
37 ATH_CHECK( m_VxPrimaryCandidateName.initialize() );
38 ATH_CHECK( m_trackContainerName.initialize() );
39 ATH_CHECK( m_refPVContainerName.initialize() );
40 ATH_CHECK( m_cascadeOutputsKeys.initialize() );
41 ATH_CHECK( m_eventInfo_key.initialize() );
42
43 ATH_CHECK( m_partPropSvc.retrieve() );
44 auto pdt = m_partPropSvc->PDT();
45
46 // retrieve particle masses
47 // https://gitlab.cern.ch/atlas/athena/-/blob/main/Generators/TruthUtils/TruthUtils/AtlasPID.h
51
58
59 return StatusCode::SUCCESS;
60 }
static double getParticleMass(const HepPDT::ParticleDataTable *pdt, int pdg)
ServiceHandle< IPartPropSvc > m_partPropSvc
static const int PSI2S
static const int MUON
static const int PIPLUS
static const int JPSI

◆ performSearch()

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

Definition at line 341 of file JpsiPlusPsiCascade.cxx.

341 {
342 ATH_MSG_DEBUG( "JpsiPlusPsiCascade::performSearch" );
343 assert(cascadeinfoContainer!=nullptr && cascadeinfoContainer_noConstr!=nullptr);
344
345 // Get TrackParticle container (for setting links to the original tracks)
346 SG::ReadHandle<xAOD::TrackParticleContainer> trackContainer(m_trackContainerName, ctx);
347 ATH_CHECK( trackContainer.isValid() );
348
349 std::vector<const xAOD::TrackParticle*> tracksJpsi;
350 std::vector<const xAOD::TrackParticle*> tracksDiTrk;
351 std::vector<const xAOD::TrackParticle*> tracksPsi;
352 std::vector<const xAOD::TrackParticle*> tracksJpsi2;
353 std::vector<double> massesPsi;
354 massesPsi.push_back(m_vtx1Daug1MassHypo);
355 massesPsi.push_back(m_vtx1Daug2MassHypo);
356 massesPsi.push_back(m_vtx1Daug3MassHypo);
357 if(m_vtx1Daug_num==4) massesPsi.push_back(m_vtx1Daug4MassHypo);
358 std::array<double,2> massesJpsi2{m_vtx2Daug1MassHypo, m_vtx2Daug2MassHypo};
359
360 // Get Psi container
361 SG::ReadHandle<xAOD::VertexContainer> psiContainer(m_vertexPsiContainerKey, ctx);
362 ATH_CHECK( psiContainer.isValid() );
363
364 // Get Jpsi container
365 SG::ReadHandle<xAOD::VertexContainer> jpsiContainer(m_vertexContainerKey, ctx);
366 ATH_CHECK( jpsiContainer.isValid() );
367
368 // Select the J/psi candidates before calling cascade fit
369 std::vector<const xAOD::Vertex*> selectedJpsiCandidates;
370 for(auto vxcItr=jpsiContainer.cptr()->cbegin(); vxcItr!=jpsiContainer.cptr()->cend(); ++vxcItr) {
371 // Check the passed flag first
372 const xAOD::Vertex* vtx = *vxcItr;
373 bool passed = false;
374 for(size_t i=0; i<m_vertexJpsiHypoNames.size(); i++) {
376 if(flagAcc.isAvailable(*vtx) && flagAcc(*vtx)) {
377 passed |= 1;
378 }
379 }
380 if(m_vertexJpsiHypoNames.size() && !passed) continue;
381
382 // Check Jpsi candidate invariant mass and skip if need be
383 double mass_jpsi2 = m_V0Tools->invariantMass(*vxcItr, massesJpsi2);
384 if (mass_jpsi2 < m_jpsi2MassLower || mass_jpsi2 > m_jpsi2MassUpper) continue;
385
386 double chi2DOF = (*vxcItr)->chiSquared()/(*vxcItr)->numberDoF();
387 if(m_chi2cut_Jpsi>0 && chi2DOF>m_chi2cut_Jpsi) continue;
388
389 selectedJpsiCandidates.push_back(*vxcItr);
390 }
391 if(selectedJpsiCandidates.size()==0) return StatusCode::SUCCESS;
392
393 // Select the Psi candidates before calling cascade fit
394 std::vector<const xAOD::Vertex*> selectedPsiCandidates;
395 for(auto vxcItr=psiContainer.cptr()->cbegin(); vxcItr!=psiContainer.cptr()->cend(); ++vxcItr) {
396 // Check the passed flag first
397 const xAOD::Vertex* vtx = *vxcItr;
398 bool passed = false;
399 for(size_t i=0; i<m_vertexPsiHypoNames.size(); i++) {
401 if(flagAcc.isAvailable(*vtx) && flagAcc(*vtx)) {
402 passed |= 1;
403 }
404 }
405 if(m_vertexPsiHypoNames.size() && !passed) continue;
406
407 // Check Psi candidate invariant mass and skip if need be
408 double mass_psi = m_V0Tools->invariantMass(*vxcItr,massesPsi);
409 if(mass_psi < m_psiMassLower || mass_psi > m_psiMassUpper) continue;
410
411 // Add loose cut on Jpsi mass from Psi -> Jpsi pi+ pi-, or on phi mass from Ds+ -> phi pi+
412 TLorentzVector p4_mu1, p4_mu2;
413 p4_mu1.SetPtEtaPhiM( (*vxcItr)->trackParticle(0)->pt(),
414 (*vxcItr)->trackParticle(0)->eta(),
415 (*vxcItr)->trackParticle(0)->phi(), m_vtx1Daug1MassHypo);
416 p4_mu2.SetPtEtaPhiM( (*vxcItr)->trackParticle(1)->pt(),
417 (*vxcItr)->trackParticle(1)->eta(),
418 (*vxcItr)->trackParticle(1)->phi(), m_vtx1Daug2MassHypo);
419 double mass_jpsi = (p4_mu1 + p4_mu2).M();
420 if (mass_jpsi < m_jpsiMassLower || mass_jpsi > m_jpsiMassUpper) continue;
421
423 TLorentzVector p4_trk1, p4_trk2;
424 p4_trk1.SetPtEtaPhiM( (*vxcItr)->trackParticle(2)->pt(),
425 (*vxcItr)->trackParticle(2)->eta(),
426 (*vxcItr)->trackParticle(2)->phi(), m_vtx1Daug3MassHypo);
427 p4_trk2.SetPtEtaPhiM( (*vxcItr)->trackParticle(3)->pt(),
428 (*vxcItr)->trackParticle(3)->eta(),
429 (*vxcItr)->trackParticle(3)->phi(), m_vtx1Daug4MassHypo);
430 double mass_diTrk = (p4_trk1 + p4_trk2).M();
431 if (mass_diTrk < m_diTrackMassLower || mass_diTrk > m_diTrackMassUpper) continue;
432 }
433
434 double chi2DOF = (*vxcItr)->chiSquared()/(*vxcItr)->numberDoF();
435 if(m_chi2cut_Psi>0 && chi2DOF>m_chi2cut_Psi) continue;
436
437 selectedPsiCandidates.push_back(*vxcItr);
438 }
439 if(selectedPsiCandidates.size()==0) return StatusCode::SUCCESS;
440
441 std::vector<std::pair<const xAOD::Vertex*, const xAOD::Vertex*> > candidatePairs;
442 for(auto jpsiItr=selectedJpsiCandidates.cbegin(); jpsiItr!=selectedJpsiCandidates.cend(); ++jpsiItr) {
443 tracksJpsi2.clear();
444 for(size_t i=0; i<(*jpsiItr)->nTrackParticles(); i++) tracksJpsi2.push_back((*jpsiItr)->trackParticle(i));
445 for(auto psiItr=selectedPsiCandidates.cbegin(); psiItr!=selectedPsiCandidates.cend(); ++psiItr) {
446 bool skip = false;
447 for(size_t j=0; j<(*psiItr)->nTrackParticles(); j++) {
448 if(std::find(tracksJpsi2.cbegin(), tracksJpsi2.cend(), (*psiItr)->trackParticle(j)) != tracksJpsi2.cend()) { skip = true; break; }
449 }
450 if(skip) continue;
451 candidatePairs.push_back(std::pair<const xAOD::Vertex*, const xAOD::Vertex*>(*jpsiItr,*psiItr));
452 }
453 }
454
455 std::sort( candidatePairs.begin(), candidatePairs.end(), [](std::pair<const xAOD::Vertex*, const xAOD::Vertex*> a, std::pair<const xAOD::Vertex*, const xAOD::Vertex*> b) { return a.first->chiSquared()/a.first->numberDoF()+a.second->chiSquared()/a.second->numberDoF() < b.first->chiSquared()/b.first->numberDoF()+b.second->chiSquared()/b.second->numberDoF(); } );
456 if(m_maxCandidates>0 && candidatePairs.size()>m_maxCandidates) {
457 candidatePairs.erase(candidatePairs.begin()+m_maxCandidates, candidatePairs.end());
458 }
459
460 for(size_t ic=0; ic<candidatePairs.size(); ic++) {
461 const xAOD::Vertex* jpsiVertex = candidatePairs[ic].first;
462 const xAOD::Vertex* psiVertex = candidatePairs[ic].second;
463
464 tracksJpsi2.clear();
465 for(size_t it=0; it<jpsiVertex->nTrackParticles(); it++) tracksJpsi2.push_back(jpsiVertex->trackParticle(it));
466 if (tracksJpsi2.size() != 2 || massesJpsi2.size() != 2) {
467 ATH_MSG_ERROR("Problems with Jpsi input: number of tracks or track mass inputs is not 2!");
468 }
469 tracksPsi.clear();
470 for(size_t it=0; it<psiVertex->nTrackParticles(); it++) tracksPsi.push_back(psiVertex->trackParticle(it));
471 if (tracksPsi.size() != massesPsi.size()) {
472 ATH_MSG_ERROR("Problems with Psi input: number of tracks or track mass inputs is not correct!");
473 }
474
475 tracksJpsi.clear();
476 tracksJpsi.push_back(psiVertex->trackParticle(0));
477 tracksJpsi.push_back(psiVertex->trackParticle(1));
478 tracksDiTrk.clear();
479 if(m_vtx1Daug_num==4) {
480 tracksDiTrk.push_back(psiVertex->trackParticle(2));
481 tracksDiTrk.push_back(psiVertex->trackParticle(3));
482 }
483
484 TLorentzVector p4_moth;
485 TLorentzVector tmp;
486 for(size_t it=0; it<jpsiVertex->nTrackParticles(); it++) {
487 tmp.SetPtEtaPhiM(jpsiVertex->trackParticle(it)->pt(),jpsiVertex->trackParticle(it)->eta(),jpsiVertex->trackParticle(it)->phi(),massesJpsi2[it]);
488 p4_moth += tmp;
489 }
490 for(size_t it=0; it<psiVertex->nTrackParticles(); it++) {
491 tmp.SetPtEtaPhiM(psiVertex->trackParticle(it)->pt(),psiVertex->trackParticle(it)->eta(),psiVertex->trackParticle(it)->phi(),massesPsi[it]);
492 p4_moth += tmp;
493 }
494 if (p4_moth.M() < m_MassLower || p4_moth.M() > m_MassUpper) continue;
495
496 // Apply the user's settings to the fitter
497 std::unique_ptr<Trk::IVKalState> state = m_iVertexFitter->makeState();
498 // Robustness: http://cdsweb.cern.ch/record/685551
499 int robustness = 0;
500 m_iVertexFitter->setRobustness(robustness, *state);
501 // Build up the topology
502 // Vertex list
503 std::vector<Trk::VertexID> vrtList;
504 // Psi vertex
505 Trk::VertexID vID1;
506 // https://gitlab.cern.ch/atlas/athena/-/blob/21.2/Tracking/TrkVertexFitter/TrkVKalVrtFitter/TrkVKalVrtFitter/IVertexCascadeFitter.h
507 if (m_constrPsi) {
508 vID1 = m_iVertexFitter->startVertex(tracksPsi,massesPsi,*state,m_mass_psi);
509 } else {
510 vID1 = m_iVertexFitter->startVertex(tracksPsi,massesPsi,*state);
511 }
512 vrtList.push_back(vID1);
513 // Jpsi vertex
514 Trk::VertexID vID2;
515 if (m_constrJpsi2) {
516 vID2 = m_iVertexFitter->nextVertex(tracksJpsi2,massesJpsi2,*state,m_mass_jpsi2);
517 } else {
518 vID2 = m_iVertexFitter->nextVertex(tracksJpsi2,massesJpsi2,*state);
519 }
520 vrtList.push_back(vID2);
521 // Mother vertex including Jpsi and Psi
522 std::vector<const xAOD::TrackParticle*> tp; tp.clear();
523 std::vector<double> tp_masses; tp_masses.clear();
524 m_iVertexFitter->nextVertex(tp,tp_masses,vrtList,*state);
525 if (m_constrJpsi) {
526 std::vector<Trk::VertexID> cnstV; cnstV.clear();
527 if ( !m_iVertexFitter->addMassConstraint(vID1,tracksJpsi,cnstV,*state,m_mass_jpsi).isSuccess() ) {
528 ATH_MSG_WARNING("addMassConstraint for Jpsi failed");
529 }
530 }
532 std::vector<Trk::VertexID> cnstV; cnstV.clear();
533 if ( !m_iVertexFitter->addMassConstraint(vID1,tracksDiTrk,cnstV,*state,m_mass_diTrk).isSuccess() ) {
534 ATH_MSG_WARNING("addMassConstraint for DiTrk failed");
535 }
536 }
537 // Do the work
538 std::unique_ptr<Trk::VxCascadeInfo> result(m_iVertexFitter->fitCascade(*state));
539
540 bool pass = false;
541 if (result != nullptr) {
542 for(auto v : result->vertices()) {
543 if(v->nTrackParticles()==0) {
544 std::vector<ElementLink<xAOD::TrackParticleContainer> > nullLinkVector;
545 v->setTrackParticleLinks(nullLinkVector);
546 }
547 }
548 // reset links to original tracks
549 BPhysPVCascadeTools::PrepareVertexLinks(result.get(), trackContainer.cptr());
550
551 // necessary to prevent memory leak
552 result->setSVOwnership(true);
553
554 // Chi2/DOF cut
555 double chi2DOF = result->fitChi2()/result->nDoF();
556 bool chi2CutPassed = (m_chi2cut <= 0.0 || chi2DOF < m_chi2cut);
557
558 if(chi2CutPassed) {
559 cascadeinfoContainer->push_back(result.release());
560 pass = true;
561 }
562 }
563
564 // do cascade fit again without any mass constraints
565 if(pass) {
567 std::unique_ptr<Trk::IVKalState> state (m_iVertexFitter->makeState());
568 m_iVertexFitter->setRobustness(robustness, *state);
569 std::vector<Trk::VertexID> vrtList_nc;
570 // Psi vertex
571 Trk::VertexID vID1_nc = m_iVertexFitter->startVertex(tracksPsi,massesPsi,*state);
572 vrtList_nc.push_back(vID1_nc);
573 Trk::VertexID vID2_nc = m_iVertexFitter->nextVertex(tracksJpsi2,massesJpsi2,*state);
574 vrtList_nc.push_back(vID2_nc);
575 // Mother vertex including Jpsi and Psi
576 std::vector<const xAOD::TrackParticle*> tp; tp.clear();
577 std::vector<double> tp_masses; tp_masses.clear();
578 m_iVertexFitter->nextVertex(tp,tp_masses,vrtList_nc,*state);
579 // Do the work
580 std::unique_ptr<Trk::VxCascadeInfo> result_nc(m_iVertexFitter->fitCascade(*state));
581
582 if (result_nc != nullptr) {
583 for(auto v : result_nc->vertices()) {
584 if(v->nTrackParticles()==0) {
585 std::vector<ElementLink<xAOD::TrackParticleContainer> > nullLinkVector;
586 v->setTrackParticleLinks(nullLinkVector);
587 }
588 }
589 // reset links to original tracks
590 BPhysPVCascadeTools::PrepareVertexLinks(result_nc.get(), trackContainer.cptr());
591
592 // necessary to prevent memory leak
593 result_nc->setSVOwnership(true);
594 cascadeinfoContainer_noConstr->push_back(result_nc.release());
595 }
596 else cascadeinfoContainer_noConstr->push_back(0);
597 }
598 else cascadeinfoContainer_noConstr->push_back(0);
599 }
600 } //Iterate over candidatePairs
601
602 return StatusCode::SUCCESS;
603 }
#define ATH_MSG_DEBUG(x)
bool passed(DecisionID id, const DecisionIDContainer &)
checks if required decision ID is in the set of IDs in the container
static Double_t a
static void PrepareVertexLinks(Trk::VxCascadeInfo *result, const xAOD::TrackParticleContainer *importedTrackCollection)
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
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.
size_t nTrackParticles() const
Get the number of tracks associated with this vertex.
const TrackParticle * trackParticle(size_t i) const
Get the pointer to a given track that was used in vertex reco.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.

Member Data Documentation

◆ m_cascadeOutputsKeys

SG::WriteHandleKeyArray<xAOD::VertexContainer> DerivationFramework::JpsiPlusPsiCascade::m_cascadeOutputsKeys
private

Definition at line 47 of file JpsiPlusPsiCascade.h.

◆ m_CascadeTools

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

Definition at line 86 of file JpsiPlusPsiCascade.h.

◆ m_chi2cut

double DerivationFramework::JpsiPlusPsiCascade::m_chi2cut
private

Definition at line 80 of file JpsiPlusPsiCascade.h.

◆ m_chi2cut_Jpsi

double DerivationFramework::JpsiPlusPsiCascade::m_chi2cut_Jpsi
private

Definition at line 79 of file JpsiPlusPsiCascade.h.

◆ m_chi2cut_Psi

double DerivationFramework::JpsiPlusPsiCascade::m_chi2cut_Psi
private

Definition at line 78 of file JpsiPlusPsiCascade.h.

◆ m_constrDiTrk

bool DerivationFramework::JpsiPlusPsiCascade::m_constrDiTrk
private

Definition at line 76 of file JpsiPlusPsiCascade.h.

◆ m_constrJpsi

bool DerivationFramework::JpsiPlusPsiCascade::m_constrJpsi
private

Definition at line 75 of file JpsiPlusPsiCascade.h.

◆ m_constrJpsi2

bool DerivationFramework::JpsiPlusPsiCascade::m_constrJpsi2
private

Definition at line 77 of file JpsiPlusPsiCascade.h.

◆ m_constrPsi

bool DerivationFramework::JpsiPlusPsiCascade::m_constrPsi
private

Definition at line 74 of file JpsiPlusPsiCascade.h.

◆ m_diTrackMassLower

double DerivationFramework::JpsiPlusPsiCascade::m_diTrackMassLower
private

Definition at line 54 of file JpsiPlusPsiCascade.h.

◆ m_diTrackMassUpper

double DerivationFramework::JpsiPlusPsiCascade::m_diTrackMassUpper
private

Definition at line 55 of file JpsiPlusPsiCascade.h.

◆ m_DoVertexType

int DerivationFramework::JpsiPlusPsiCascade::m_DoVertexType
private

Definition at line 93 of file JpsiPlusPsiCascade.h.

◆ m_eventInfo_key

SG::ReadHandleKey<xAOD::EventInfo> DerivationFramework::JpsiPlusPsiCascade::m_eventInfo_key
private

Definition at line 50 of file JpsiPlusPsiCascade.h.

◆ m_hypoName

std::string DerivationFramework::JpsiPlusPsiCascade::m_hypoName
private

Definition at line 91 of file JpsiPlusPsiCascade.h.

◆ m_iVertexFitter

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

Definition at line 83 of file JpsiPlusPsiCascade.h.

◆ m_jpsi2MassLower

double DerivationFramework::JpsiPlusPsiCascade::m_jpsi2MassLower
private

Definition at line 58 of file JpsiPlusPsiCascade.h.

◆ m_jpsi2MassUpper

double DerivationFramework::JpsiPlusPsiCascade::m_jpsi2MassUpper
private

Definition at line 59 of file JpsiPlusPsiCascade.h.

◆ m_jpsiMassLower

double DerivationFramework::JpsiPlusPsiCascade::m_jpsiMassLower
private

Definition at line 52 of file JpsiPlusPsiCascade.h.

◆ m_jpsiMassUpper

double DerivationFramework::JpsiPlusPsiCascade::m_jpsiMassUpper
private

Definition at line 53 of file JpsiPlusPsiCascade.h.

◆ m_mass_diTrk

double DerivationFramework::JpsiPlusPsiCascade::m_mass_diTrk
private

Definition at line 71 of file JpsiPlusPsiCascade.h.

◆ m_mass_jpsi

double DerivationFramework::JpsiPlusPsiCascade::m_mass_jpsi
private

Definition at line 70 of file JpsiPlusPsiCascade.h.

◆ m_mass_jpsi2

double DerivationFramework::JpsiPlusPsiCascade::m_mass_jpsi2
private

Definition at line 73 of file JpsiPlusPsiCascade.h.

◆ m_mass_psi

double DerivationFramework::JpsiPlusPsiCascade::m_mass_psi
private

Definition at line 72 of file JpsiPlusPsiCascade.h.

◆ m_MassLower

double DerivationFramework::JpsiPlusPsiCascade::m_MassLower
private

Definition at line 60 of file JpsiPlusPsiCascade.h.

◆ m_MassUpper

double DerivationFramework::JpsiPlusPsiCascade::m_MassUpper
private

Definition at line 61 of file JpsiPlusPsiCascade.h.

◆ m_maxCandidates

unsigned int DerivationFramework::JpsiPlusPsiCascade::m_maxCandidates
private

Definition at line 81 of file JpsiPlusPsiCascade.h.

◆ m_partPropSvc

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

Definition at line 87 of file JpsiPlusPsiCascade.h.

87{this, "PartPropSvc", "PartPropSvc"};

◆ m_psiMassLower

double DerivationFramework::JpsiPlusPsiCascade::m_psiMassLower
private

Definition at line 56 of file JpsiPlusPsiCascade.h.

◆ m_psiMassUpper

double DerivationFramework::JpsiPlusPsiCascade::m_psiMassUpper
private

Definition at line 57 of file JpsiPlusPsiCascade.h.

◆ m_PV_max

int DerivationFramework::JpsiPlusPsiCascade::m_PV_max
private

Definition at line 92 of file JpsiPlusPsiCascade.h.

◆ m_PV_minNTracks

size_t DerivationFramework::JpsiPlusPsiCascade::m_PV_minNTracks
private

Definition at line 94 of file JpsiPlusPsiCascade.h.

◆ m_pvRefitter

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

Definition at line 84 of file JpsiPlusPsiCascade.h.

◆ m_refitPV

bool DerivationFramework::JpsiPlusPsiCascade::m_refitPV
private

Definition at line 89 of file JpsiPlusPsiCascade.h.

◆ m_refPVContainerName

SG::WriteHandleKey<xAOD::VertexContainer> DerivationFramework::JpsiPlusPsiCascade::m_refPVContainerName
private

Definition at line 90 of file JpsiPlusPsiCascade.h.

◆ m_trackContainerName

SG::ReadHandleKey<xAOD::TrackParticleContainer> DerivationFramework::JpsiPlusPsiCascade::m_trackContainerName
private

Definition at line 49 of file JpsiPlusPsiCascade.h.

◆ m_V0Tools

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

Definition at line 85 of file JpsiPlusPsiCascade.h.

◆ m_vertexContainerKey

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::JpsiPlusPsiCascade::m_vertexContainerKey
private

Definition at line 43 of file JpsiPlusPsiCascade.h.

◆ m_vertexJpsiHypoNames

std::vector<std::string> DerivationFramework::JpsiPlusPsiCascade::m_vertexJpsiHypoNames
private

Definition at line 45 of file JpsiPlusPsiCascade.h.

◆ m_vertexPsiContainerKey

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::JpsiPlusPsiCascade::m_vertexPsiContainerKey
private

Definition at line 44 of file JpsiPlusPsiCascade.h.

◆ m_vertexPsiHypoNames

std::vector<std::string> DerivationFramework::JpsiPlusPsiCascade::m_vertexPsiHypoNames
private

Definition at line 46 of file JpsiPlusPsiCascade.h.

◆ m_vtx1Daug1MassHypo

double DerivationFramework::JpsiPlusPsiCascade::m_vtx1Daug1MassHypo
private

Definition at line 63 of file JpsiPlusPsiCascade.h.

◆ m_vtx1Daug2MassHypo

double DerivationFramework::JpsiPlusPsiCascade::m_vtx1Daug2MassHypo
private

Definition at line 64 of file JpsiPlusPsiCascade.h.

◆ m_vtx1Daug3MassHypo

double DerivationFramework::JpsiPlusPsiCascade::m_vtx1Daug3MassHypo
private

Definition at line 65 of file JpsiPlusPsiCascade.h.

◆ m_vtx1Daug4MassHypo

double DerivationFramework::JpsiPlusPsiCascade::m_vtx1Daug4MassHypo
private

Definition at line 66 of file JpsiPlusPsiCascade.h.

◆ m_vtx1Daug_num

int DerivationFramework::JpsiPlusPsiCascade::m_vtx1Daug_num
private

Definition at line 62 of file JpsiPlusPsiCascade.h.

◆ m_vtx2Daug1MassHypo

double DerivationFramework::JpsiPlusPsiCascade::m_vtx2Daug1MassHypo
private

Definition at line 67 of file JpsiPlusPsiCascade.h.

◆ m_vtx2Daug2MassHypo

double DerivationFramework::JpsiPlusPsiCascade::m_vtx2Daug2MassHypo
private

Definition at line 68 of file JpsiPlusPsiCascade.h.

◆ m_VxPrimaryCandidateName

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::JpsiPlusPsiCascade::m_VxPrimaryCandidateName
private

Name of primary vertex container.

Definition at line 48 of file JpsiPlusPsiCascade.h.


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