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

#include <ReVertex.h>

Inheritance diagram for DerivationFramework::ReVertex:

Public Member Functions

 ReVertex (const std::string &t, const std::string &n, const IInterface *p)
virtual StatusCode initialize () override
virtual StatusCode addBranches (const EventContext &ctx) const override
void fitAndStore (const EventContext &ctx, xAOD::VertexContainer *vtxContainer, const xAOD::Vertex *v, const xAOD::VertexContainer *InVtxContainer, const std::vector< const xAOD::TrackParticle * > &inputTracks, const xAOD::TrackParticleContainer *importedTrackCollection, const xAOD::VertexContainer *pvContainer) const
std::unique_ptr< xAOD::Vertexfit (const EventContext &ctx, const std::vector< const xAOD::TrackParticle * > &inputTracks, const xAOD::TrackParticleContainer *importedTrackCollection, const xAOD::Vertex *pv) const

Private Attributes

std::vector< int > m_TrackIndices
ToolHandle< InDet::VertexPointEstimatorm_vertexEstimator
ToolHandle< Trk::IVertexFitterm_iVertexFitter
Trk::TrkVKalVrtFitterm_VKVFitter {}
SG::WriteHandleKey< xAOD::VertexContainerm_OutputContainerName
SG::ReadHandleKey< xAOD::VertexContainerm_inputContainerName
SG::ReadHandleKey< xAOD::TrackParticleContainerm_trackContainer
SG::WriteHandleKey< xAOD::VertexContainerm_refPVContainerName
SG::ReadHandleKey< xAOD::VertexContainerm_pvContainerName
SG::ReadHandleKey< xAOD::VertexContainerm_defaultPVContainerName
std::vector< double > m_trkMasses
std::vector< int > m_indices
double m_massConst {}
double m_totalMassConst {}
std::vector< std::string > m_hypoNames
PublicToolHandle< Trk::V0Toolsm_v0Tools
PublicToolHandle< Analysis::PrimaryVertexRefitterm_pvRefitter
SG::ReadHandleKey< xAOD::EventInfom_eventInfo_key {this, "EventInfo", "EventInfo", "Input event information"}
int m_PV_max {}
int m_DoVertexType {}
size_t m_PV_minNTracks {}
bool m_do3d {}
bool m_AddPVData {}
bool m_refitPV {}
bool m_doMassConst {}
bool m_startingpoint0 {}
bool m_vertexFittingWithPV {}
double m_BMassUpper {}
double m_BMassLower {}
double m_chi2cut {}
double m_trkDeltaZ {}
bool m_useAdditionalTrack {}
SG::ReadHandleKeyArray< xAOD::VertexContainerm_CollectionsToCheck {this, "CheckVertexContainers", {}}
SG::ReadHandleKeyArray< xAOD::TrackParticleContainerm_RelinkContainers {this, "RelinkTracks", {}, "Track Containers if they need to be relinked through indirect use" }

Detailed Description

Definition at line 42 of file ReVertex.h.

Constructor & Destructor Documentation

◆ ReVertex()

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

Definition at line 27 of file ReVertex.cxx.

