ATLAS Offline Software
Loading...
Searching...
No Matches
DStarSelectionTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5//==================================================
6// Selection of D*+ -> pi+ + D0
7//==================================================
8
12#include "GaudiKernel/ThreadLocalContext.h"
13#include "Gaudi/Property.h"
14
17
18namespace{
19 static const SG::AuxElement::ConstAccessor<Char_t> flag_D0("passed_D0");
20 static const SG::AuxElement::ConstAccessor<Char_t> flag_D0b("passed_D0b");
21 static const SG::Accessor<std::vector<float>> acc_reFit_Px("RefTrackPx");
22 static const SG::Accessor<std::vector<float>> acc_reFit_Py("RefTrackPy");
23 static const SG::Accessor<std::vector<float>> acc_reFit_Pz("RefTrackPz");
24}
25
26
27namespace DerivationFramework {
28
30 ATH_MSG_DEBUG("in initialize()");
31 ATH_CHECK(m_trackKey.initialize());
33 ATH_CHECK(m_trackDecoKey.initialize());
34 ATH_CHECK(m_vertexDecoKey.initialize());
35 return StatusCode::SUCCESS;
36
37}
38
39StatusCode DStarSelectionTool::addBranches(const EventContext& ctx) const {
40
41 // Track container
42 const xAOD::TrackParticleContainer* trackParticleContainer{nullptr};
43 ATH_CHECK(SG::get(trackParticleContainer, m_trackKey, ctx));
44
45 // Vertex container
46 const xAOD::VertexContainer* D0Container{nullptr};
47 ATH_CHECK(SG::get(D0Container, m_inputVtxContainerName, ctx));
48
49 // Mask for track / vertex skimming
50 std::vector<char> track_pass_map(trackParticleContainer->size());
51 std::vector<char> vertex_pass_map(D0Container->size());
52
53 // Loop over D0 candidate vertices
54 for(const xAOD::Vertex* vertex: *D0Container) {
55
56 // Check vertex passes loose D0 cuts defined in python
57
58 bool passed_D0_or_D0b = (flag_D0(*vertex) || flag_D0b(*vertex));
59
60 // Only consider vertices which have passed D0 mass window for
61 // matching to soft pion canddiates
62 if(!passed_D0_or_D0b) continue;
63
64 const xAOD::TrackParticle* track1 = vertex->trackParticle(0);
65 const xAOD::TrackParticle* track2 = vertex->trackParticle(1);
66
67 if(!track1) { ATH_MSG_WARNING("Could not find track at D0 vertex (index 0)"); continue; }
68 if(!track2) { ATH_MSG_WARNING("Could not find track at D0 vertex (index 1)"); continue; }
69
70 // Re-fitted D0 track 4-vectors
71
72 const std::vector<float>& reFit_Px = acc_reFit_Px(*vertex);
73 const std::vector<float>& reFit_Py = acc_reFit_Py(*vertex);
74 const std::vector<float>& reFit_Pz = acc_reFit_Pz(*vertex);
75
76 TLorentzVector track1_pion, track1_kaon;
77 track1_pion.SetXYZM(reFit_Px.at(0),reFit_Py.at(0),reFit_Pz.at(0),m_pionMass);
78 track1_kaon.SetXYZM(reFit_Px.at(0),reFit_Py.at(0),reFit_Pz.at(0),m_kaonMass);
79
80 TLorentzVector track2_pion, track2_kaon;
81 track2_pion.SetXYZM(reFit_Px.at(1),reFit_Py.at(1),reFit_Pz.at(1),m_pionMass);
82 track2_kaon.SetXYZM(reFit_Px.at(1),reFit_Py.at(1),reFit_Pz.at(1),m_kaonMass);
83
84 // Consider both D0 and D0bar hypotheses (irrespective of charges)
85 TLorentzVector D0_hypo1 = track1_pion + track2_kaon;
86 TLorentzVector D0_hypo2 = track1_kaon + track2_pion;
87
88 bool passed_Dstar = false;
89
90 // Loop over tracks
91 for(const xAOD::TrackParticle* track3: *trackParticleContainer) {
92
93 // Avoid duplicates
94 if( track1->index() == track3->index() ) continue;
95 if( track2->index() == track3->index() ) continue;
96
97 TLorentzVector track3_pion= track3->p4();
98
99 const double deltaM_hypo1 = (D0_hypo1 + track3_pion).M() - D0_hypo1.M();
100 const double deltaM_hypo2 = (D0_hypo2 + track3_pion).M() - D0_hypo2.M();
101
102 // Selection of loose D*+ ->pi+ + D0 candidates based on mass difference
103 // Keep both "right" and "wrong" charge combinations
104 if( deltaM_hypo1 < m_deltaMassMax || deltaM_hypo2 < m_deltaMassMax ) {
105
106 passed_Dstar = true;
107
108 // Mark all tracks as being part of a D*+ -> D0 + pi+ (+c.c.) candidate
109 track_pass_map[track1->index()] = true;
110 track_pass_map[track2->index()] = true;
111 track_pass_map[track3->index()] = true;
112
113 }
114
115 } // Loop over tracks
116
117 // Adjust mask if we'd like to retain this vertex
118 if(passed_Dstar) {
119 vertex_pass_map[vertex->index()] = true;
120 }
121
122 } // Loop over D0 candidate vertices
123
124 // Decorate tracks with skim decision flag
125
127 // Loop over tracks
128 for(const xAOD::TrackParticle* track:*trackParticleContainer) {
129
130 flagTrackPass(*track) = track_pass_map[track->index()] ; // Mark track to keep/drop during skimming
131
132 } // Loop over tracks
133
134 // Decorate vertex for skimming decision (Char used for compatability for BPhys skim/slim tools)
136
137 // Loop over D0 candidate vertices
138 for(const xAOD::Vertex* vertex: *D0Container) {
139
140 // Mark vertex to keep during skimming, based on earlier decision
141 flagVertexPass(*vertex) = vertex_pass_map[vertex->index()];
142
143 } // Loop over D0 candidate vertices
144
145 return StatusCode::SUCCESS;
146}
147
148
149
150}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Handle class for reading from StoreGate.
Handle class for adding a decoration to an object.
Handle for requesting thinning for a data object.
size_type size() const noexcept
Returns the number of elements in the collection.
SG::WriteDecorHandleKey< xAOD::VertexContainer > m_vertexDecoKey
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trackKey
SG::WriteDecorHandleKey< xAOD::TrackParticleContainer > m_trackDecoKey
virtual StatusCode addBranches(const EventContext &ctx) const override
SG::ReadHandleKey< xAOD::VertexContainer > m_inputVtxContainerName
Helper class to provide type-safe access to aux data.
SG::ConstAccessor< T, ALLOC > ConstAccessor
Definition AuxElement.h:569
size_t index() const
Return the index of this element within its container.
Handle class for adding a decoration to an object.
THE reconstruction tool.
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".