ATLAS Offline Software
Loading...
Searching...
No Matches
InDet::InDetConversionFinderTools Class Referencefinal

#include <InDetConversionFinderTools.h>

Inheritance diagram for InDet::InDetConversionFinderTools:
Collaboration diagram for InDet::InDetConversionFinderTools:

Public Member Functions

 InDetConversionFinderTools (const std::string &t, const std::string &n, const IInterface *p)
 ~InDetConversionFinderTools ()
virtual StatusCode initialize () override
virtual StatusCode finalize () override
virtual std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const EventContext &ctx, const TrackCollection *trk_coll) const override
 Find vertex from Trk::TrackCollection.
virtual std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const EventContext &ctx, const xAOD::TrackParticleContainer *trk_coll) const override
 Conversion candidate reconstruction for Trk::TrackParticle (default)
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Static Public Member Functions

static const InterfaceID & interfaceID ()

Public Attributes

BooleanProperty m_removeTrt {this, "RemoveTrtTracks", false, "Remove standalone TRT tracks"}
BooleanProperty m_isConversion {this, "IsConversion", true, "Conversions or V0s"}
BooleanProperty m_decorateVertices {this, "DecorateVertices", true, "Decorate vertices with values used for vertex selection"}
 Conversion candidate reconstruction for Trk::Tracks.

Protected Member Functions

bool passPreSelection (TrackPairsSelector::Cache &cache, const xAOD::TrackParticle *track_pos, const xAOD::TrackParticle *track_neg, std::vector< Amg::Vector3D > &trackList, Amg::Vector3D &initPos, int &flag, std::map< std::string, float > &intersectionDecors) const
void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Protected Attributes

ToolHandle< Trk::IVertexFitterm_iVertexFitter
 Vertex fitter interface.
ToolHandle< InDet::TrackPairsSelectorm_trackPairsSelector
 Initial conversion vertex estimator tool.
ToolHandle< InDet::VertexPointEstimatorm_vertexEstimator
 Initial conversion vertex estimator tool.
ToolHandle< InDet::ConversionPostSelectorm_postSelector
 Conversion post-fit selector tool.
ToolHandle< InDet::SingleTrackConversionToolm_singleTrkConvTool
 Single track conversion tool.
ToolHandle< Trk::ITrackSelectorToolm_trkSelector
 Track Selector Tool.
DoubleProperty m_mindR
 Cuts.
DoubleProperty m_maxdR
DoubleProperty m_MinInitVtxR
DoubleProperty m_MinFlightAngle

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 48 of file InDetConversionFinderTools.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ InDetConversionFinderTools()

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

Definition at line 25 of file InDetConversionFinderTools.cxx.

28 : AthAlgTool(t, n, p)
29{
30 declareInterface<IVertexFinder>(this);
31}
AthAlgTool()
Default constructor:

◆ ~InDetConversionFinderTools()

InDet::InDetConversionFinderTools::~InDetConversionFinderTools ( )
default

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ finalize()

StatusCode InDet::InDetConversionFinderTools::finalize ( )
overridevirtual

Definition at line 88 of file InDetConversionFinderTools.cxx.

89 {
90 return StatusCode::SUCCESS;
91 }

◆ findVertex() [1/2]

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > InDet::InDetConversionFinderTools::findVertex ( const EventContext & ctx,
const TrackCollection * trackTES ) const
overridevirtual

Find vertex from Trk::TrackCollection.

Parameters
EventContext
inputtrack container
Returns
a pair of newly created container and auxiliary store

Implements InDet::IVertexFinder.

Definition at line 95 of file InDetConversionFinderTools.cxx.

97 {
98
99 ATH_MSG_ERROR("Using Track Container not currently supported returning an empty conatiner");
100
101 // Make collection for conversions.
102 xAOD::VertexContainer* InDetConversionContainer = new xAOD::VertexContainer();
103 xAOD::VertexAuxContainer* InDetConversionContainerAux = new xAOD::VertexAuxContainer();
104 InDetConversionContainer->setStore( InDetConversionContainerAux );
105
106 return std::make_pair(InDetConversionContainer,InDetConversionContainerAux);
107 }
#define ATH_MSG_ERROR(x)
VertexAuxContainer_v1 VertexAuxContainer
Definition of the current jet auxiliary container.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".