29 :
30 base_class(t,n,p), m_vertexEstimator("InDet::VertexPointEstimator"), m_iVertexFitter("Trk::TrkVKalVrtFitter"),
31 m_massConst(0.),
33 m_v0Tools("Trk::V0Tools"),
34 m_pvRefitter("Analysis::PrimaryVertexRefitter"),
35 m_doMassConst(false),
37 m_chi2cut(-1.0),
38 m_trkDeltaZ(-1.0),
40{
41
42 declareProperty("TrackIndices", m_TrackIndices);
43 declareProperty("TrkVertexFitterTool", m_iVertexFitter);
44 declareProperty("VertexPointEstimator",m_vertexEstimator);
45
46 declareProperty("OutputVtxContainerName", m_OutputContainerName);
47 declareProperty("InputVtxContainerName", m_inputContainerName);
48 declareProperty("TrackContainerName", m_trackContainer = "InDetTrackParticles");
49 declareProperty("UseVertexFittingWithPV", m_vertexFittingWithPV);
50
51 declareProperty("HypothesisNames",m_hypoNames);
52
53 declareProperty("V0Tools" , m_v0Tools);
54 declareProperty("PVRefitter" , m_pvRefitter);
55 declareProperty("DefaultPVContainerName", m_defaultPVContainerName = "PrimaryVertices");
56 declareProperty("PVContainerName" , m_pvContainerName = "PrimaryVertices");
57 declareProperty("RefPVContainerName" , m_refPVContainerName = "RefittedPrimaryVertices");
58
59 declareProperty("UseMassConstraint", m_doMassConst);
60 declareProperty("VertexMass", m_totalMassConst);
61 declareProperty("SubVertexMass", m_massConst);
62 declareProperty("MassInputParticles", m_trkMasses);
63 declareProperty("SubVertexTrackIndices", m_indices);
64
65 declareProperty("UseAdditionalTrack", m_useAdditionalTrack);
66
67 declareProperty("RefitPV" , m_refitPV = false);
68 //This parameter will allow us to optimize the number of PVs under consideration as the probability
69 //of a useful primary vertex drops significantly the higher you go
70 declareProperty("MaxPVrefit" , m_PV_max = 1000);
71 declareProperty("DoVertexType" , m_DoVertexType = 7);
72 // minimum number of tracks for PV to be considered for PV association
73 declareProperty("MinNTracksInPV" , m_PV_minNTracks = 0);
74 declareProperty("Do3d" , m_do3d = false);
75 declareProperty("AddPVData" , m_AddPVData = true);
76 declareProperty("StartingPoint0" , m_startingpoint0 = false);
77 declareProperty("BMassUpper",m_BMassUpper = std::numeric_limits<double>::max() );
78 declareProperty("BMassLower",m_BMassLower = std::numeric_limits<double>::min() );
79 declareProperty("Chi2Cut",m_chi2cut = std::numeric_limits<double>::max() );
80 declareProperty("TrkDeltaZ",m_trkDeltaZ);
81
82
83}
PublicToolHandle< Analysis::PrimaryVertexRefitter > m_pvRefitter
Definition ReVertex.h:81
SG::ReadHandleKey< xAOD::VertexContainer > m_defaultPVContainerName
Definition ReVertex.h:72
std::vector< int > m_indices
Definition ReVertex.h:75
ToolHandle< InDet::VertexPointEstimator > m_vertexEstimator
Definition ReVertex.h:64
std::vector< std::string > m_hypoNames
Definition ReVertex.h:78
PublicToolHandle< Trk::V0Tools > m_v0Tools
Definition ReVertex.h:80
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trackContainer
Definition ReVertex.h:69
ToolHandle< Trk::IVertexFitter > m_iVertexFitter
Definition ReVertex.h:65
std::vector< double > m_trkMasses
Definition ReVertex.h:74
SG::WriteHandleKey< xAOD::VertexContainer > m_refPVContainerName
Definition ReVertex.h:70
std::vector< int > m_TrackIndices
Definition ReVertex.h:63
SG::WriteHandleKey< xAOD::VertexContainer > m_OutputContainerName
Definition ReVertex.h:67
SG::ReadHandleKey< xAOD::VertexContainer > m_pvContainerName
Definition ReVertex.h:71
SG::ReadHandleKey< xAOD::VertexContainer > m_inputContainerName
Definition ReVertex.h:68

Member Function Documentation

◆ addBranches()

StatusCode ReVertex::addBranches ( const EventContext & ctx) const
overridevirtual

Definition at line 110 of file ReVertex.cxx.

