ATLAS Offline Software
Loading...
Searching...
No Matches
xAOD::BPhysHelper Class Reference

#include <BPhysHelper.h>

Inheritance diagram for xAOD::BPhysHelper:
Collaboration diagram for xAOD::BPhysHelper:

Public Types

enum  pv_type { PV_MAX_SUM_PT2 , PV_MIN_A0 , PV_MIN_Z0 , PV_MIN_Z0_BA }
 : Enum type of the PV More...

Public Member Functions

 BPhysHelper (const xAOD::Vertex *b)
 : Main constructor
const xAOD::Vertexvtx () const
 Getter method for the cached xAOD::Vertex.
const TMatrixTSym< double > & covariance ()
 : Returns full covariance matrix
int nRefTrks ()
 Returns number of stored refitted track momenta.
TVector3 refTrk (const size_t index)
 Returns i-th refitted track 3-momentum.
const std::vector< TVector3 > & refTrks ()
 Returns refitted track momenta.
TLorentzVector refTrk (const size_t index, const float mass)
 Returns i-th refitted track as a 4-vector given the mass hypothesis.
const xAOD::IParticlerefTrkOrigin (const size_t index) const
 : Returns the original track (charged or neutral) corresponding to the i-th refitted track
TVector3 refTrkOriginP (const size_t index) const
 Returns perigee 3-momentum of the original track corresponding i-th refitted track.
TLorentzVector refTrkOriginP (const size_t index, const float mass) const
 : Returns lorentz vector build from the original track perigee momentum and given mass hypothesis
float refTrkCharge (const size_t index) const
 Returns charge of the i-th track.
bool setRefTrks (std::vector< float > px, std::vector< float > py, std::vector< float > pz)
 Sets refitted track momenta.
bool setRefTrks (const std::vector< TVector3 > &refTrks)
 Sets refitted track momenta.
bool setRefTrks ()
 : Sets refitted track momenta
TVector3 totalP ()
 : Returns total 3-momentum calculated from the refitted tracks
TLorentzVector totalP (std::span< const double > masses)
 Returns total 4-momentum calculated from the refitted tracks given mass hypotheses for individual tracks.
float ptErr ()
 Returns pT error.
bool setPtErr (const float val)
 Set pT error.
int nMuons ()
 : Methods providing access to the linked muons
int nElectrons ()
const xAOD::Muonmuon (const size_t index)
 Returns pointer to the i-th linked muon.
const xAOD::Electronelectron (const size_t index)
const std::vector< const xAOD::Muon * > & muons ()
 Returns linked muons.
const std::vector< const xAOD::Electron * > & electrons ()
bool setMuons (const std::vector< const xAOD::Muon * > &muons, const xAOD::MuonContainer *muonContainer)
 Set links to muons.
bool setElectrons (const std::vector< const xAOD::Electron * > &electrons, const xAOD::ElectronContainer *electronContainer)
int nPrecedingVertices ()
 : Links to preceding vertices
const xAOD::VertexprecedingVertex (const size_t index)
 Returns pointer to a preceding vertex.
const std::vector< const xAOD::Vertex * > & precedingVertices ()
 Returns vector of pointers to preceding vertices.
bool setPrecedingVertices (const std::vector< const xAOD::Vertex * > &vertices, const xAOD::VertexContainer *vertexContainer)
 Sets links to preceding vertices.
int nCascadeVertices ()
 : Links to cascade vertices
const xAOD::VertexcascadeVertex (const size_t index)
 Returns pointer to a cascade vertex.
const std::vector< const xAOD::Vertex * > & cascadeVertices ()
 Returns vector of pointers to cascade vertices.
bool setCascadeVertices (const std::vector< const xAOD::Vertex * > &vertices, const xAOD::VertexContainer *vertexContainer)
 Sets links to cascade vertices.
int nRefTrksCascade ()
 Returns number of stored refitted tracks INCLUDING those from the linked cascade vertices.