◆ findVertex() [2/2]

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > InDet::InDetConversionFinderTools::findVertex ( const EventContext & ctx,
const xAOD::TrackParticleContainer * trk_coll ) const
overridevirtual

Conversion candidate reconstruction for Trk::TrackParticle (default)

Implements InDet::IVertexFinder.

Definition at line 111 of file InDetConversionFinderTools.cxx.

113 {
114 // Make collection for conversions.
115 xAOD::VertexContainer* InDetConversionContainer =
117 xAOD::VertexAuxContainer* InDetConversionContainerAux =
119 InDetConversionContainer->setStore(InDetConversionContainerAux);
120 // Count for number of successful conversions
121 int numConversions = 0;
122
123 // have to be used for the vertex fit
124 Amg::Vector3D pos(0., 0., 0.);
125 Amg::Vector3D initPos(0., 0., 0.);
126
127 // Make seperate lists of positive and negative tracks (after presection
128 // cuts)
129 std::vector<const xAOD::TrackParticle*> negSelectedTracks;
130 negSelectedTracks.clear();
131 std::vector<const xAOD::TrackParticle*> posSelectedTracks;
132 posSelectedTracks.clear();
133 std::vector<int> negIndx;
134 std::vector<int> posIndx;
135
136 TrackPairsSelector::Cache cache{};
137 // track preselection -->pt-cut
139 for (iter = (*trk_coll).begin(); iter != (*trk_coll).end(); ++iter) {
140 if (m_trkSelector->decision(
141 **iter,
142 nullptr)) { // Only save if track passes the pt, d0, z0 and TRT PID cuts
143 ATH_MSG_DEBUG("Track passed preselection");
144 if ((*iter)->charge() < 0) {
145 negSelectedTracks.push_back(*iter);
146 negIndx.push_back(0);
147 } else {
148 posSelectedTracks.push_back(*iter);
149 posIndx.push_back(0);
150 }
151 } else
152 ATH_MSG_DEBUG("Track failed preselection");
153 } // end pt,d0.z0-cuts
154
155 // Make track pairs. To be used for double leg conversions
156 std::vector<const xAOD::TrackParticle*>::const_iterator iter_pos;
157 std::vector<const xAOD::TrackParticle*>::const_iterator iter_neg;
158 std::vector<Amg::Vector3D> positionList;
159 positionList.clear();
160 std::vector<const xAOD::TrackParticle*> trackParticleList;
161 trackParticleList.clear();
162 std::vector<const xAOD::TrackParticle*> singleTrackConvList;
163 singleTrackConvList.clear();
164 int ipos = -1;
165 int ineg = -1;
166 // Outer loop: Loop over positive tracks
167 for (iter_pos = posSelectedTracks.begin();
168 iter_pos != posSelectedTracks.end();
169 ++iter_pos) {
170 ineg = -1;
171 ipos++;
172 // Inner loop: Loop over negative tracks
173 for (iter_neg = negSelectedTracks.begin();
174 iter_neg != negSelectedTracks.end();
175 ++iter_neg) {
176 ineg++;
177 int flag = 0;
178
179 std::map<std::string, float> intersectionDecors;
180 if (!passPreSelection(cache,
181 *iter_pos,
182 *iter_neg,
183 positionList,
184 initPos,
185 flag,
186 intersectionDecors)) {
187 positionList.clear();
188 continue;
189 }
190
191 // Do the fit
192 if (positionList.size() < 2) {
193 ATH_MSG_DEBUG("No tracks to fit ");
194 positionList.clear();
195 continue;
196 }
197
198 trackParticleList.push_back(*iter_pos);
199 trackParticleList.push_back(*iter_neg);
200
201 std::unique_ptr<xAOD::Vertex> myVertex =
202 m_iVertexFitter->fit(ctx, trackParticleList, initPos);
203 trackParticleList.clear();
204
205 // We have a new vertex
206 if (myVertex) {
207 ATH_MSG_DEBUG("VertexFit successful!");
208 int type = -1;
209 if ((m_isConversion && m_postSelector->selectConversionCandidate(
210 myVertex.get(), flag, positionList)) ||
211 (!m_isConversion &&
212 m_postSelector->selectSecVtxCandidate(
213 myVertex.get(), flag, positionList, type))) {
214
215 ATH_MSG_DEBUG(" Conversion passed postselection cuts");
216 // Remove old element links
217 myVertex->clearTracks();
218
219 // If we do not have a valid type just reset
220 if (!m_isConversion && !(type == 101) && !(type == 110) &&
221 !(type == 11)) {
222 myVertex.reset();
223 }
224
225 // If we have the right type (not reset above) lets fill information
226 // and then push to the containers
227 if (myVertex) {
228 if (m_decorateVertices) {
230 "Decorating vertex with values used in track pair selector");
231 for (const auto& kv :
233 SG::Accessor<float> acc (kv.first);
234 acc (*myVertex) = kv.second;
235 }
236 ATH_MSG_DEBUG("Decorating vertex with values used in vertex "
237 "point estimator");
238 for (const auto& kv : intersectionDecors) {
239 SG::Accessor<float> acc (kv.first);
240 acc (*myVertex) = kv.second;
241 }
242 }
243
244 ElementLink<xAOD::TrackParticleContainer> newLinkPos(*iter_pos,
245 *trk_coll);
246 ElementLink<xAOD::TrackParticleContainer> newLinkNeg(*iter_neg,
247 *trk_coll);
248 myVertex->addTrackAtVertex(newLinkPos);
249 myVertex->addTrackAtVertex(newLinkNeg);
250
251 // Now fill in the containers depending on the 2 possible
252 // cases
253 if (m_isConversion) {
254 myVertex->setVertexType(xAOD::VxType::ConvVtx);
255 InDetConversionContainer->push_back(std::move(myVertex));
256 }
257 else if (type == 101 || type == 110 || type == 11) { // V0
258 myVertex->setVertexType(xAOD::VxType::V0Vtx);
259 InDetConversionContainer->push_back(std::move(myVertex));
260 }
261
262 } // End if on right type
263
264 negIndx[ineg] = 1;
265 posIndx[ipos] = 1;
266 ++numConversions;
267 } else {
268 ATH_MSG_DEBUG("VxCandidate failed the post selection cuts!");
269 myVertex.reset();
270 }
271 } else {
272 ATH_MSG_DEBUG("VertexFit was NOT successful!");
273 }
274 positionList.clear();
275 } // neg loop
276 } // pos loop
277 ATH_MSG_DEBUG("Number of conversions found passing post selection cuts: "
278 << numConversions);
279
280 // single track conversions
281 if (m_isConversion) {
282 for (int ip = 0; ip < int(posIndx.size()); ++ip) {
283 if (posIndx[ip] == 0)
284 singleTrackConvList.push_back(posSelectedTracks[ip]);
285 }
286 for (int in = 0; in < int(negIndx.size()); ++in) {
287 if (negIndx[in] == 0)
288 singleTrackConvList.push_back(negSelectedTracks[in]);
289 }
290
291 std::vector<const xAOD::TrackParticle*>::iterator itk;
292 std::vector<const xAOD::TrackParticle*>::iterator itke = singleTrackConvList.end();
293 int numSingle = 0;
294 for (itk = singleTrackConvList.begin(); itk != itke; ++itk) {
295 if (!m_singleTrkConvTool->selectSingleTrackParticleConversion((*itk)))
296 ATH_MSG_DEBUG("Track failed single track conversion selection");
297 else {
298 xAOD::Vertex* sConver(nullptr);
299 sConver = m_singleTrkConvTool->buildSingleTrackParticleConversion(
300 (*itk), InDetConversionContainer);
301 if (sConver) {
302 sConver->clearTracks();
303
304 ElementLink<xAOD::TrackParticleContainer> newLink;
305 newLink.toContainedElement(*trk_coll, *itk);
306 sConver->addTrackAtVertex(newLink);
307 sConver->setVertexType(xAOD::VxType::ConvVtx);
308 numSingle++;
309
310 if (m_decorateVertices) {
311 ATH_MSG_DEBUG("Decorating single track vertex with dummy values "
312 "used in track pair selector");
313 for (const auto& kv : InDet::TrackPairsSelector::getLastValues(cache)) {
314 SG::Accessor<float> acc (kv.first);
315 acc (*sConver) = 0.;
316 }
317
318 ATH_MSG_DEBUG("Decorating single track vertex with dummy values "
319 "used in vertex point estimator");
320 for (const std::string& k : InDet::VertexPointEstimator::decorKeys()) {
321 SG::Accessor<float> acc (k);
322 acc (*sConver) = 0.;
323 }
324
325 ATH_MSG_DEBUG("Decorating single track vertex with dummy values "
326 "used in post selector");
327 InDet::ConversionPostSelector::decorateVertex(*sConver, 0., 0., 0., 0., 0.);
328 }
329 }
330 }
331 }
332 ATH_MSG_DEBUG("Number successful reconstructed single track conversion: "
333 << numSingle);
334 }
335
336 return std::make_pair(InDetConversionContainer,
337 InDetConversionContainerAux);
338 }
#define ATH_MSG_DEBUG(x)
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
value_type push_back(value_type pElem)
Add an element to the end of the collection.
static void decorateVertex(xAOD::Vertex &vertex, float inv_mass, float pt1, float pt2, float fR, float deltaPhiVtxTrk)
Decorate vertices with values used in post selector.
ToolHandle< Trk::ITrackSelectorTool > m_trkSelector
Track Selector Tool.
BooleanProperty m_decorateVertices
Conversion candidate reconstruction for Trk::Tracks.
bool passPreSelection(TrackPairsSelector::Cache &cache, const xAOD::TrackParticle *track_pos, const xAOD::TrackParticle *track_neg, std::vector< Amg::Vector3D > &trackList, Amg::Vector3D &initPos, int &flag, std::map< std::string, float > &intersectionDecors) const
ToolHandle< Trk::IVertexFitter > m_iVertexFitter
Vertex fitter interface.
ToolHandle< InDet::ConversionPostSelector > m_postSelector
Conversion post-fit selector tool.
ToolHandle< InDet::SingleTrackConversionTool > m_singleTrkConvTool
Single track conversion tool.
static std::map< std::string, float > getLastValues(const Cache &cache)
Return a map with the values calculated for the last pair to decorate the vertex once it is created.
static std::vector< std::string > decorKeys()
Return list of keys used for decorations.
Eigen::Matrix< double, 3, 1 > Vector3D
bool flag
Definition master.py:29
@ V0Vtx
Vertex from V0 decay.
@ ConvVtx
Conversion vertex.
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ initialize()