110 {
111 SG::WriteHandle<xAOD::VertexContainer> vtxContainer(m_OutputContainerName, ctx);
112 ATH_CHECK(vtxContainer.record(std::make_unique<xAOD::VertexContainer>(), std::make_unique<xAOD::VertexAuxContainer>()));
113
114 const size_t Ntracks = m_TrackIndices.size();
115
116 SG::ReadHandle<xAOD::VertexContainer> InVtxContainer(m_inputContainerName, ctx);
117 SG::ReadHandle<xAOD::TrackParticleContainer> importedTrackCollection(m_trackContainer, ctx);
118 ATH_CHECK(InVtxContainer.isValid());
119 ATH_CHECK(importedTrackCollection.isValid());
120 //----------------------------------------------------
121 // retrieve primary vertices
122 //----------------------------------------------------
123 SG::ReadHandle<xAOD::VertexContainer> defaultPVContainer(m_defaultPVContainerName, ctx);
124 ATH_CHECK(defaultPVContainer.isValid());
125
126 SG::ReadHandle<xAOD::VertexContainer> pvContainer(m_pvContainerName, ctx);
127 ATH_CHECK(pvContainer.isValid());
128
129 std::vector<const xAOD::TrackParticle*> fitpair(Ntracks + m_useAdditionalTrack);
130 for(const xAOD::Vertex* v : *InVtxContainer)
131 {
132
133 bool passed = false;
134 for(size_t i=0;i<m_hypoNames.size();i++) {
135 xAOD::BPhysHypoHelper onia(m_hypoNames.at(i), v);
136 passed |= onia.pass();
137 }
138 if (!passed && m_hypoNames.size()) continue;
139
140 for(size_t i =0; i<Ntracks; i++)
141 {
142 size_t trackN = m_TrackIndices[i];
143 if(trackN >= v->nTrackParticles())
144 {
145 ATH_MSG_FATAL("Indices exceeds limit in particle");
146 return StatusCode::FAILURE;
147 }
148 fitpair[i] = v->trackParticle(trackN);
149 }
150
152 {
153 // Loop over ID tracks, call vertexing
154 for (auto trkItr=importedTrackCollection->cbegin(); trkItr!=importedTrackCollection->cend(); ++trkItr) {
155 const xAOD::TrackParticle* tp (*trkItr);
156 fitpair.back() = nullptr;
157 if (Analysis::JpsiUpsilonCommon::isContainedIn(tp,fitpair)) continue; // remove tracks which were used to build J/psi+2Tracks
158 fitpair.back() = tp;
159
160 // Daniel Scheirich: remove track too far from the Jpsi+2Tracks vertex (DeltaZ cut)
161 if(m_trkDeltaZ>0 &&
162 std::abs((tp)->z0() + (tp)->vz() - v->z()) > m_trkDeltaZ )
163 continue;
164
165 fitAndStore(ctx,vtxContainer.ptr(),v,InVtxContainer.cptr(),fitpair,importedTrackCollection.cptr(),pvContainer.cptr());
166 }
167 }
168 else
169 {
170 fitAndStore(ctx,vtxContainer.ptr(),v,InVtxContainer.cptr(),fitpair,importedTrackCollection.cptr(),pvContainer.cptr());
171 }
172 }
173
174 if(m_AddPVData){
175 // Give the helper class the ptr to v0tools and beamSpotsSvc to use
176 SG::ReadHandle<xAOD::EventInfo> evt(m_eventInfo_key, ctx);
177 if(not evt.isValid()) ATH_MSG_ERROR("Cannot Retrieve " << evt.key() );
178 BPhysPVTools helper(&(*m_v0Tools), evt.cptr());
179 helper.SetMinNTracksInPV(m_PV_minNTracks);
180 helper.SetSave3d(m_do3d);
181
182 if(m_refitPV) {
183 //----------------------------------------------------
184 // Try to retrieve refitted primary vertices
185 //----------------------------------------------------
186 SG::WriteHandle<xAOD::VertexContainer> refPvContainer(m_refPVContainerName, ctx);
187 ATH_CHECK(refPvContainer.record(std::make_unique<xAOD::VertexContainer>(), std::make_unique<xAOD::VertexAuxContainer>()));
188
189 if(vtxContainer->size() >0){
190 ATH_CHECK(helper.FillCandwithRefittedVertices(vtxContainer.ptr(), pvContainer.cptr(), refPvContainer.ptr(), &(*m_pvRefitter) , m_PV_max, m_DoVertexType));
191 }
192 }else{
193 if(pvContainer->size()==0) {
194 if(vtxContainer->size() >0) ATH_CHECK(helper.FillCandExistingVertices(vtxContainer.ptr(), defaultPVContainer.cptr(), m_DoVertexType));
195 }
196 else {
197 if(vtxContainer->size() >0) ATH_CHECK(helper.FillCandExistingVertices(vtxContainer.ptr(), pvContainer.cptr(), m_DoVertexType));
198 }
199 }
200 }
201
202 using Analysis::JpsiUpsilonCommon;
203
204 std::vector<const xAOD::TrackParticleContainer*> trackCols;
205 for(const auto &str : m_RelinkContainers){
206 SG::ReadHandle<xAOD::TrackParticleContainer> handle(str,ctx);
207 trackCols.push_back(handle.cptr());
208 }
209 if(not trackCols.empty()){
210 for(xAOD::Vertex* vtx : *vtxContainer){
211 try{
213 }catch(std::runtime_error const& e){
214 ATH_MSG_ERROR(e.what());
215 return StatusCode::FAILURE;
216 }
217 }
218 }
219 return StatusCode::SUCCESS;
220}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
bool passed(DecisionID id, const DecisionIDContainer &)
checks if required decision ID is in the set of IDs in the container
static void RelinkVertexTracks(std::span< const xAOD::TrackParticleContainer *const > trkcols, xAOD::Vertex *vtx)
static bool isContainedIn(const xAOD::TrackParticle *, std::span< const xAOD::TrackParticle *const >) noexcept
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfo_key
Definition ReVertex.h:82
void fitAndStore(const EventContext &ctx, xAOD::VertexContainer *vtxContainer, const xAOD::Vertex *v, const xAOD::VertexContainer *InVtxContainer, const std::vector< const xAOD::TrackParticle * > &inputTracks, const xAOD::TrackParticleContainer *importedTrackCollection, const xAOD::VertexContainer *pvContainer) const
Definition ReVertex.cxx:222
SG::ReadHandleKeyArray< xAOD::TrackParticleContainer > m_RelinkContainers
Definition ReVertex.h:102
virtual void handle(const Incident &inc)
Handle end of run incidents to save the metadata at that point.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ fit()

