ATLAS Offline Software
Loading...
Searching...
No Matches
Select_onia2mumu.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// Select_onia2mumu.cxx
8// Author: Daniel Scheirich <daniel.scheirich@cern.ch>
9// Based on the Integrated Simulation Framework
10//
11// Basic Jpsi->mu mu derivation example
12
13#include "Select_onia2mumu.h"
14
20#include <vector>
21#include <string>
22
23namespace DerivationFramework {
24
26 const std::string& n,
27 const IInterface* p) :
28 base_class(t,n,p),
29 m_v0Tools("Trk::V0Tools")
30 {
31
32 // Declare tools
33 declareProperty("V0Tools", m_v0Tools);
34
35 // Declare user-defined properties
36
37 declareProperty("HypothesisName" , m_hypoName = "A");
38 declareProperty("InputVtxContainerName", m_inputVtxContainerName = "JpsiCandidates");
39 declareProperty("TrkMasses" , m_trkMasses = std::vector<double>(2, ParticleConstants::muonMassInMeV) );
40 declareProperty("VtxMassHypo" , m_massHypo = ParticleConstants::JpsiMassInMeV );
41 declareProperty("MassMax" , m_massMax = 6000);
42 declareProperty("MassMin" , m_massMin = 2000);
43 declareProperty("Chi2Max" , m_chi2Max = 200);
44 declareProperty("DoVertexType" , m_DoVertexType = 7);
45 declareProperty("LxyMin" , m_lxyMin = std::numeric_limits<double>::lowest());
46 declareProperty("Do3d" , m_do3d = false);
47
48 }
49
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
51
53 {
54
55 ATH_MSG_DEBUG("in initialize()");
56
57 // retrieve V0 tools
58 CHECK( m_v0Tools.retrieve() );
60
61 return StatusCode::SUCCESS;
62
63 }
64
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
66
68 constexpr float errConst = -9999999;
69 const xAOD::Vertex* pv = onia.pv(pv_t);
70 if(pv) {
71 // decorate the vertex.
72 // Proper decay time assuming constant mass hypothesis m_massHypo
73 BPHYS_CHECK( onia.setTau( m_v0Tools->tau(onia.vtx(), pv, m_massHypo),
74 pv_t,
76 // Proper decay time assuming error constant mass hypothesis m_massHypo
77 BPHYS_CHECK( onia.setTauErr( m_v0Tools->tauError(onia.vtx(), pv, m_massHypo),
78 pv_t,
80
81 BPHYS_CHECK( onia.setTau( m_v0Tools->tau(onia.vtx(), pv, m_trkMasses),
82 pv_t,
84
85 BPHYS_CHECK( onia.setTauErr( m_v0Tools->tauError(onia.vtx(), pv, m_trkMasses),
86 pv_t,
88
89 //enum pv_type {PV_MAX_SUM_PT2, PV_MIN_A0, PV_MIN_Z0, PV_MIN_Z0_BA};
90 }else{
91
92
93 BPHYS_CHECK( onia.setTau(errConst, pv_t,
95 // Proper decay time assuming error constant mass hypothesis m_massHypo
96 BPHYS_CHECK( onia.setTauErr( errConst,
97 pv_t,
99
100 BPHYS_CHECK( onia.setTau( errConst,
101 pv_t,
103
104 BPHYS_CHECK( onia.setTauErr( errConst,
105 pv_t,
107 }
108
109 if(m_do3d){
110 BPHYS_CHECK( onia.setTau3d( pv ? m_v0Tools->tau3D(onia.vtx(), pv, m_massHypo) : errConst,
111 pv_t,
113 // Proper decay time assuming error constant mass hypothesis m_massHypo
114 BPHYS_CHECK( onia.setTau3dErr( pv ? m_v0Tools->tau3DError(onia.vtx(), pv, m_massHypo) : errConst,
115 pv_t,
117
118 BPHYS_CHECK( onia.setTau3d( pv ? m_v0Tools->tau3D(onia.vtx(), pv, m_trkMasses) : errConst,
119 pv_t,
121
122 BPHYS_CHECK( onia.setTau3dErr( pv ? m_v0Tools->tau3DError(onia.vtx(), pv, m_trkMasses) : errConst,
123 pv_t,
125
126 }
127
128 }
129
130
131 StatusCode Select_onia2mumu::addBranches(const EventContext& ctx) const
132 {
134 SG::auxid_set_t decor_auxids;
135
136 bool doPt = (m_DoVertexType & 1) != 0;
137 bool doA0 = (m_DoVertexType & 2) != 0;
138 bool doZ0 = (m_DoVertexType & 4) != 0;
139 bool doZ0BA = (m_DoVertexType & 8) != 0;
140 // loop over onia candidates and perform selection and augmentation
141 xAOD::VertexContainer::const_iterator oniaItr = oniaContainer->begin();
142 for(; oniaItr!=oniaContainer->end(); ++oniaItr) {
143 // create BPhysHypoHelper
144 xAOD::BPhysHypoHelper onia(m_hypoName, *oniaItr, &decor_auxids);
145 if((*oniaItr)->nTrackParticles() != m_trkMasses.size())
146 ATH_MSG_WARNING("Vertex has " << (*oniaItr)->nTrackParticles() << " while provided masses " << m_trkMasses.size());
147 //----------------------------------------------------
148 // decorate the vertex
149 //----------------------------------------------------
150 // a) invariant mass and error
151 if( !onia.setMass(m_trkMasses) ) ATH_MSG_WARNING("Decoration onia.setMass failed");
152
153 double massErr = m_v0Tools->invariantMassError(onia.vtx(), m_trkMasses);
154 if( !onia.setMassErr(massErr) ) ATH_MSG_WARNING("Decoration onia.setMassErr failed");
155
156 // b) proper decay time and error:
157 // retrieve the refitted PV (or the original one, if the PV refitting was turned off)
162
163 //----------------------------------------------------
164 // perform the selection (i.e. flag the vertex)
165 //----------------------------------------------------
166 // flag the vertex indicating that it is selected by this selector
167 onia.setPass(true);
168
169 // now we check othe cuts. if one of them didn't pass, set the flag to 0
170 // and continue to the next candidate:
171
172 // 1) invariant mass cut
173 if( onia.mass() < m_massMin || onia.mass() > m_massMax) {
174 onia.setPass(false); // flag as failed
175 continue;
176 }
177
178 // 2) chi2 cut
179 if( onia.vtx()->chiSquared() > m_chi2Max) {
180 onia.setPass(false);; // flag as failed
181 continue;
182 }
183 // 3) lxy cut
185 onia.setPass(false);; // flag as failed
186 continue;
187 }
188
189 } // end of loop over onia candidates
190
191 // Lock the decorations we just produced.
193 const_cast<xAOD::VertexContainer*> (oniaContainer.cptr());
194 for (SG::auxid_t auxid : decor_auxids) {
195 onia_nc->lockDecoration (auxid);
196 }
197
198 // all OK
199 return StatusCode::SUCCESS;
200 }
201
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
203
204
205}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define BPHYS_CHECK(EXP)
Useful CHECK macro.
: B-physics xAOD helpers.
#define CHECK(...)
Evaluate an expression and check for errors.
A number of constexpr particle constants to avoid hardcoding them directly in various places.
Define macros for attributes used to control the static checker.
#define ATLAS_THREAD_SAFE
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
SG::ReadHandleKey< xAOD::VertexContainer > m_inputVtxContainerName
name of the input container name
virtual StatusCode addBranches(const EventContext &ctx) const override
: augmentation and selection Retrieved vertices are augmented with usual information.
StatusCode initialize() override
inirialization and finalization
Select_onia2mumu(const std::string &t, const std::string &n, const IInterface *p)
double m_massMin
invariant mass range
int m_DoVertexType
Allows user to skip certain vertexes - bitwise test 7==all(111)
ToolHandle< Trk::V0Tools > m_v0Tools
tools
std::vector< double > m_trkMasses
track mass hypotheses
double m_massHypo
vertex mass hypothesis
double m_massMax
invariant mass range
void ProcessVertex(xAOD::BPhysHypoHelper &, xAOD::BPhysHelper::pv_type) const
const_pointer_type cptr()
Dereference the pointer.
A set of aux data identifiers.
Definition AuxTypes.h:47
float lxy(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Get the transverse decay distance and its error measured between the refitted primary vertex of type ...
const xAOD::Vertex * vtx() const
Getter method for the cached xAOD::Vertex.
const xAOD::Vertex * pv(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Get the refitted collision vertex of type pv_type.
pv_type
: Enum type of the PV
bool setTau3d(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0, const tau_type tauType=BPhysHypoHelper::TAU_CONST_MASS)
proper decay time
bool setTau3dErr(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0, const tau_type tauType=BPhysHypoHelper::TAU_CONST_MASS)
proper decay time error
float mass() const
Get invariant mass and its error.
bool setMass(const float val)
Set given invariant mass and its error.
bool setPass(bool passVal)
get the pass flag for this hypothesis
bool setTau(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0, const tau_type tauType=BPhysHypoHelper::TAU_CONST_MASS)
: Set the proper decay time and error.
bool setMassErr(const float val)
invariant mass error
bool setTauErr(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0, const tau_type tauType=BPhysHypoHelper::TAU_CONST_MASS)
proper decay time error
float chiSquared() const
Returns the of the vertex fit as float.
THE reconstruction tool.
constexpr double muonMassInMeV
the mass of the muon (in MeV)
size_t auxid_t
Identifier for a particular aux data item.
Definition AuxTypes.h:27
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.