StatusCode InDet::InDetConversionFinderTools::initialize ( )
overridevirtual

Definition at line 35 of file InDetConversionFinderTools.cxx.

36 {
37 StatusCode sc = AthAlgTool::initialize();
38 if (sc.isFailure()) {
39 ATH_MSG_FATAL("Unable to initialize InDetConversionFinderTools");
40 return StatusCode::FAILURE;
41 }
42 /* Get the right vertex fitting tool from ToolSvc */
43 if (m_iVertexFitter.retrieve().isFailure()) {
44 ATH_MSG_FATAL("Failed to retrieve tool " << m_iVertexFitter);
45 return StatusCode::FAILURE;
46 }
47 ATH_MSG_DEBUG("Retrieved tool " << m_iVertexFitter);
48
49 /* Get the track pairs selector tool from ToolSvc */
50 if (m_trackPairsSelector.retrieve().isFailure()) {
51 ATH_MSG_FATAL("Failed to retrieve tool " << m_trackPairsSelector);
52 return StatusCode::FAILURE;
53 }
54 ATH_MSG_DEBUG("Retrieved tool " << m_trackPairsSelector);
55
56 /* Get the vertex point estimator tool from ToolSvc */
57 if (m_vertexEstimator.retrieve().isFailure()) {
58 ATH_MSG_FATAL("Failed to retrieve tool " << m_vertexEstimator);
59 return StatusCode::FAILURE;
60 }
61 ATH_MSG_DEBUG("Retrieved tool " << m_vertexEstimator);
62
63 /* Get the postselector tool from ToolSvc */
64 if (m_postSelector.retrieve().isFailure()) {
65 ATH_MSG_FATAL("Failed to retrieve tool " << m_postSelector);
66 return StatusCode::FAILURE;
67 }
68 ATH_MSG_DEBUG("Retrieved tool " << m_postSelector);
69
70 /* Get the single track conversion tool from ToolSvc */
71 if (m_singleTrkConvTool.retrieve().isFailure()) {
72 ATH_MSG_FATAL("Failed to retrieve tool " << m_singleTrkConvTool);
73 return StatusCode::FAILURE;
74 }
75 ATH_MSG_DEBUG("Retrieved tool " << m_singleTrkConvTool);
76
77 /* Get the track selector tool from ToolSvc */
78 if (m_trkSelector.retrieve().isFailure()) {
79 ATH_MSG_FATAL("Failed to retrieve tool " << m_trkSelector);
80 return StatusCode::FAILURE;
81 }
82 ATH_MSG_DEBUG("Retrieved tool " << m_trkSelector);
83
84 ATH_MSG_DEBUG("Initialization successful");
85 return StatusCode::SUCCESS;
86 }
#define ATH_MSG_FATAL(x)
static Double_t sc
ToolHandle< InDet::VertexPointEstimator > m_vertexEstimator
Initial conversion vertex estimator tool.
ToolHandle< InDet::TrackPairsSelector > m_trackPairsSelector
Initial conversion vertex estimator tool.
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ interfaceID()