std::unique_ptr< xAOD::Vertex > ReVertex::fit ( const EventContext & ctx,
const std::vector< const xAOD::TrackParticle * > & inputTracks,
const xAOD::TrackParticleContainer * importedTrackCollection,
const xAOD::Vertex * pv ) const

Definition at line 283 of file ReVertex.cxx.

287{
288 std::unique_ptr<Trk::IVKalState> state = m_VKVFitter->makeState(ctx);
289 if (m_doMassConst && (m_trkMasses.size()==inputTracks.size())) {
290 m_VKVFitter->setMassInputParticles(m_trkMasses, *state);
291 if (m_totalMassConst) m_VKVFitter->setMassForConstraint(m_totalMassConst, *state);
292 if (m_massConst) m_VKVFitter->setMassForConstraint(m_massConst, m_indices, *state);
293 }
294 if (pv) {
295 m_VKVFitter->setCnstType(8, *state);
296 m_VKVFitter->setVertexForConstraint(pv->position().x(),
297 pv->position().y(),
298 pv->position().z(), *state);
299 m_VKVFitter->setCovVrtForConstraint(pv->covariancePosition()(Trk::x,Trk::x),
300 pv->covariancePosition()(Trk::y,Trk::x),
301 pv->covariancePosition()(Trk::y,Trk::y),
302 pv->covariancePosition()(Trk::z,Trk::x),
303 pv->covariancePosition()(Trk::z,Trk::y),
304 pv->covariancePosition()(Trk::z,Trk::z), *state );
305 }
306
307 // Do the fit itself.......
308 // Starting point (use the J/psi position)
309 const Trk::Perigee& aPerigee1 = inputTracks[0]->perigeeParameters();
310 const Trk::Perigee& aPerigee2 = inputTracks[1]->perigeeParameters();
311 int sflag = 0;
312 int errorcode = 0;
313 Amg::Vector3D startingPoint = m_vertexEstimator->getCirclesIntersectionPoint(&aPerigee1,&aPerigee2,sflag,errorcode);
314 if (errorcode != 0) {startingPoint(0) = 0.0; startingPoint(1) = 0.0; startingPoint(2) = 0.0;}
315 std::unique_ptr<xAOD::Vertex> theResult = m_VKVFitter->fit(inputTracks, startingPoint, *state);
316
317 // Added by ASC
318 if(theResult){
319 std::vector<ElementLink<DataVector<xAOD::TrackParticle> > > newLinkVector;
320 for(unsigned int i=0; i< theResult->trackParticleLinks().size(); i++)
321 {
322 ElementLink<DataVector<xAOD::TrackParticle> > mylink=theResult->trackParticleLinks()[i]; //makes a copy (non-const)
323 mylink.setStorableObject( *importedTrackCollection, true);
324 newLinkVector.push_back( mylink );
325 }
326 theResult->clearTracks();
327 theResult->setTrackParticleLinks( newLinkVector );
328 }
329
330 return theResult;
331}
size_t size() const
Number of registered mappings.
Trk::TrkVKalVrtFitter * m_VKVFitter
Definition ReVertex.h:66
Eigen::Matrix< double, 3, 1 > Vector3D
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ x
Definition ParamDefs.h:55
@ z
global position (cartesian)
Definition ParamDefs.h:57
@ y
Definition ParamDefs.h:56