const xAOD::Vertexpv (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Get the refitted collision vertex of type pv_type.
const xAOD::VertexorigPv (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 original PV
bool setOrigPv (const xAOD::Vertex *pv, const xAOD::VertexContainer *vertexContainer, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Set the original collision vertex of type pv_type.
bool setPv (const xAOD::Vertex *pv, const xAOD::VertexContainer *vertexContainer, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Set the refitted collision vertex of type pv_type.
bool setRefitPVStatus (int code, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Set the exitCode of the refitter for vertex of type pv_type.
int RefitPVStatus (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Get the exitCode of the refitter for vertex of type pv_type.
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 pv_type and the B vertex.
float lxyErr (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 its error
bool setLxy (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Set the transverse decay distance and its error measured between the refitted primary vertex of type pv_type and the B vertex.
bool setLxyErr (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 its error
float lxyz (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 decay distance
float lxyzErr (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 its error
bool setLxyz (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 decay distance
bool setLxyzErr (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 its error
float a0 (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Get the 3D, transverse, and longitudinal impact parameters and their error.
float a0Err (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 3D impact parameter error
float a0xy (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 transverse impact parameter
float a0xyErr (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 rtansverse impact parameter error
float z0 (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 longitudinal impact parameter
float z0Err (const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 longitudinal impact parameter error
float setA0 (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 Set the 3D and transverse impact parameters and their error.
float setA0Err (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 3D impact parameter error
float setA0xy (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 transverse impact parameter
float setA0xyErr (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 transverse impact parameter error
float setZ0 (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 longitudinal impact parameter
float setZ0Err (const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
 longitudinal impact parameter error

Static Public Attributes

static const std::string pv_type_str []
static const unsigned int n_pv_types = 4

Protected Member Functions

bool cacheCov ()
 : Cache covariance matrix
bool cacheRefTracks ()
 : Cache refitted tracks
bool cacheMuons ()
 : Cache linked muons
bool cacheElectrons ()
bool cachePrecedingVertices ()
 : Cache preceding vertices
bool cacheCascadeVertices ()
 : Cache cascade vertices

Protected Attributes

const xAOD::Vertexm_b
 Cached B decay xAOD vertex.
bool m_covCached
TMatrixTSym< doublem_cachedCov
bool m_refTracksCached
std::vector< TVector3 > m_cachedRefTracks
bool m_muonsCached
bool m_electronsCached
std::vector< const xAOD::Muon * > m_cachedMuons
std::vector< const xAOD::Electron * > m_cachedElectrons
bool m_precedingVerticesCached
std::vector< const xAOD::Vertex * > m_cachedPrecedingVertices
bool m_cascadeVerticesCached
std::vector< const xAOD::Vertex * > m_cachedCascadeVertices

Static Protected Attributes

static const std::vector< TVector3 > s_emptyVectorOfTVector3
static const std::vector< const xAOD::Muon * > s_emptyVectorOfMuons
static const std::vector< const xAOD::Electron * > s_emptyVectorOfElectrons
static const TMatrixTSym< doubles_emptyMatrix
static const std::vector< const xAOD::Vertex * > s_emptyVectorOfVertices

Detailed Description

Definition at line 71 of file BPhysHelper.h.

Member Enumeration Documentation

◆ pv_type

: Enum type of the PV

This enum type is used to identify type of the PV and quantities calculated using the PV. For instance, the transverse distance can be calculated wrt the default PV, the closest PV etc.

@value: PV_MAX_SUM_PT2 vertex with the largest sum Pt^2 (the default one) @value: PV_MIN_A0 collision vertex closest in 3D to the particle's trajectory, i.e. the vertex with the smalles 3D impact parameter a0 @value: PV_MIN_Z0 collision vertex closest in delta(z0) = pv_z - z0, where z0 is the z coordinate of the intesection of the particle's trajectory with the beam axis in the Z-rho plane. @value: PV_MIN_Z0_BA collision vertex closest in delta(z0) = pv_z - z0, where z0 is the z coordinate of the point of closest aproach (in 2D) to the particle's trajectory with the beam axis in the Z-rho plane.

Enumerator
PV_MAX_SUM_PT2 
PV_MIN_A0 
PV_MIN_Z0 
PV_MIN_Z0_BA 

Definition at line 475 of file BPhysHelper.h.

Constructor & Destructor Documentation

◆ BPhysHelper()

xAOD::BPhysHelper::BPhysHelper ( const xAOD::Vertex * b)
inline

: Main constructor

The Main constructor. BPhysHelper does NOT take ownership of the class "b"

Parameters
[in]bPointer to the xAOD::Vertex

Definition at line 85 of file BPhysHelper.h.

85 :
86 m_b(b),
87 m_covCached(false),
88 m_cachedCov(0),
89 m_refTracksCached(false),
91 m_muonsCached(false),
92 m_electronsCached(false),
99 {
100 assert(m_b!=0); // sanity check: m_b must not be null
101 }
std::vector< const xAOD::Vertex * > m_cachedCascadeVertices
const xAOD::Vertex * m_b
Cached B decay xAOD vertex.
TMatrixTSym< double > m_cachedCov
std::vector< const xAOD::Electron * > m_cachedElectrons
bool m_precedingVerticesCached
std::vector< TVector3 > m_cachedRefTracks
std::vector< const xAOD::Vertex * > m_cachedPrecedingVertices
std::vector< const xAOD::Muon * > m_cachedMuons

Member Function Documentation

◆ a0()

float xAOD::BPhysHelper::a0 ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

Get the 3D, transverse, and longitudinal impact parameters and their error.

Impact parameters and their error

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Returns
: impact parameter or its error, -9999999. in case augmentation is not available 3D impact parameter

Definition at line 957 of file BPhysHelper.cxx.

958{
959 switch(vertexType) {
960 case PV_MAX_SUM_PT2 : GET_FLOAT("A0MaxSumPt2");
961 case PV_MIN_A0 : GET_FLOAT("A0MinA0");
962 case PV_MIN_Z0 : GET_FLOAT("A0MinZ0");
963 case PV_MIN_Z0_BA : GET_FLOAT("A0MinZ0BA");
964 default: return -9999999.;
965 }
966}
#define GET_FLOAT(name)

◆ a0Err()

float xAOD::BPhysHelper::a0Err ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

3D impact parameter error

Definition at line 968 of file BPhysHelper.cxx.

969{
970 switch(vertexType) {
971 case PV_MAX_SUM_PT2 : GET_FLOAT("A0ErrMaxSumPt2");
972 case PV_MIN_A0 : GET_FLOAT("A0ErrMinA0");
973 case PV_MIN_Z0 : GET_FLOAT("A0ErrMinZ0");
974 case PV_MIN_Z0_BA : GET_FLOAT("A0ErrMinZ0BA");
975 default: return -9999999.;
976 }
977}

◆ a0xy()

float xAOD::BPhysHelper::a0xy ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

transverse impact parameter

Definition at line 979 of file BPhysHelper.cxx.

980{
981 switch(vertexType) {
982 case PV_MAX_SUM_PT2 : GET_FLOAT("A0xyMaxSumPt2");
983 case PV_MIN_A0 : GET_FLOAT("A0xyMinA0");
984 case PV_MIN_Z0 : GET_FLOAT("A0xyMinZ0");
985 case PV_MIN_Z0_BA : GET_FLOAT("A0xyMinZ0BA");
986 default: return -9999999.;
987 }
988}

◆ a0xyErr()

float xAOD::BPhysHelper::a0xyErr ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

rtansverse impact parameter error

Definition at line 990 of file BPhysHelper.cxx.

991{
992 switch(vertexType) {
993 case PV_MAX_SUM_PT2 : GET_FLOAT("A0xyErrMaxSumPt2");
994 case PV_MIN_A0 : GET_FLOAT("A0xyErrMinA0");
995 case PV_MIN_Z0 : GET_FLOAT("A0xyErrMinZ0");
996 case PV_MIN_Z0_BA : GET_FLOAT("A0xyErrMinZ0BA");
997 default: return -9999999.;
998 }
999}

◆ cacheCascadeVertices()

bool xAOD::BPhysHelper::cacheCascadeVertices ( )
protected

: Cache cascade vertices

To speed up access to linked cascade vertices, m_cachedCascadeVertices vector is created first time linked cascade vertices accessor methods are called. In subsequent calls, cached vector is used.

Returns
: true on success

Definition at line 1316 of file BPhysHelper.cxx.

1317{
1318 // if linked cascade vertices are already cached do nothing and return true
1320 return true;
1321
1322 // if linked cascade vertices are cached but the cache is empty, return false (error)
1324 return false;
1325
1326 // Linked cascade vertices are not cached, yet. Retrieve them from the aux store:
1329
1330 // Create auxiliary branches accessors
1331 static const SG::AuxElement::Accessor<VertexLinkVector> cascadeVertexLinksAcc("CascadeVertexLinks");
1332
1333 // check if branch exists
1334 if(!cascadeVertexLinksAcc.isAvailable(*m_b)) {
1335 return false;
1336 }
1337
1338 // retrieve the cascadeVertex links...
1339 const VertexLinkVector& cascadeVertexLinks = cascadeVertexLinksAcc(*m_b);
1340
1341 // ... and check if they are all valid
1342 VertexLinkVector::const_iterator cascadeVertexLinksItr = cascadeVertexLinks.begin();
1343 for(; cascadeVertexLinksItr!=cascadeVertexLinks.end(); ++cascadeVertexLinksItr) {
1344 // check if links are valid
1345 if(!(*cascadeVertexLinksItr).isValid()) {
1346 return false;
1347 }
1348 }
1349
1350 // all OK, cache the cascadeVertices
1351 cascadeVertexLinksItr = cascadeVertexLinks.begin();
1352 for(; cascadeVertexLinksItr!=cascadeVertexLinks.end(); ++cascadeVertexLinksItr) {
1353 m_cachedCascadeVertices.push_back(*(*cascadeVertexLinksItr));
1354 }
1355
1356 return true;
1357
1358}
std::vector< VertexLink > VertexLinkVector
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572

◆ cacheCov()

bool xAOD::BPhysHelper::cacheCov ( )
protected

: Cache covariance matrix

To speed up access to covariance matrix, m_cachedCov attribute is created first time covariance matrix accessor method is called. In subsequent calls, cached matrix is used.

Returns
: true on success

Definition at line 1089 of file BPhysHelper.cxx.

1090{
1091 // if cov matrix is already cached, do nothing and return true
1092 if(m_covCached && m_cachedCov.GetNrows()!=0)
1093 return true;
1094
1095 // if cov matrix is cached but the cache is empty, return false (error)
1096 if(m_covCached && m_cachedCov.GetNrows()==0)
1097 return false;
1098
1099 // get number of tracks
1100 // if ref tracks are stripped from the vertex, this will not work and
1101 // error will be returned
1102 const int nTrk = nRefTrks();
1103 if(nTrk<0)
1104 return false;
1105
1106 // Covariance matrix not cached, yet. Retrieve them from
1107 // the vertex and convert:
1108 m_covCached = true;
1109
1110 // retrieve native covariance matrix
1111 const std::vector<float>& cov = m_b->covariance();
1112
1113 // convert covariance matrix
1114 if(int(cov.size())==(3*nTrk+3)*(3*nTrk+3+1)/2) {
1115 // this is a VKalVrtFitter covariance matrix
1116 m_cachedCov.ResizeTo(3+3*nTrk, 3+3*nTrk);
1117
1118 long int ij=0;
1119 for(int i=0; i<(3+3*nTrk); i++) {
1120 for(int j=0; j<=i; j++) {
1121 m_cachedCov(i,j) = cov[ij];
1122 ij++;
1123 }
1124 }
1125
1126 return true;
1127 }
1128
1129 // FIXME: so far we can only conver VKalVrt cov matrix
1130
1131 // error
1132 return false;
1133
1134}
int nRefTrks()
Returns number of stored refitted track momenta.

◆ cacheElectrons()

bool xAOD::BPhysHelper::cacheElectrons ( )
protected

Definition at line 1226 of file BPhysHelper.cxx.

1227{
1228 // if linked electrons are already cached do nothing and return true
1229 if(m_electronsCached && !m_cachedElectrons.empty())
1230 return true;
1231
1232 // if linked electrons are cached but the cache is empty, return false (error)
1234 return false;
1235
1236 // Linked electrons are not cached, yet. Retrieve them from the aux store:
1237 m_electronsCached = true;
1238 m_cachedElectrons.clear();
1239
1240 // Create auxiliary branches accessors
1241 static const SG::AuxElement::Accessor<ElectronLinkVector> electronLinksAcc("ElectronLinks");
1242
1243 // check if branch exists
1244 if(!electronLinksAcc.isAvailable(*m_b)) {
1245 return false;
1246 }
1247
1248 // retrieve the electron links...
1249 const ElectronLinkVector& electronLinks = electronLinksAcc(*m_b);
1250
1251 // ... and check if they are all valid
1252 ElectronLinkVector::const_iterator electronLinksItr = electronLinks.begin();
1253 for(; electronLinksItr!=electronLinks.end(); ++electronLinksItr) {
1254 // check if links are valid
1255 if(!(*electronLinksItr).isValid()) {
1256 return false;
1257 }
1258 }
1259
1260 // all OK, cache the electrons
1261 electronLinksItr = electronLinks.begin();
1262 for(; electronLinksItr!=electronLinks.end(); ++electronLinksItr) {
1263 m_cachedElectrons.push_back(*(*electronLinksItr));
1264 }
1265
1266 return true;
1267
1268}
std::vector< ElectronLink > ElectronLinkVector

◆ cacheMuons()

bool xAOD::BPhysHelper::cacheMuons ( )
protected

: Cache linked muons

To speed up access to linked muons, m_cachedMuons vector is created first time linked muon accessor methods are called. In subsequent calls, cached vector is used.

Returns
: true on success

Definition at line 1181 of file BPhysHelper.cxx.

1182{
1183 // if linked muons are already cached do nothing and return true
1184 if(m_muonsCached && !m_cachedMuons.empty())
1185 return true;
1186
1187 // if linked muons are cached but the cache is empty, return false (error)
1188 if(m_muonsCached && m_cachedMuons.empty())
1189 return false;
1190
1191 // Linked muons are not cached, yet. Retrieve them from the aux store:
1192 m_muonsCached = true;
1193 m_cachedMuons.clear();
1194
1195 // Create auxiliary branches accessors
1196 static const SG::AuxElement::Accessor<MuonLinkVector> muonLinksAcc("MuonLinks");
1197
1198 // check if branch exists
1199 if(!muonLinksAcc.isAvailable(*m_b)) {
1200 return false;
1201 }
1202
1203 // retrieve the muon links...
1204 const MuonLinkVector& muonLinks = muonLinksAcc(*m_b);
1205
1206 // ... and check if they are all valid
1207 MuonLinkVector::const_iterator muonLinksItr = muonLinks.begin();
1208 for(; muonLinksItr!=muonLinks.end(); ++muonLinksItr) {
1209 // check if links are valid
1210 if(!(*muonLinksItr).isValid()) {
1211 return false;
1212 }
1213 }
1214
1215 // all OK, cache the muons
1216 muonLinksItr = muonLinks.begin();
1217 for(; muonLinksItr!=muonLinks.end(); ++muonLinksItr) {
1218 m_cachedMuons.push_back(*(*muonLinksItr));
1219 }
1220
1221 return true;
1222
1223}
std::vector< MuonLink > MuonLinkVector

◆ cachePrecedingVertices()

bool xAOD::BPhysHelper::cachePrecedingVertices ( )
protected

: Cache preceding vertices

To speed up access to linked preceding vertices, m_cachedPrecedingVertices vector is created first time linked preceding vertices accessor methods are called. In subsequent calls, cached vector is used.

Returns
: true on success

Definition at line 1271 of file BPhysHelper.cxx.

1272{
1273 // if linked preceding vertices are already cached do nothing and return true
1275 return true;
1276
1277 // if linked preceding vertices are cached but the cache is empty, return false (error)
1279 return false;
1280
1281 // Linked preceding vertices are not cached, yet. Retrieve them from the aux store:
1284
1285 // Create auxiliary branches accessors
1286 static const SG::AuxElement::Accessor<VertexLinkVector> precedingVertexLinksAcc("PrecedingVertexLinks");
1287
1288 // check if branch exists
1289 if(!precedingVertexLinksAcc.isAvailable(*m_b)) {
1290 return false;
1291 }
1292
1293 // retrieve the precedingVertex links...
1294 const VertexLinkVector& precedingVertexLinks = precedingVertexLinksAcc(*m_b);
1295
1296 // ... and check if they are all valid
1297 VertexLinkVector::const_iterator precedingVertexLinksItr = precedingVertexLinks.begin();
1298 for(; precedingVertexLinksItr!=precedingVertexLinks.end(); ++precedingVertexLinksItr) {
1299 // check if links are valid
1300 if(!(*precedingVertexLinksItr).isValid()) {
1301 return false;
1302 }
1303 }
1304
1305 // all OK, cache the precedingVertices
1306 precedingVertexLinksItr = precedingVertexLinks.begin();
1307 for(; precedingVertexLinksItr!=precedingVertexLinks.end(); ++precedingVertexLinksItr) {
1308 m_cachedPrecedingVertices.push_back(*(*precedingVertexLinksItr));
1309 }
1310
1311 return true;
1312
1313}

◆ cacheRefTracks()

bool xAOD::BPhysHelper::cacheRefTracks ( )
protected

: Cache refitted tracks

To speed up access to refitted track momenta, m_cachedRefTracks vector is created first time refitted track accessor methods are called. In subsequent calls, cached vector is used.

Returns
: true on success

Definition at line 1137 of file BPhysHelper.cxx.

1138{
1139 // if refitted tracks are already cached, do nothing and return true
1140 if(m_refTracksCached && !m_cachedRefTracks.empty())
1141 return true;
1142
1143 // if refitted tracks are cached but the cache is empty, return false (error)
1145 return false;
1146
1147 // Refitted tracks are not cached, yet. Retrieve them from the aux store:
1148 m_refTracksCached = true;
1149 m_cachedRefTracks.clear();
1150
1151 // Create auxiliary branches accessors
1152 static const SG::AuxElement::Accessor< std::vector<float> > refTrackPxAcc("RefTrackPx");
1153 static const SG::AuxElement::Accessor< std::vector<float> > refTrackPyAcc("RefTrackPy");
1154 static const SG::AuxElement::Accessor< std::vector<float> > refTrackPzAcc("RefTrackPz");
1155
1156 // check if branches are available:
1157 if(!refTrackPxAcc.isAvailable(*m_b) ||
1158 !refTrackPyAcc.isAvailable(*m_b) ||
1159 !refTrackPzAcc.isAvailable(*m_b) )
1160 return false;
1161
1162 // retrieve branches:
1163 const std::vector<float>& refTrackPx = refTrackPxAcc(*m_b);
1164 const std::vector<float>& refTrackPy = refTrackPyAcc(*m_b);
1165 const std::vector<float>& refTrackPz = refTrackPzAcc(*m_b);
1166
1167 // all must have the same lenght
1168 if(refTrackPx.size()!=refTrackPy.size() || refTrackPx.size()!=refTrackPz.size())
1169 return false;
1170
1171 // all OK: store refitted tracks in the cache
1172 for(uint i=0; i<refTrackPx.size(); ++i) {
1173 m_cachedRefTracks.emplace_back(refTrackPx[i], refTrackPy[i], refTrackPz[i] );
1174 }
1175
1176 return true;
1177
1178}
unsigned int uint

◆ cascadeVertex()

const xAOD::Vertex * xAOD::BPhysHelper::cascadeVertex ( const size_t index)

Returns pointer to a cascade vertex.

Parameters
[in]indexindex of the cascade vertex
Returns
: pointer to the vertex, 0 on error

Definition at line 694 of file BPhysHelper.cxx.

695{
696 // cache linked cascadeVertices
698 return nullptr;
699
700 // range check
701 if(index>=m_cachedCascadeVertices.size())
702 return nullptr;
703
704 // all OK:
706
707}
bool cacheCascadeVertices()
: Cache cascade vertices
str index
Definition DeMoScan.py:362

◆ cascadeVertices()

const std::vector< const xAOD::Vertex * > & xAOD::BPhysHelper::cascadeVertices ( )

Returns vector of pointers to cascade vertices.

Returns
: vector of pointers to vertices, empty vector on error

Definition at line 710 of file BPhysHelper.cxx.

711{
712 // cache linked cascadeVertices
715
716 // all OK:
718
719}
static const std::vector< const xAOD::Vertex * > s_emptyVectorOfVertices

◆ covariance()

const TMatrixTSym< double > & xAOD::BPhysHelper::covariance ( )

: Returns full covariance matrix

Converts native covariance matrix stored as a std::vector<float> into a symmetric matrix-form represented by TMatrixTSym<double> class.

Note
This method needs to know the number of tracks. If refitted tracks are stripped from the verrtex, this method will no longer work.
Returns
: covariance matrix, 0x0 matrix in case of an error

Definition at line 103 of file BPhysHelper.cxx.

104{
105
106 // cache covariance matrix
107 if(!cacheCov())
109
110 // all OK:
111 return m_cachedCov;
112}
static const TMatrixTSym< double > s_emptyMatrix
bool cacheCov()
: Cache covariance matrix

◆ electron()

const xAOD::Electron * xAOD::BPhysHelper::electron ( const size_t index)

Definition at line 532 of file BPhysHelper.cxx.

533{
534 // cache linked electrons
535 if(!cacheElectrons())
536 return nullptr;
537
538 // range check
539 if(index>=m_cachedElectrons.size())
540 return nullptr;
541
542 // all OK:
543 return m_cachedElectrons[index];
544
545}

◆ electrons()

const std::vector< const xAOD::Electron * > & xAOD::BPhysHelper::electrons ( )

Definition at line 548 of file BPhysHelper.cxx.

549{
550 // cache linked electrons
551 if(!cacheElectrons())
553
554 // all OK:
555 return m_cachedElectrons;
556}
static const std::vector< const xAOD::Electron * > s_emptyVectorOfElectrons

◆ lxy()

float xAOD::BPhysHelper::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 pv_type and the B vertex.

Transverse decay distance and its error

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Returns
: lxy, -9999999. in case augmentation is not available decay distance

Definition at line 866 of file BPhysHelper.cxx.

867{
868 switch(vertexType) {
869 case PV_MAX_SUM_PT2 : GET_FLOAT("LxyMaxSumPt2");
870 case PV_MIN_A0 : GET_FLOAT("LxyMinA0");
871 case PV_MIN_Z0 : GET_FLOAT("LxyMinZ0");
872 case PV_MIN_Z0_BA : GET_FLOAT("LxyMinZ0BA");
873 default: return -9999999.;
874 }
875}

◆ lxyErr()

float xAOD::BPhysHelper::lxyErr ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

its error

Definition at line 877 of file BPhysHelper.cxx.

878{
879 switch(vertexType) {
880 case PV_MAX_SUM_PT2 : GET_FLOAT("LxyErrMaxSumPt2");
881 case PV_MIN_A0 : GET_FLOAT("LxyErrMinA0");
882 case PV_MIN_Z0 : GET_FLOAT("LxyErrMinZ0");
883 case PV_MIN_Z0_BA : GET_FLOAT("LxyErrMinZ0BA");
884 default: return -9999999.;
885 }
886}

◆ lxyz()

float xAOD::BPhysHelper::lxyz ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

decay distance

Definition at line 912 of file BPhysHelper.cxx.

913{
914 switch(vertexType) {
915 case PV_MAX_SUM_PT2 : GET_FLOAT("LxyzMaxSumPt2");
916 case PV_MIN_A0 : GET_FLOAT("LxyzMinA0");
917 case PV_MIN_Z0 : GET_FLOAT("LxyzMinZ0");
918 case PV_MIN_Z0_BA : GET_FLOAT("LxyzMinZ0BA");
919 default: return -9999999.;
920 }
921}

◆ lxyzErr()

float xAOD::BPhysHelper::lxyzErr ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

its error

Definition at line 923 of file BPhysHelper.cxx.

924{
925 switch(vertexType) {
926 case PV_MAX_SUM_PT2 : GET_FLOAT("LxyzErrMaxSumPt2");
927 case PV_MIN_A0 : GET_FLOAT("LxyzErrMinA0");
928 case PV_MIN_Z0 : GET_FLOAT("LxyzErrMinZ0");
929 case PV_MIN_Z0_BA : GET_FLOAT("LxyzErrMinZ0BA");
930 default: return -9999999.;
931 }
932}

◆ muon()

const xAOD::Muon * xAOD::BPhysHelper::muon ( const size_t index)

Returns pointer to the i-th linked muon.

Parameters
[in]indexindex of the linked muon
Returns
: Pointer to the muon, 0 in case of an error

Definition at line 452 of file BPhysHelper.cxx.

453{
454 // cache linked muons
455 if(!cacheMuons())
456 return nullptr;
457
458 // range check
459 if(index>=m_cachedMuons.size())
460 return nullptr;
461
462 // all OK:
463 return m_cachedMuons[index];
464
465}
bool cacheMuons()
: Cache linked muons

◆ muons()

const std::vector< const xAOD::Muon * > & xAOD::BPhysHelper::muons ( )

Returns linked muons.

Returns
: linked muons, empty vector in case of an error

Definition at line 468 of file BPhysHelper.cxx.

469{
470 // cache linked muons
471 if(!cacheMuons())
473
474 // all OK:
475 return m_cachedMuons;
476}
static const std::vector< const xAOD::Muon * > s_emptyVectorOfMuons

◆ nCascadeVertices()

int xAOD::BPhysHelper::nCascadeVertices ( )

: Links to cascade vertices

Cascade decays, e.g. Bd->Jpsi(mumu)+Ks(pipi), consist of more than one reconstructed vertex. for instance in case of the Bd decay, 2 vertices are required. Cascade vertices links are used to store relations between mother and daughter vertices in the cascade decays. Returns number of cascade vertices

Returns
: number of cascade vertices, -1 on error

Definition at line 682 of file BPhysHelper.cxx.

683{
684 // cache linked cascade vertices
686 return -1;
687
688 // all OK:
689 return m_cachedCascadeVertices.size();
690
691}

◆ nElectrons()

int xAOD::BPhysHelper::nElectrons ( )

Definition at line 521 of file BPhysHelper.cxx.

522{
523 // cache linked electrons
524 if(!cacheElectrons())
525 return -1;
526
527 // all OK:
528 return m_cachedElectrons.size();
529}

◆ nMuons()

int xAOD::BPhysHelper::nMuons ( )

: Methods providing access to the linked muons

By default, xAOD::Vertex provides only links to the TrackParticles used to build the vertex. Links to muons (e.g. for Jpsi->mu+mu) are stored using auxiliary branches. Returns number of stored muon links

Returns
: number of muons, -1 in case of an error

Definition at line 441 of file BPhysHelper.cxx.

442{
443 // cache linked muons
444 if(!cacheMuons())
445 return -1;
446
447 // all OK:
448 return m_cachedMuons.size();
449}

◆ nPrecedingVertices()

int xAOD::BPhysHelper::nPrecedingVertices ( )

: Links to preceding vertices

The B vertices are often build from tracks that have already been pre-selected by requirements based on some preceding vertex fit. For instance, B --> J/psi(mumu)+K vertex is only build using muon tracks that have been successfully refitted in a J/psi->mumu vertex.

The following methods enable to retrieve/store links to the vertices that were preceded this one in the reconstruction flow. Returns number of preceding vertices

Returns
: number of preceding vertices, -1 on error

Definition at line 601 of file BPhysHelper.cxx.

602{
603 // cache linked preceding vertices
605 return -1;
606
607 // all OK:
608 return m_cachedPrecedingVertices.size();
609
610}
bool cachePrecedingVertices()
: Cache preceding vertices

◆ nRefTrks()

int xAOD::BPhysHelper::nRefTrks ( )

Returns number of stored refitted track momenta.

Methods providing access to the refitted tracks momenta and related quantities:

Returns
: number of refitted tracks, -1 in case of an error

Definition at line 115 of file BPhysHelper.cxx.

116{
117 // cache refitted tracks
118 if( !cacheRefTracks() ) return -1;
119
120 // all OK:
121 return m_cachedRefTracks.size();
122
123}
bool cacheRefTracks()
: Cache refitted tracks

◆ nRefTrksCascade()

int xAOD::BPhysHelper::nRefTrksCascade ( )

Returns number of stored refitted tracks INCLUDING those from the linked cascade vertices.

Returns
: number of refitted tracks, -1 in case of an error

Definition at line 764 of file BPhysHelper.cxx.

765{
766 // get number of refitted tracks in this vertex
767 int n = nRefTrks();
768
769 // check for errors. If there are cascade vertices stored in this vertex,
770 // it is possible there are no refitted tracks, so report error only in case
771 // when this veretx has no cascade vertices and no refitted tracks:
772 if(n<0 && nCascadeVertices()<=0) {
773 return -1;
774 } else if(n<0 && nCascadeVertices()>0) {
775 n = 0;
776 }
777
778 // loop over cascade vertices, if they exist:
779 for(int i=0; i<nCascadeVertices(); ++i) {
780 // create helper class for the daughter vertex
781 xAOD::BPhysHelper daughterVtx(cascadeVertex(i));
782
783 // check for errors
784 if(daughterVtx.nRefTrksCascade()<0)
785 return -1;
786
787 // OK, add the numbers
788 n += daughterVtx.nRefTrksCascade();
789 }
790
791
792 // all OK, return the number of tracks
793 return n;
794}
int nCascadeVertices()
: Links to cascade vertices
const xAOD::Vertex * cascadeVertex(const size_t index)
Returns pointer to a cascade vertex.

◆ origPv()

const xAOD::Vertex * xAOD::BPhysHelper::origPv ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

original PV

Definition at line 807 of file BPhysHelper.cxx.

808{
809 switch(vertexType) {
810 case PV_MAX_SUM_PT2 : GET_PV("OrigPvMaxSumPt2Link");
811 case PV_MIN_A0 : GET_PV("OrigPvMinA0Link");
812 case PV_MIN_Z0 : GET_PV("OrigPvMinZ0Link");
813 case PV_MIN_Z0_BA : GET_PV("OrigPvMinZ0BALink");
814 default: return nullptr;
815 }
816}
#define GET_PV(name)

◆ precedingVertex()

const xAOD::Vertex * xAOD::BPhysHelper::precedingVertex ( const size_t index)

Returns pointer to a preceding vertex.

Parameters
[in]indexindex of the preceding vertex
Returns
: pointer to the vertex, 0 on error

Definition at line 613 of file BPhysHelper.cxx.

614{
615 // cache linked precedingVertices
617 return nullptr;
618
619 // range check
620 if(index>=m_cachedPrecedingVertices.size())
621 return nullptr;
622
623 // all OK:
625
626}

◆ precedingVertices()

const std::vector< const xAOD::Vertex * > & xAOD::BPhysHelper::precedingVertices ( )

Returns vector of pointers to preceding vertices.

Returns
: vector of pointers to vertices, empty vector on error

Definition at line 629 of file BPhysHelper.cxx.

630{
631 // cache linked precedingVertices
634
635 // all OK:
637
638}

◆ ptErr()

float xAOD::BPhysHelper::ptErr ( )

Returns pT error.

Returns
: pT error, -9999999. if augmentation is not available

Definition at line 431 of file BPhysHelper.cxx.

432{
433 GET_FLOAT("PtErr");
434}

◆ pv()

const xAOD::Vertex * xAOD::BPhysHelper::pv ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

Get the refitted collision vertex of type pv_type.

Links to the refitted primary vertices

Todo
: In case the refitted vertices are augmented with links to the original primary vertices, method origPv will resolve the links and return the original PVs.
Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Returns
: pointer to the PV, 0 if not stored refitted PV

Definition at line 796 of file BPhysHelper.cxx.

797{
798 switch(vertexType) {
799 case PV_MAX_SUM_PT2 : GET_PV("PvMaxSumPt2Link");
800 case PV_MIN_A0 : GET_PV("PvMinA0Link");
801 case PV_MIN_Z0 : GET_PV("PvMinZ0Link");
802 case PV_MIN_Z0_BA : GET_PV("PvMinZ0BALink");
803 default: return nullptr;
804 }
805}

◆ RefitPVStatus()

int xAOD::BPhysHelper::RefitPVStatus ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

Get the exitCode of the refitter for vertex of type pv_type.

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Returns
: the code 0=refit not request, 1=completed sucessfully, >1 not refitted for good reasons, <0 not refitted for bad reasons (see PrimaryVertexRefitter.cxx for more details)

Definition at line 855 of file BPhysHelper.cxx.

856{
857 switch(vertexType) {
858 case PV_MAX_SUM_PT2 : GET_INT("PvMaxSumPt2Status");
859 case PV_MIN_A0 : GET_INT("PvMinA0Status" );
860 case PV_MIN_Z0 : GET_INT("PvMinZ0Status" );
861 case PV_MIN_Z0_BA : GET_INT("PvMinZ0BAStatus" );
862 default: return -999999;
863 }
864}
#define GET_INT(name)

◆ refTrk() [1/2]

TVector3 xAOD::BPhysHelper::refTrk ( const size_t index)

Returns i-th refitted track 3-momentum.

Parameters
[in]indexindex of the refitted track
Returns
: refitted track momentum, (0,0,0) in case of an error

Definition at line 126 of file BPhysHelper.cxx.

127{
128 // cache refitted tracks
129 if( !cacheRefTracks() ) return
130 TVector3(0,0,0);
131
132 // range check:
133 if(index>=m_cachedRefTracks.size())
134 return TVector3(0,0,0);
135
136 // all OK:
137 return m_cachedRefTracks[index];
138
139}

◆ refTrk() [2/2]

TLorentzVector xAOD::BPhysHelper::refTrk ( const size_t index,
const float mass )

Returns i-th refitted track as a 4-vector given the mass hypothesis.

Parameters
[in]indexindex of the refitted track
[in]massmass hypothesis assigned to the refitted track
Returns
: refitted track 4-momentum, (0,0,0,0) in case of an error

Definition at line 153 of file BPhysHelper.cxx.

154{
155 // cache refitted tracks
156 if( !cacheRefTracks() ) return
157 TLorentzVector(0,0,0,0);
158
159 // range check:
160 if(index>=m_cachedRefTracks.size())
161 return TLorentzVector(0,0,0,0);
162
163 // all OK, create the 4-momentum
164 TLorentzVector mom;
165 mom.SetVectM(m_cachedRefTracks[index], mass);
166 return mom;
167
168}

◆ refTrkCharge()

float xAOD::BPhysHelper::refTrkCharge ( const size_t index) const

Returns charge of the i-th track.

Note
This method needs access to original TrackParticles (and/or NeutralParticles). If the collection(s) of original tracks has(have) been stripped, this method will return the error value.
Parameters
[in]indexindex of the track
Returns
: 1. for positive track, -1. for negative, and -9999. on error

Definition at line 255 of file BPhysHelper.cxx.

256{
257 // range check:
258 // FIXME: add neutral particles once they are available
259 //if(index>=m_b->nTrackParticles()+m_b->nNeutralParticles())
260 // return 0;
261
262 // FIXME:supports charged particles only for now
263 if(index>=m_b->nTrackParticles())
264 return -9999.;
265
266 if(index<m_b->nTrackParticles()) {
267 // charged track
268 const xAOD::TrackParticle* trk = dynamic_cast<const xAOD::TrackParticle*>(refTrkOrigin(index));
269
270 // null-pointer protection
271 if(!trk)
272 return -9999.;
273
274 // all OK
275 return trk->charge();
276
277 }else{
278 // we do not need to do anything here, since we know this track must be neutral
279 return 0.;
280 }
281
282
283}
const xAOD::IParticle * refTrkOrigin(const size_t index) const
: Returns the original track (charged or neutral) corresponding to the i-th refitted track
float charge() const
Returns the charge.
TrackParticle_v1 TrackParticle
Reference the current persistent version:

◆ refTrkOrigin()

const xAOD::IParticle * xAOD::BPhysHelper::refTrkOrigin ( const size_t index) const

: Returns the original track (charged or neutral) corresponding to the i-th refitted track

In xAOD model, charged and neutral particles are represented by two distict classes: xAOD::TrackParticle and xAOD::NeutralParticle. Vertex fitters can use both as input and store charded tracks first (needs checking!!!) in the list of refitted momenta.

This method returns pointer to the original track corresponding to the i-th refitted track, i.e. either TrackParticle (if index<nTrackParticles) or NeutralParticle (if index>=nTrackParticles).

One has to dynamicly cast the returned value to either TrackParticle or NeutralParticle. One can use refTrkCharge to determine if the refitted track is charged or neutral.

Note
This method needs access to original TrackParticles (and/or NeutralParticles). If the collection(s) of original tracks has(have) been stripped, this method will return the error value.
Parameters
[in]indexindex of the refitted track
Returns
: pointer to the original TrackParticle (if index<nTrackParticles), pointer to the original NeutralParticle (if index>=nNeutralParticles), 0 in case of an error

Definition at line 171 of file BPhysHelper.cxx.

172{
173 // range check:
174 // FIXME: add neutral particles once they are available
175 //if(index>=m_b->nTrackParticles()+m_b->nNeutralParticles())
176 // return 0;
177
178 // FIXME:supports charged particles only for now
179 if(index>=m_b->nTrackParticles())
180 return nullptr;
181
182 // charged tracks
183 if(index<m_b->nTrackParticles()) {
184 return m_b->trackParticle(index);
185 } else {
186 // FIXME: add neutral particles once they are available
187 //return m_b->neutralParticle(index - m_b->nNeutralParticles());
188 return nullptr;
189 }
190
191}

◆ refTrkOriginP() [1/2]

TVector3 xAOD::BPhysHelper::refTrkOriginP ( const size_t index) const

Returns perigee 3-momentum of the original track corresponding i-th refitted track.

Note
This method needs access to original TrackParticles (and/or NeutralParticles). If the collection(s) of original tracks has(have) been stripped, this method will return the error value.
Parameters
[in]indexindex of the refitted track
Returns
: Original perigee momentum, (0,0,0) in case of an error

Definition at line 194 of file BPhysHelper.cxx.

195{
196 // range check:
197 // FIXME: add neutral particles once they are available
198 //if(index>=m_b->nTrackParticles()+m_b->nNeutralParticles())
199 // return TVector3(0,0,0);
200
201 // FIXME:supports charged particles only for now
202 if(index>=m_b->nTrackParticles())
203 return TVector3(0,0,0);
204
205 if(index<m_b->nTrackParticles()) {
206 // charged track
207 const xAOD::TrackParticle* trk = dynamic_cast<const xAOD::TrackParticle*>(refTrkOrigin(index));
208
209 // null-pointer protection
210 if(!trk) return TVector3(0,0,0);
211
212 // all OK
213 TVector3 tmp;
214 tmp.SetPtEtaPhi(trk->pt(), trk->eta(), trk->phi());
215 return tmp;
216
217 }else{
218 // // neutral track
219 // const xAOD::NeutralParticle* trk = dynamic_cast<const xAOD::NeutralParticle*>(refTrkOrigin(index));
220
221 // // null-pointer protection
222 // if(!trk) return TVector3(0,0,0);
223
224 // // all OK
225 // TVector3 tmp;
226 // tmp.SetPtEtaPhi(trk->pt(), trk->eta(), trk->phi());
227 // return tmp;
228
229 // FIXME neutral tracks not implemented yet:
230 return TVector3(0,0,0);
231 }
232
233}
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.

◆ refTrkOriginP() [2/2]

TLorentzVector xAOD::BPhysHelper::refTrkOriginP ( const size_t index,
const float mass ) const

: Returns lorentz vector build from the original track perigee momentum and given mass hypothesis

The lorentz vector is created from the perigee parameters of the original track (neutral or charged) which corresponds to the i-th refitted track and a given mass hypothesis

Note
This method needs access to original TrackParticles (and/or NeutralParticles). If the collection(s) of original tracks has(have) been stripped, this method will return the error value.
Parameters
[in]indexindex of the refitted track
[in]massmass hypothesis of the track
Returns
: Original perigee momentum, (0,0,0,0) in case of an error

Definition at line 236 of file BPhysHelper.cxx.

237{
238 // range check:
239 // FIXME: add neutral particles once they are available
240 //if(index>=m_b->nTrackParticles()+m_b->nNeutralParticles())
241 // return TVector3(0,0,0);
242
243 // FIXME:supports charged particles only for now
244 if(index>=m_b->nTrackParticles())
245 return TLorentzVector(0,0,0,0);
246
247 // create TLorentzVector
248 TLorentzVector tmp;
249 tmp.SetVectM(refTrkOriginP(index), mass);
250 return tmp;
251
252}
TVector3 refTrkOriginP(const size_t index) const
Returns perigee 3-momentum of the original track corresponding i-th refitted track.

◆ refTrks()

const std::vector< TVector3 > & xAOD::BPhysHelper::refTrks ( )

Returns refitted track momenta.

Note
The std::vector is returned as a reference! It will cease to exist when BPhysHelper class is destroyed.
Returns
: refitted track, empty vector in case of an error

Definition at line 142 of file BPhysHelper.cxx.

143{
144 // cache refitted tracks
145 if( !cacheRefTracks() )
147
148 // return
149 return m_cachedRefTracks;
150}
static const std::vector< TVector3 > s_emptyVectorOfTVector3

◆ setA0()

float xAOD::BPhysHelper::setA0 ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

Set the 3D and transverse impact parameters and their error.

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Parameters
[in]valinput value
Returns
: true on success 3D impact parameter

Definition at line 1023 of file BPhysHelper.cxx.

1024{
1025 switch(vertexType) {
1026 case PV_MAX_SUM_PT2 : SET_FLOAT("A0MaxSumPt2", val);
1027 case PV_MIN_A0 : SET_FLOAT("A0MinA0" , val);
1028 case PV_MIN_Z0 : SET_FLOAT("A0MinZ0" , val);
1029 case PV_MIN_Z0_BA : SET_FLOAT("A0MinZ0BA" , val);
1030 default: return false;
1031 }
1032}
#define SET_FLOAT(name, val)

◆ setA0Err()

float xAOD::BPhysHelper::setA0Err ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

3D impact parameter error

Definition at line 1034 of file BPhysHelper.cxx.

1035{
1036 switch(vertexType) {
1037 case PV_MAX_SUM_PT2 : SET_FLOAT("A0ErrMaxSumPt2", val);
1038 case PV_MIN_A0 : SET_FLOAT("A0ErrMinA0" , val);
1039 case PV_MIN_Z0 : SET_FLOAT("A0ErrMinZ0" , val);
1040 case PV_MIN_Z0_BA : SET_FLOAT("A0ErrMinZ0BA" , val);
1041 default: return false;
1042 }
1043}

◆ setA0xy()

float xAOD::BPhysHelper::setA0xy ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

transverse impact parameter

Definition at line 1045 of file BPhysHelper.cxx.

1046{
1047 switch(vertexType) {
1048 case PV_MAX_SUM_PT2 : SET_FLOAT("A0xyMaxSumPt2", val);
1049 case PV_MIN_A0 : SET_FLOAT("A0xyMinA0" , val);
1050 case PV_MIN_Z0 : SET_FLOAT("A0xyMinZ0" , val);
1051 case PV_MIN_Z0_BA : SET_FLOAT("A0xyMinZ0BA" , val);
1052 default: return false;
1053 }
1054}

◆ setA0xyErr()

float xAOD::BPhysHelper::setA0xyErr ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

transverse impact parameter error

Definition at line 1056 of file BPhysHelper.cxx.

1057{
1058 switch(vertexType) {
1059 case PV_MAX_SUM_PT2 : SET_FLOAT("A0xyErrMaxSumPt2", val);
1060 case PV_MIN_A0 : SET_FLOAT("A0xyErrMinA0" , val);
1061 case PV_MIN_Z0 : SET_FLOAT("A0xyErrMinZ0" , val);
1062 case PV_MIN_Z0_BA : SET_FLOAT("A0xyErrMinZ0BA" , val);
1063 default: return false;
1064 }
1065}

◆ setCascadeVertices()

bool xAOD::BPhysHelper::setCascadeVertices ( const std::vector< const xAOD::Vertex * > & vertices,
const xAOD::VertexContainer * vertexContainer )

Sets links to cascade vertices.

Parameters
[in]indexindex of the cascade vertex
Returns
: pointer to the vertex, 0 on error

Definition at line 722 of file BPhysHelper.cxx.

724{
725 // erase cascadeVertex cache to preserve consistency
728
729 // Create cascade vertex links decorator
730 static const SG::AuxElement::Decorator<VertexLinkVector> cascadeVertexLinksDecor("CascadeVertexLinks");
731
732 // create tmp vector of cascade vertex links
733 VertexLinkVector cascadeVertexLinks;
734
735 // loop over input cascadeVertices
736 std::vector<const xAOD::Vertex*>::const_iterator cascadeVerticesItr = vertices.begin();
737 for(; cascadeVerticesItr!=vertices.end(); ++cascadeVerticesItr) {
738 // sanity check 1: protect against null pointers
739 if( !(*cascadeVerticesItr) )
740 return false;
741
742 // create element link
743 VertexLink vertexLink;
744 vertexLink.setElement(*cascadeVerticesItr);
745 vertexLink.setStorableObject(*vertexContainer);
746
747 // sanity check 2: is the link valid?
748 if( !vertexLink.isValid() )
749 return false;
750
751 // link is OK, store it in the tmp vector
752 cascadeVertexLinks.push_back( vertexLink );
753
754 } // end of loop over cascade vertices
755
756 // all OK: store cascade vertex links in the aux store
757 cascadeVertexLinksDecor(*m_b) = std::move(cascadeVertexLinks);
758
759 return true;
760
761}
ElementLink< xAOD::VertexContainer > VertexLink
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575

◆ setElectrons()

bool xAOD::BPhysHelper::setElectrons ( const std::vector< const xAOD::Electron * > & electrons,
const xAOD::ElectronContainer * electronContainer )

Definition at line 559 of file BPhysHelper.cxx.

561{
562 // erase electron cache to preserve cosistency
563 m_electronsCached = false;
564 m_cachedElectrons.clear();
565
566 // Create electron links decorator
567 static const SG::AuxElement::Decorator<ElectronLinkVector> electronLinksDecor("ElectronLinks");
568
569 // create tmp vector of electron links
570 ElectronLinkVector electronLinks;
571
572 // loop over input electrons
573 std::vector<const xAOD::Electron*>::const_iterator electronsItr = electrons.begin();
574 for(; electronsItr!=electrons.end(); ++electronsItr) {
575 // sanity check 1: protect against null pointers
576 if( !(*electronsItr) )
577 return false;
578
579 // create element link
580 ElementLink<xAOD::ElectronContainer> elLink;
581 elLink.setElement(*electronsItr);
583
584 // sanity check 2: is the link valid?
585 if( !elLink.isValid() )
586 return false;
587
588 // link is OK, store it in the tmp vector
589 electronLinks.push_back( elLink );
590
591 } // end of loop over electrons
592
593 // all OK: store electron links in the aux store
594 electronLinksDecor(*m_b) = std::move(electronLinks);
595
596 return true;
597
598}
xAOD::ElectronContainer * electronContainer
const std::vector< const xAOD::Electron * > & electrons()

◆ setLxy()

bool xAOD::BPhysHelper::setLxy ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

Set the transverse decay distance and its error measured between the refitted primary vertex of type pv_type and the B vertex.

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Parameters
[in]valinput value
Returns
: true on success decay distance

Definition at line 888 of file BPhysHelper.cxx.

889{
890 switch(vertexType) {
891 case PV_MAX_SUM_PT2 : SET_FLOAT("LxyMaxSumPt2", val);
892 case PV_MIN_A0 : SET_FLOAT("LxyMinA0" , val);
893 case PV_MIN_Z0 : SET_FLOAT("LxyMinZ0" , val);
894 case PV_MIN_Z0_BA : SET_FLOAT("LxyMinZ0BA" , val);
895 default: return false;
896 }
897}

◆ setLxyErr()

bool xAOD::BPhysHelper::setLxyErr ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

its error

Definition at line 899 of file BPhysHelper.cxx.

900{
901 switch(vertexType) {
902 case PV_MAX_SUM_PT2 : SET_FLOAT("LxyErrMaxSumPt2", val);
903 case PV_MIN_A0 : SET_FLOAT("LxyErrMinA0" , val);
904 case PV_MIN_Z0 : SET_FLOAT("LxyErrMinZ0" , val);
905 case PV_MIN_Z0_BA : SET_FLOAT("LxyErrMinZ0BA" , val);
906 default: return false;
907 }
908}

◆ setLxyz()

bool xAOD::BPhysHelper::setLxyz ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

decay distance

Definition at line 934 of file BPhysHelper.cxx.

935{
936 switch(vertexType) {
937 case PV_MAX_SUM_PT2 : SET_FLOAT("LxyzMaxSumPt2", val);
938 case PV_MIN_A0 : SET_FLOAT("LxyzMinA0" , val);
939 case PV_MIN_Z0 : SET_FLOAT("LxyzMinZ0" , val);
940 case PV_MIN_Z0_BA : SET_FLOAT("LxyzMinZ0BA" , val);
941 default: return false;
942 }
943}

◆ setLxyzErr()

bool xAOD::BPhysHelper::setLxyzErr ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

its error

Definition at line 945 of file BPhysHelper.cxx.

946{
947 switch(vertexType) {
948 case PV_MAX_SUM_PT2 : SET_FLOAT("LxyzErrMaxSumPt2", val);
949 case PV_MIN_A0 : SET_FLOAT("LxyzErrMinA0" , val);
950 case PV_MIN_Z0 : SET_FLOAT("LxyzErrMinZ0" , val);
951 case PV_MIN_Z0_BA : SET_FLOAT("LxyzErrMinZ0BA" , val);
952 default: return false;
953 }
954}

◆ setMuons()

bool xAOD::BPhysHelper::setMuons ( const std::vector< const xAOD::Muon * > & muons,
const xAOD::MuonContainer * muonContainer )

Set links to muons.

Parameters
[in]muonsstd::vector of muons to be linked to this vertex
[in]muonContainercollection these muons belong to
Returns
: true on success

Definition at line 479 of file BPhysHelper.cxx.

481{
482 // erase muon cache to preserve cosistency
483 m_muonsCached = false;
484 m_cachedMuons.clear();
485
486 // Create muon links decorator
487 static const SG::AuxElement::Decorator<MuonLinkVector> muonLinksDecor("MuonLinks");
488
489 // create tmp vector of muon links
490 MuonLinkVector muonLinks;
491
492 // loop over input muons
493 std::vector<const xAOD::Muon*>::const_iterator muonsItr = muons.begin();
494 for(; muonsItr!=muons.end(); ++muonsItr) {
495 // sanity check 1: protect against null pointers
496 if( !(*muonsItr) )
497 return false;
498
499 // create element link
500 ElementLink<xAOD::MuonContainer> muLink;
501 muLink.setElement(*muonsItr);
503
504 // sanity check 2: is the link valid?
505 if( !muLink.isValid() )
506 return false;
507
508 // link is OK, store it in the tmp vector
509 muonLinks.push_back( muLink );
510
511 } // end of loop over muons
512
513 // all OK: store muon links in the aux store
514 muonLinksDecor(*m_b) = std::move(muonLinks);
515
516 return true;
517
518}
xAOD::MuonContainer * muonContainer
const std::vector< const xAOD::Muon * > & muons()
Returns linked muons.

◆ setOrigPv()

bool xAOD::BPhysHelper::setOrigPv ( const xAOD::Vertex * pv,
const xAOD::VertexContainer * vertexContainer,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

Set the original collision vertex of type pv_type.

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Parameters
[in]pvpointer to the PV
[in]vertexContainercollection in which the PV is stored
Returns
: true on success, false on failure

Definition at line 818 of file BPhysHelper.cxx.

821{
822 switch(vertexType) {
823 case PV_MAX_SUM_PT2 : SET_PV("OrigPvMaxSumPt2Link", pv, vertexContainer);
824 case PV_MIN_A0 : SET_PV("OrigPvMinA0Link" , pv, vertexContainer);
825 case PV_MIN_Z0 : SET_PV("OrigPvMinZ0Link" , pv, vertexContainer);
826 case PV_MIN_Z0_BA : SET_PV("OrigPvMinZ0BALink" , pv, vertexContainer);
827 default: return 0;
828 }
829}
#define SET_PV(name, pv, vertexContainer)
const xAOD::Vertex * pv(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Get the refitted collision vertex of type pv_type.

◆ setPrecedingVertices()

bool xAOD::BPhysHelper::setPrecedingVertices ( const std::vector< const xAOD::Vertex * > & vertices,
const xAOD::VertexContainer * vertexContainer )

Sets links to preceding vertices.

Parameters
[in]indexindex of the preceding vertex
Returns
: pointer to the vertex, 0 on error

Definition at line 641 of file BPhysHelper.cxx.

643{
644 // erase precedingVertex cache to preserve consistency
647
648 // Create preceding vertex links decorator
649 static const SG::AuxElement::Decorator<VertexLinkVector> precedingVertexLinksDecor("PrecedingVertexLinks");
650
651 // create tmp vector of preceding vertex links
652 VertexLinkVector precedingVertexLinks;
653
654 // loop over input precedingVertices
655 std::vector<const xAOD::Vertex*>::const_iterator precedingVerticesItr = vertices.begin();
656 for(; precedingVerticesItr!=vertices.end(); ++precedingVerticesItr) {
657 // sanity check 1: protect against null pointers
658 if( !(*precedingVerticesItr) )
659 return false;
660
661 // create element link
662 VertexLink vertexLink;
663 vertexLink.setElement(*precedingVerticesItr);
664 vertexLink.setStorableObject(*vertexContainer);
665
666 // sanity check 2: is the link valid?
667 if( !vertexLink.isValid() )
668 return false;
669
670 // link is OK, store it in the tmp vector
671 precedingVertexLinks.push_back( vertexLink );
672
673 } // end of loop over preceding vertices
674
675 // all OK: store preceding vertex links in the aux store
676 precedingVertexLinksDecor(*m_b) = std::move(precedingVertexLinks);
677
678 return true;
679
680}

◆ setPtErr()

bool xAOD::BPhysHelper::setPtErr ( const float val)

Set pT error.

Parameters
[in]inputvalue
Returns
: true on success

Definition at line 436 of file BPhysHelper.cxx.

437{
438 SET_FLOAT("PtErr", val);
439}

◆ setPv()

bool xAOD::BPhysHelper::setPv ( const xAOD::Vertex * pv,
const xAOD::VertexContainer * vertexContainer,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

Set the refitted collision vertex of type pv_type.

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Parameters
[in]pvpointer to the PV
[in]vertexContainercollection in which the PV is stored
Returns
: true on success, false on failure

Definition at line 831 of file BPhysHelper.cxx.

834{
835 switch(vertexType) {
836 case PV_MAX_SUM_PT2 : SET_PV("PvMaxSumPt2Link", pv, vertexContainer);
837 case PV_MIN_A0 : SET_PV("PvMinA0Link" , pv, vertexContainer);
838 case PV_MIN_Z0 : SET_PV("PvMinZ0Link" , pv, vertexContainer);
839 case PV_MIN_Z0_BA : SET_PV("PvMinZ0BALink" , pv, vertexContainer);
840 default: return false;
841 }
842}

◆ setRefitPVStatus()

bool xAOD::BPhysHelper::setRefitPVStatus ( int code,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

Set the exitCode of the refitter for vertex of type pv_type.

Parameters
[in]vertexTypetype (selection) of the PV (
See also
: pv_type)
Parameters
[in]codeint containing the code
Returns
: true on success, false on failure

Definition at line 844 of file BPhysHelper.cxx.

845{
846 switch(vertexType) {
847 case PV_MAX_SUM_PT2 : SET_INT("PvMaxSumPt2Status", code);
848 case PV_MIN_A0 : SET_INT("PvMinA0Status" , code);
849 case PV_MIN_Z0 : SET_INT("PvMinZ0Status" , code);
850 case PV_MIN_Z0_BA : SET_INT("PvMinZ0BAStatus" , code);
851 default: return 0;
852 }
853}
#define SET_INT(name, val)

◆ setRefTrks() [1/3]

bool xAOD::BPhysHelper::setRefTrks ( )

: Sets refitted track momenta

This method uses as input the tracks-at-vertex branch (std::vector< Trk::VxTrackAtVertex >) in xAOD::Vertex. This method only works in ATHENA and only for transient xAOD::Vertex instances, which still hold the refitted tracks as vxTrackAtVertex branch. Once the vertex class is stored in the xAOD file, the vxTrackAtVertex branch is lost.

This method is meant to be called by the vertex fitters immediately after the vertex is created to convert vxTrackAtVertex branch into the persistifiable form.

Returns
: true on success

Definition at line 341 of file BPhysHelper.cxx.

342{
343 // refitted momenta
344 std::vector<float> px;
345 std::vector<float> py;
346 std::vector<float> pz;
347 const auto N = vtx()->vxTrackAtVertex().size();
348 px.reserve(N);
349 py.reserve(N);
350 pz.reserve(N);
351 // loop over refitted tracks at vertex
352 for(uint i=0; i<N; ++i) {
353 const Trk::TrackParameters* aPerigee = vtx()->vxTrackAtVertex()[i].perigeeAtVertex();
354 //sanity check
355 if(!aPerigee)
356 return false;
357
358 // store individual momentum components
359 px.push_back( aPerigee->momentum()[Trk::px] );
360 py.push_back( aPerigee->momentum()[Trk::py] );
361 pz.push_back( aPerigee->momentum()[Trk::pz] );
362 }
363
364 // store as augmentation:
365 setRefTrks(std::move(px), std::move(py), std::move(pz));
366
367 // all OK
368 return true;
369
370}
const Amg::Vector3D & momentum() const
Access method for the momentum.
bool setRefTrks()
: Sets refitted track momenta
const xAOD::Vertex * vtx() const
Getter method for the cached xAOD::Vertex.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
@ pz
global momentum (cartesian)
Definition ParamDefs.h:61
@ px
Definition ParamDefs.h:59
@ py
Definition ParamDefs.h:60
ParametersBase< TrackParametersDim, Charged > TrackParameters

◆ setRefTrks() [2/3]

bool xAOD::BPhysHelper::setRefTrks ( const std::vector< TVector3 > & refTrks)

Sets refitted track momenta.

Parameters
[in]refTrksstd::vector of refitted momenta as TVector3
Returns
: true on success

Definition at line 312 of file BPhysHelper.cxx.

313{
314 // sanity check
315 if(refTrks.empty())
316 return false;
317
318 // tmp vectors of px,py,pz components:
319 std::vector<float> px;
320 std::vector<float> py;
321 std::vector<float> pz;
322 px.reserve(refTrks.size());
323 py.reserve(refTrks.size());
324 pz.reserve(refTrks.size());
325 // loop over refitted track momenta and store the components
326 std::vector<TVector3>::const_iterator refTrksItr = refTrks.begin();
327 for(; refTrksItr!=refTrks.end(); ++refTrksItr) {
328 px.push_back( refTrksItr->Px() );
329 py.push_back( refTrksItr->Py() );
330 pz.push_back( refTrksItr->Pz() );
331 }
332
333 // call overloaded method:
334 return setRefTrks(std::move(px),std::move(py),std::move(pz));
335
336}
const std::vector< TVector3 > & refTrks()
Returns refitted track momenta.

◆ setRefTrks() [3/3]

bool xAOD::BPhysHelper::setRefTrks ( std::vector< float > px,
std::vector< float > py,
std::vector< float > pz )

Sets refitted track momenta.

Parameters
[in]px,py,pzcomponents of the refitted momenta
Returns
: true on success

Definition at line 286 of file BPhysHelper.cxx.

289{
290 // sanity check:
291 if(px.size()!=py.size() || px.size()!=pz.size())
292 return false;
293
294 // erase refitted track cache to preserve consistency
295 m_refTracksCached = false;
296 m_cachedRefTracks.clear();
297
298 // create decorators
299 static const SG::AuxElement::Decorator< std::vector<float> > refTrackPxDeco("RefTrackPx");
300 static const SG::AuxElement::Decorator< std::vector<float> > refTrackPyDeco("RefTrackPy");
301 static const SG::AuxElement::Decorator< std::vector<float> > refTrackPzDeco("RefTrackPz");
302
303 // store the elements:
304 refTrackPxDeco(*m_b) = std::move(px);
305 refTrackPyDeco(*m_b) = std::move(py);
306 refTrackPzDeco(*m_b) = std::move(pz);
307
308 return true;
309}

◆ setZ0()

float xAOD::BPhysHelper::setZ0 ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

longitudinal impact parameter

Definition at line 1067 of file BPhysHelper.cxx.

1068{
1069 switch(vertexType) {
1070 case PV_MAX_SUM_PT2 : SET_FLOAT("Z0MaxSumPt2", val);
1071 case PV_MIN_A0 : SET_FLOAT("Z0MinA0" , val);
1072 case PV_MIN_Z0 : SET_FLOAT("Z0MinZ0" , val);
1073 case PV_MIN_Z0_BA : SET_FLOAT("Z0MinZ0BA" , val);
1074 default: return false;
1075 }
1076}

◆ setZ0Err()

float xAOD::BPhysHelper::setZ0Err ( const float val,
const pv_type vertexType = BPhysHelper::PV_MIN_A0 )

longitudinal impact parameter error

Definition at line 1078 of file BPhysHelper.cxx.

1079{
1080 switch(vertexType) {
1081 case PV_MAX_SUM_PT2 : SET_FLOAT("Z0ErrMaxSumPt2", val);
1082 case PV_MIN_A0 : SET_FLOAT("Z0ErrMinA0" , val);
1083 case PV_MIN_Z0 : SET_FLOAT("Z0ErrMinZ0" , val);
1084 case PV_MIN_Z0_BA : SET_FLOAT("Z0ErrMinZ0BA" , val);
1085 default: return false;
1086 }
1087}

◆ totalP() [1/2]

TVector3 xAOD::BPhysHelper::totalP ( )

: Returns total 3-momentum calculated from the refitted tracks

Total momentum calculated from refitted tracks

Returns total 3-momentum calculated from the refitted tracks.

Note
In case that this vertex has some linked cascade vertices, it adds the total momentum in the cascade vertices.
Returns
: Total momentum, (0,0,0) on error

Definition at line 374 of file BPhysHelper.cxx.

375{
376 // cache refitted tracks
377 if( !cacheRefTracks() && nCascadeVertices()<=0) {
378 // In case of cascade decays, vector with no refitted tracks can still be correct.
379 // Therefore, return error only is this is not the case
380 return TVector3(0,0,0);
381 }
382
383 // sum the momenta
384 TVector3 sum(0,0,0);
385 int n = nRefTrks();
386 for(int i=0; i<n; ++i) {
387 sum += refTrk(i);
388 }
389
390 // treat cascade vertices, if they exist
391 if(nCascadeVertices()>0) {
392 for(int i=0; i<nCascadeVertices(); ++i) {
393 // create helper class for the daughter vertex
394 xAOD::BPhysHelper daughterVtx( cascadeVertex(i) );
395 TVector3 daughterP = daughterVtx.totalP();
396
397 // check that the daughter total momentum was retrieved correctly
398 if(daughterP==TVector3(0,0,0))
399 return TVector3(0,0,0);
400
401 //OK, sum the momenta
402 sum += daughterP;
403 }
404 }
405
406 return sum;
407}
TVector3 refTrk(const size_t index)
Returns i-th refitted track 3-momentum.

◆ totalP() [2/2]

TLorentzVector xAOD::BPhysHelper::totalP ( std::span< const double > masses)

Returns total 4-momentum calculated from the refitted tracks given mass hypotheses for individual tracks.

Parameters
[in]massesinvariant mass hypotheses
Returns
: Total 4-momentum, (0,0,0,0) on error

Definition at line 410 of file BPhysHelper.cxx.

411{
412 // cache refitted tracks
413 if( !cacheRefTracks() )
414 return TLorentzVector(0,0,0,0);
415 int n = nRefTrks();
416 // input check
417 if(int(masses.size()) != n)
418 return TLorentzVector(0,0,0,0);
419
420 // sum the 4-momenta
421 TLorentzVector sum;
422 for(int i=0; i<n; ++i) {
423 sum += refTrk(i, masses[i]);
424 }
425
426 return sum;
427
428}

◆ vtx()

const xAOD::Vertex * xAOD::BPhysHelper::vtx ( ) const
inline

Getter method for the cached xAOD::Vertex.

Returns
: cached xAOD::Vertex

Definition at line 108 of file BPhysHelper.h.

108{ return m_b; }

◆ z0()

float xAOD::BPhysHelper::z0 ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

longitudinal impact parameter

Definition at line 1001 of file BPhysHelper.cxx.

1002{
1003 switch(vertexType) {
1004 case PV_MAX_SUM_PT2 : GET_FLOAT("Z0MaxSumPt2");
1005 case PV_MIN_A0 : GET_FLOAT("Z0MinA0");
1006 case PV_MIN_Z0 : GET_FLOAT("Z0MinZ0");
1007 case PV_MIN_Z0_BA : GET_FLOAT("Z0MinZ0BA");
1008 default: return -9999999.;
1009 }
1010}

◆ z0Err()

float xAOD::BPhysHelper::z0Err ( const pv_type vertexType = BPhysHelper::PV_MIN_A0)

longitudinal impact parameter error

Definition at line 1012 of file BPhysHelper.cxx.

1013{
1014 switch(vertexType) {
1015 case PV_MAX_SUM_PT2 : GET_FLOAT("Z0ErrMaxSumPt2");
1016 case PV_MIN_A0 : GET_FLOAT("Z0ErrMinA0");
1017 case PV_MIN_Z0 : GET_FLOAT("Z0ErrMinZ0");
1018 case PV_MIN_Z0_BA : GET_FLOAT("Z0ErrMinZ0BA");
1019 default: return -9999999.;
1020 }
1021}

Member Data Documentation

◆ m_b

const xAOD::Vertex* xAOD::BPhysHelper::m_b
protected

Cached B decay xAOD vertex.

Definition at line 674 of file BPhysHelper.h.

◆ m_cachedCascadeVertices

std::vector<const xAOD::Vertex*> xAOD::BPhysHelper::m_cachedCascadeVertices
protected

Definition at line 716 of file BPhysHelper.h.

◆ m_cachedCov

TMatrixTSym<double> xAOD::BPhysHelper::m_cachedCov
protected

Definition at line 682 of file BPhysHelper.h.

◆ m_cachedElectrons

std::vector<const xAOD::Electron*> xAOD::BPhysHelper::m_cachedElectrons
protected

Definition at line 704 of file BPhysHelper.h.

◆ m_cachedMuons

std::vector<const xAOD::Muon*> xAOD::BPhysHelper::m_cachedMuons
protected

Definition at line 703 of file BPhysHelper.h.

◆ m_cachedPrecedingVertices

std::vector<const xAOD::Vertex*> xAOD::BPhysHelper::m_cachedPrecedingVertices
protected

Definition at line 714 of file BPhysHelper.h.

◆ m_cachedRefTracks

std::vector<TVector3> xAOD::BPhysHelper::m_cachedRefTracks
protected

Definition at line 692 of file BPhysHelper.h.

◆ m_cascadeVerticesCached

bool xAOD::BPhysHelper::m_cascadeVerticesCached
protected

Definition at line 715 of file BPhysHelper.h.

◆ m_covCached

bool xAOD::BPhysHelper::m_covCached
protected

cached covariance matrix

Definition at line 681 of file BPhysHelper.h.

◆ m_electronsCached

bool xAOD::BPhysHelper::m_electronsCached
protected

Definition at line 702 of file BPhysHelper.h.

◆ m_muonsCached

bool xAOD::BPhysHelper::m_muonsCached
protected

cached linked muons

Definition at line 701 of file BPhysHelper.h.

◆ m_precedingVerticesCached

bool xAOD::BPhysHelper::m_precedingVerticesCached
protected

cached linked vertices

Definition at line 713 of file BPhysHelper.h.

◆ m_refTracksCached

bool xAOD::BPhysHelper::m_refTracksCached
protected

cached refitted track momenta

Definition at line 691 of file BPhysHelper.h.

◆ n_pv_types

const unsigned int xAOD::BPhysHelper::n_pv_types = 4
static

Return number of vertex association types (useful for loops).

Some useful static consts

Definition at line 613 of file BPhysHelper.h.

◆ pv_type_str

const std::string xAOD::BPhysHelper::pv_type_str
static
Initial value:
=
{"PV_MAX_SUM_PT2", "PV_MIN_A0", "PV_MIN_Z0", "PV_MIN_Z0_BA"}

Return string names for vertex association types.

Definition at line 35 of file BPhysHelper.h.

36 {

◆ s_emptyMatrix

const TMatrixTSym< double > xAOD::BPhysHelper::s_emptyMatrix
staticprotected

Definition at line 728 of file BPhysHelper.h.

◆ s_emptyVectorOfElectrons

const std::vector< const xAOD::Electron * > xAOD::BPhysHelper::s_emptyVectorOfElectrons
staticprotected

Definition at line 727 of file BPhysHelper.h.

◆ s_emptyVectorOfMuons

const std::vector< const xAOD::Muon * > xAOD::BPhysHelper::s_emptyVectorOfMuons
staticprotected

Definition at line 726 of file BPhysHelper.h.

◆ s_emptyVectorOfTVector3

const std::vector< TVector3 > xAOD::BPhysHelper::s_emptyVectorOfTVector3
staticprotected

Return error values for vector quantities:

Definition at line 725 of file BPhysHelper.h.

◆ s_emptyVectorOfVertices

const std::vector< const xAOD::Vertex * > xAOD::BPhysHelper::s_emptyVectorOfVertices
staticprotected

Definition at line 729 of file BPhysHelper.h.


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