const InterfaceID & InDet::IVertexFinder::interfaceID ( )
inlinestaticinherited

Definition at line 56 of file IVertexFinder.h.

56{ return IID_IVertexFinder; }
static const InterfaceID IID_IVertexFinder("IVertexFinder", 1, 0)

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ passPreSelection()

bool InDet::InDetConversionFinderTools::passPreSelection ( TrackPairsSelector::Cache & cache,
const xAOD::TrackParticle * track_pos,
const xAOD::TrackParticle * track_neg,
std::vector< Amg::Vector3D > & trackList,
Amg::Vector3D & initPos,
int & flag,
std::map< std::string, float > & intersectionDecors ) const
protected

Definition at line 341 of file InDetConversionFinderTools.cxx.

349 {
350 //Track summary information
351
352 uint8_t nclusPos(0);
353 uint8_t dummy(0);
354 if(track_pos->summaryValue(dummy,xAOD::numberOfSCTHits)){
355 nclusPos += dummy;
356 }
357 if(track_pos->summaryValue(dummy,xAOD::numberOfPixelHits)){
358 nclusPos += dummy;
359 }
360
361 uint8_t nclusNeg(0);
362 if(track_neg->summaryValue(dummy,xAOD::numberOfSCTHits)){
363 nclusNeg += dummy;
364 }
365 if(track_neg->summaryValue(dummy,xAOD::numberOfPixelHits)){
366 nclusNeg += dummy;
367 }
368
369
370 if(nclusNeg>0 && nclusPos>0) flag= 0;
371 if((nclusNeg>0 && nclusPos==0) || (nclusNeg==0 && nclusPos>0)) flag = 1;
372 if(nclusNeg==0 && nclusPos==0) flag = 2;
373 if(m_removeTrt && (flag==1 || flag==2)) return false;
374
375 if (m_trackPairsSelector->selectTrackParticlePair( track_pos,track_neg,cache)){
376
377 const Trk::Perigee& perPos = track_pos->perigeeParameters();
378 const Trk::Perigee& perNeg = track_neg->perigeeParameters();
379 int errorcode = 0;
380 Amg::Vector3D startingPoint(
381 m_vertexEstimator->getCirclesIntersectionPoint(
382 &perPos, &perNeg, flag, errorcode, intersectionDecors));
383 if(m_isConversion && errorcode != 0) return false;
384 if(!m_isConversion){
385 Amg::Vector3D v_direction = perPos.momentum() + perNeg.momentum();
386 double d_alpha = (startingPoint.adjoint()*v_direction)[0]/(startingPoint.mag()*v_direction.mag());
387 if(d_alpha<m_MinFlightAngle) return false;
388 }
389 initPos = startingPoint;
390 double startingPointR = startingPoint.perp();
391 if(startingPointR<800.) {
392 //Get measured track parameters at first track measurement.
393
394 //Position of first hit in track particle
395 Trk::CurvilinearParameters trkPar_pos;
396 Trk::CurvilinearParameters trkPar_neg;
397
398 int index(-1);
399 for(unsigned int i(0); i< track_pos->numberOfParameters() ; ++i ){
400 if( xAOD::FirstMeasurement == track_pos->parameterPosition(i) ){
401 index = i;
402 break;
403 }
404 }
405 if(index!=-1){
406 trkPar_pos = track_pos->curvilinearParameters(index);
407 } else {
408 ATH_MSG_WARNING("Track Particle does not contain first Measurement track parameters");
409 return false;
410 }
411
412 index = -1;
413 for(unsigned int i(0); i< track_neg->numberOfParameters() ; ++i ){
414 if( xAOD::FirstMeasurement == track_neg->parameterPosition(i) ){
415 index = i;
416 break;
417 }
418 }
419 if(index!=-1){
420 trkPar_neg = track_neg->curvilinearParameters(index);
421 } else {
422 ATH_MSG_WARNING("Track Particle does not contain first Measurement track parameters");
423 return false;
424 }
425
426
427 if(!m_isConversion){
428 double posR = trkPar_pos.position().perp();
429 double negR = trkPar_neg.position().perp();
430 double diffR = 1000.;
431 if((startingPointR-posR)<(startingPointR-negR)) diffR = startingPointR-posR;
432 else diffR = startingPointR-negR;
433 if(startingPointR<m_MinInitVtxR) return false;
434 if(diffR<m_mindR || diffR>m_maxdR) return false;
435 }
436
437 //Alternative way to estimate the guess vertex
438 initPos = trkPar_pos.position() + trkPar_neg.position();
439 initPos *= 0.5;
440 // Not correct but will do for now
441 trackList.push_back(trkPar_pos.position());
442 trackList.push_back(trkPar_neg.position());
443
444 }
445 }
446 return true;
447 }
#define ATH_MSG_WARNING(x)
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
const Trk::CurvilinearParameters curvilinearParameters(unsigned int index) const
Returns a curvilinear representation of the parameters at 'index'.
xAOD::ParameterPosition parameterPosition(unsigned int index) const
Return the ParameterPosition of the parameters at 'index'.
bool summaryValue(uint8_t &value, const SummaryType &information) const
Accessor for TrackSummary values.
size_t numberOfParameters() const
Returns the number of additional parameters stored in the TrackParticle.
str index
Definition DeMoScan.py:362
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
CurvilinearParametersT< TrackParametersDim, Charged, PlaneSurface > CurvilinearParameters
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ FirstMeasurement
Parameter defined at the position of the 1st measurement.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_decorateVertices