◆ fitAndStore()

void ReVertex::fitAndStore ( const EventContext & ctx,
xAOD::VertexContainer * vtxContainer,
const xAOD::Vertex * v,
const xAOD::VertexContainer * InVtxContainer,
const std::vector< const xAOD::TrackParticle * > & inputTracks,
const xAOD::TrackParticleContainer * importedTrackCollection,
const xAOD::VertexContainer * pvContainer ) const

Definition at line 222 of file ReVertex.cxx.

229{
230 std::unique_ptr<xAOD::Vertex> ptr(fit(ctx, inputTracks, importedTrackCollection, nullptr));
231 if(!ptr)return;
232
233 double chi2DOF = ptr->chiSquared()/ptr->numberDoF();
234 ATH_MSG_DEBUG("Candidate chi2/DOF is " << chi2DOF);
235 bool chi2CutPassed = (m_chi2cut <= 0.0 || chi2DOF < m_chi2cut);
236 if(!chi2CutPassed) { ATH_MSG_DEBUG("Chi Cut failed!"); return; }
237 xAOD::BPhysHelper bHelper(ptr.get());//"get" does not "release" still automatically deleted
238 bHelper.setRefTrks();
239 if (m_trkMasses.size()==inputTracks.size()) {
240 TLorentzVector bMomentum = bHelper.totalP(m_trkMasses);
241 double bMass = bMomentum.M();
242 bool passesCuts = (m_BMassUpper > bMass && bMass > m_BMassLower);
243 if(!passesCuts)return;
244 }
245
246 DerivationFramework::BPhysPVTools::PrepareVertexLinks( ptr.get(), importedTrackCollection );
247 std::vector<const xAOD::Vertex*> thePreceding;
248 thePreceding.push_back(v);
250 //
251 Analysis::CleanUpVertex closestRefPV = Analysis::JpsiUpsilonCommon::ClosestRefPV(bHelper, pvContainer, &(*m_pvRefitter));
252 if (!closestRefPV.get()) return;
253 std::unique_ptr<xAOD::Vertex> ptrPV(fit(ctx, inputTracks, importedTrackCollection, closestRefPV.get()));
254 if(!ptrPV) return;
255
256 double chi2DOFPV = ptrPV->chiSquared()/ptrPV->numberDoF();
257 ATH_MSG_DEBUG("CandidatePV chi2/DOF is " << chi2DOFPV);
258 bool chi2CutPassed = (m_chi2cut <= 0.0 || chi2DOFPV < m_chi2cut);
259 if(!chi2CutPassed) { ATH_MSG_DEBUG("Chi Cut failed!"); return; }
260 xAOD::BPhysHelper bHelperPV(ptrPV.get());//"get" does not "release" still automatically deleted
261 bHelperPV.setRefTrks();
262 if (m_trkMasses.size()==inputTracks.size()) {
263 TLorentzVector bMomentumPV = bHelperPV.totalP(m_trkMasses);
264 double bMass = bMomentumPV.M();
265 bool passesCuts = (m_BMassUpper > bMass && bMass > m_BMassLower);
266 if(!passesCuts)return;
267 }
268
269 bHelperPV.setPrecedingVertices(thePreceding, InVtxContainer);
270 vtxContainer->push_back(ptrPV.release());
271 return; //Don't store other vertex
272 }
273 bHelper.setPrecedingVertices(thePreceding, InVtxContainer);
274 vtxContainer->push_back(ptr.release());
275}
#define ATH_MSG_DEBUG(x)
const xAOD::Vertex * get() const
static Analysis::CleanUpVertex ClosestRefPV(xAOD::BPhysHelper &, const xAOD::VertexContainer *, const Analysis::PrimaryVertexRefitter *)
value_type push_back(value_type pElem)
Add an element to the end of the collection.
static void PrepareVertexLinks(xAOD::Vertex *theResult, const xAOD::TrackParticleContainer *importedTrackCollection)
std::unique_ptr< xAOD::Vertex > fit(const EventContext &ctx, const std::vector< const xAOD::TrackParticle * > &inputTracks, const xAOD::TrackParticleContainer *importedTrackCollection, const xAOD::Vertex *pv) const
Definition ReVertex.cxx:283
void * ptr(T *p)
Definition SGImplSvc.cxx:74

◆ initialize()

StatusCode ReVertex::initialize ( )
overridevirtual

Definition at line 85 of file ReVertex.cxx.

85 {
86 ATH_MSG_DEBUG("in initialize()");
87 if(m_TrackIndices.empty()) {
88 ATH_MSG_FATAL("No Indices provided");
89 return StatusCode::FAILURE;
90 }
91 ATH_CHECK(m_iVertexFitter.retrieve());
92 ATH_CHECK(m_v0Tools.retrieve());
93 ATH_CHECK(m_pvRefitter.retrieve());
94 ATH_CHECK(m_vertexEstimator.retrieve());
95 m_VKVFitter = dynamic_cast<Trk::TrkVKalVrtFitter*>(&(*m_iVertexFitter));
96 if(m_VKVFitter==nullptr) return StatusCode::FAILURE;
97 ATH_CHECK(m_OutputContainerName.initialize());
98 ATH_CHECK(m_inputContainerName.initialize());
99 ATH_CHECK(m_trackContainer.initialize());
101 ATH_CHECK(m_pvContainerName.initialize());
102 ATH_CHECK(m_refPVContainerName.initialize());
103 ATH_CHECK(m_eventInfo_key.initialize());
104 ATH_CHECK(m_RelinkContainers.initialize());
105 ATH_CHECK(m_CollectionsToCheck.initialize());
106 return StatusCode::SUCCESS;
107}
SG::ReadHandleKeyArray< xAOD::VertexContainer > m_CollectionsToCheck
Definition ReVertex.h:101

Member Data Documentation

◆ m_AddPVData

bool DerivationFramework::ReVertex::m_AddPVData {}
private

Definition at line 87 of file ReVertex.h.

87{};

◆ m_BMassLower

double DerivationFramework::ReVertex::m_BMassLower {}
private

Definition at line 95 of file ReVertex.h.

95{};

◆ m_BMassUpper

double DerivationFramework::ReVertex::m_BMassUpper {}
private

Definition at line 94 of file ReVertex.h.

94{};

◆ m_chi2cut

double DerivationFramework::ReVertex::m_chi2cut {}
private

Definition at line 96 of file ReVertex.h.

96{}; // chi2/Ndof of the final veretx

◆ m_CollectionsToCheck

SG::ReadHandleKeyArray<xAOD::VertexContainer> DerivationFramework::ReVertex::m_CollectionsToCheck {this, "CheckVertexContainers", {}}
private

Definition at line 101 of file ReVertex.h.

101{this, "CheckVertexContainers", {}};

◆ m_defaultPVContainerName

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::ReVertex::m_defaultPVContainerName
private

Definition at line 72 of file ReVertex.h.

◆ m_do3d

bool DerivationFramework::ReVertex::m_do3d {}
private

Definition at line 86 of file ReVertex.h.

86{};

◆ m_doMassConst

bool DerivationFramework::ReVertex::m_doMassConst {}
private

Definition at line 89 of file ReVertex.h.

89{};

◆ m_DoVertexType

int DerivationFramework::ReVertex::m_DoVertexType {}
private

Definition at line 84 of file ReVertex.h.

84{};

◆ m_eventInfo_key

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

Definition at line 82 of file ReVertex.h.

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

◆ m_hypoNames

std::vector<std::string> DerivationFramework::ReVertex::m_hypoNames
private