BooleanProperty InDet::InDetConversionFinderTools::m_decorateVertices {this, "DecorateVertices", true, "Decorate vertices with values used for vertex selection"}

Conversion candidate reconstruction for Trk::Tracks.

Definition at line 65 of file InDetConversionFinderTools.h.

66{this, "DecorateVertices", true, "Decorate vertices with values used for vertex selection"};

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_isConversion

BooleanProperty InDet::InDetConversionFinderTools::m_isConversion {this, "IsConversion", true, "Conversions or V0s"}

Definition at line 63 of file InDetConversionFinderTools.h.

64{this, "IsConversion", true, "Conversions or V0s"};

◆ m_iVertexFitter

ToolHandle<Trk::IVertexFitter> InDet::InDetConversionFinderTools::m_iVertexFitter
protected
Initial value:
{ this,
"VertexFitterTool",
"Trk::TrkVKalVrtFitter",
"Vertex fitter Tool" }

Vertex fitter interface.

Track pair selection tool.

Definition at line 82 of file InDetConversionFinderTools.h.

82 { this,
83 "VertexFitterTool",
84 "Trk::TrkVKalVrtFitter",
85 "Vertex fitter Tool" };

◆ m_maxdR

DoubleProperty InDet::InDetConversionFinderTools::m_maxdR
protected
Initial value:
{this, "MaxDistVtxHit", 250.,
"Maximum allowed radial distance between guess vertex and closest 1st hit of participating track"}

Definition at line 134 of file InDetConversionFinderTools.h.

134 {this, "MaxDistVtxHit", 250.,
135 "Maximum allowed radial distance between guess vertex and closest 1st hit of participating track"};

◆ m_mindR

DoubleProperty InDet::InDetConversionFinderTools::m_mindR
protected
Initial value:
{this, "MinDistVtxHit", -350.,
"Minimum allowed radial distance beteeen guess vertex and closest 1st hit of participating track"}

Cuts.

Definition at line 132 of file InDetConversionFinderTools.h.

132 {this, "MinDistVtxHit", -350.,
133 "Minimum allowed radial distance beteeen guess vertex and closest 1st hit of participating track"};

◆ m_MinFlightAngle

DoubleProperty InDet::InDetConversionFinderTools::m_MinFlightAngle
protected
Initial value:
{this, "MinFlightAngle", 0.,
"Minimum allowed angular difference between V0 and children direction. Used only in V0 reconstruction."}

Definition at line 138 of file InDetConversionFinderTools.h.

138 {this, "MinFlightAngle", 0.,
139 "Minimum allowed angular difference between V0 and children direction. Used only in V0 reconstruction."};

◆ m_MinInitVtxR

DoubleProperty InDet::InDetConversionFinderTools::m_MinInitVtxR
protected
Initial value:
{this, "MinInitVtxR", 0.,
"Minimum allowed radial position for initial guess vertex. Used only in V0 reconstruction."}