Definition at line 78 of file ReVertex.h.

◆ m_indices

std::vector<int> DerivationFramework::ReVertex::m_indices
private

Definition at line 75 of file ReVertex.h.

◆ m_inputContainerName

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::ReVertex::m_inputContainerName
private

Definition at line 68 of file ReVertex.h.

◆ m_iVertexFitter

ToolHandle< Trk::IVertexFitter > DerivationFramework::ReVertex::m_iVertexFitter
private

Definition at line 65 of file ReVertex.h.

◆ m_massConst

double DerivationFramework::ReVertex::m_massConst {}
private

Definition at line 76 of file ReVertex.h.

76{};

◆ m_OutputContainerName

SG::WriteHandleKey<xAOD::VertexContainer> DerivationFramework::ReVertex::m_OutputContainerName
private

Definition at line 67 of file ReVertex.h.

◆ m_PV_max

int DerivationFramework::ReVertex::m_PV_max {}
private

Definition at line 83 of file ReVertex.h.

83{};

◆ m_PV_minNTracks

size_t DerivationFramework::ReVertex::m_PV_minNTracks {}
private

Definition at line 85 of file ReVertex.h.

85{};

◆ m_pvContainerName

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::ReVertex::m_pvContainerName
private

Definition at line 71 of file ReVertex.h.

◆ m_pvRefitter

PublicToolHandle<Analysis::PrimaryVertexRefitter> DerivationFramework::ReVertex::m_pvRefitter
private

Definition at line 81 of file ReVertex.h.

◆ m_refitPV

bool DerivationFramework::ReVertex::m_refitPV {}
private

Definition at line 88 of file ReVertex.h.

88{};

◆ m_refPVContainerName

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

Definition at line 70 of file ReVertex.h.

◆ m_RelinkContainers

SG::ReadHandleKeyArray<xAOD::TrackParticleContainer> DerivationFramework::ReVertex::m_RelinkContainers {this, "RelinkTracks", {}, "Track Containers if they need to be relinked through indirect use" }
private

Definition at line 102 of file ReVertex.h.

102{this, "RelinkTracks", {}, "Track Containers if they need to be relinked through indirect use" };

◆ m_startingpoint0

bool DerivationFramework::ReVertex::m_startingpoint0 {}
private

Definition at line 90 of file ReVertex.h.

90{};

◆ m_totalMassConst

double DerivationFramework::ReVertex::m_totalMassConst {}
private

Definition at line 77 of file ReVertex.h.

77{};

◆ m_trackContainer

SG::ReadHandleKey<xAOD::TrackParticleContainer> DerivationFramework::ReVertex::m_trackContainer
private

Definition at line 69 of file ReVertex.h.

◆ m_TrackIndices

std::vector<int> DerivationFramework::ReVertex::m_TrackIndices
private

Definition at line 63 of file ReVertex.h.

◆ m_trkDeltaZ

double DerivationFramework::ReVertex::m_trkDeltaZ {}
private

Definition at line 97 of file ReVertex.h.

97{}; // DeltaZ between the JPsi vertex and hadronic tracks Z0

◆ m_trkMasses

std::vector<double> DerivationFramework::ReVertex::m_trkMasses
private

Definition at line 74 of file ReVertex.h.

◆ m_useAdditionalTrack

bool DerivationFramework::ReVertex::m_useAdditionalTrack {}
private

Definition at line 99 of file ReVertex.h.

99{};

◆ m_v0Tools

PublicToolHandle<Trk::V0Tools> DerivationFramework::ReVertex::m_v0Tools
private

Definition at line 80 of file ReVertex.h.

◆ m_vertexEstimator

ToolHandle< InDet::VertexPointEstimator > DerivationFramework::ReVertex::m_vertexEstimator
private

Definition at line 64 of file ReVertex.h.

◆ m_vertexFittingWithPV

bool DerivationFramework::ReVertex::m_vertexFittingWithPV {}
private

Definition at line 92 of file ReVertex.h.

92{};

◆ m_VKVFitter

Trk::TrkVKalVrtFitter* DerivationFramework::ReVertex::m_VKVFitter {}
private

Definition at line 66 of file ReVertex.h.

66{};

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