Definition at line 136 of file InDetConversionFinderTools.h.

136 {this, "MinInitVtxR", 0.,
137 "Minimum allowed radial position for initial guess vertex. Used only in V0 reconstruction."};

◆ m_postSelector

ToolHandle<InDet::ConversionPostSelector> InDet::InDetConversionFinderTools::m_postSelector
protected
Initial value:
{
this,
"PostSelector",
"InDet::ConversionPostSelector",
"Tool for post selection of conversion candidates"
}

Conversion post-fit selector tool.

Definition at line 101 of file InDetConversionFinderTools.h.

101 {
102 this,
103 "PostSelector",
104 "InDet::ConversionPostSelector",
105 "Tool for post selection of conversion candidates"
106 };

◆ m_removeTrt

BooleanProperty InDet::InDetConversionFinderTools::m_removeTrt {this, "RemoveTrtTracks", false, "Remove standalone TRT tracks"}

Definition at line 61 of file InDetConversionFinderTools.h.

62{this, "RemoveTrtTracks", false, "Remove standalone TRT tracks"};

◆ m_singleTrkConvTool

ToolHandle<InDet::SingleTrackConversionTool> InDet::InDetConversionFinderTools::m_singleTrkConvTool
protected
Initial value:
{
this,
"SingleTrackConversionTool",
"InDet::SingleTrackConversionTool",
"Tool for single track conversions"
}

Single track conversion tool.

Definition at line 108 of file InDetConversionFinderTools.h.

108 {
109 this,
110 "SingleTrackConversionTool",
111 "InDet::SingleTrackConversionTool",
112 "Tool for single track conversions"
113 };

◆ m_trackPairsSelector

ToolHandle<InDet::TrackPairsSelector> InDet::InDetConversionFinderTools::m_trackPairsSelector
protected
Initial value:
{
this,
"TrackPairsSelector",
"InDet::TrackPairsSelector",
"Track Pair Selector Tool"
}

Initial conversion vertex estimator tool.

Definition at line 87 of file InDetConversionFinderTools.h.

87 {
88 this,
89 "TrackPairsSelector",
90 "InDet::TrackPairsSelector",
91 "Track Pair Selector Tool"
92 };

◆ m_trkSelector

ToolHandle<Trk::ITrackSelectorTool> InDet::InDetConversionFinderTools::m_trkSelector
protected
Initial value:
{
this,
"TrackSelectorTool",
"InDet::TrackSelectorTool",
"Tool for track Selection"
}

Track Selector Tool.

Definition at line 116 of file InDetConversionFinderTools.h.

116 {
117 this,
118 "TrackSelectorTool",
119 "InDet::TrackSelectorTool",
120 "Tool for track Selection"
121 };

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vertexEstimator

ToolHandle<InDet::VertexPointEstimator> InDet::InDetConversionFinderTools::m_vertexEstimator
protected
Initial value:
{
this,
"VertexPointEstimator",
"InDet::VertexPointEstimator",
"Vertex point estimator"
}

Initial conversion vertex estimator tool.

Definition at line 94 of file InDetConversionFinderTools.h.

94 {
95 this,
96 "VertexPointEstimator",
97 "InDet::VertexPointEstimator",
98 "Vertex point estimator"
99 };

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


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