ATLAS Offline Software
Loading...
Searching...
No Matches
Analysis::JpsiPlus1Track Class Reference

#include <JpsiPlus1Track.h>

Inheritance diagram for Analysis::JpsiPlus1Track:

Public Member Functions

 JpsiPlus1Track (const std::string &t, const std::string &n, const IInterface *p)
 ~JpsiPlus1Track ()
virtual StatusCode initialize () override
virtual StatusCode performSearch (const EventContext &ctx, xAOD::VertexContainer &) const override
xAOD::Vertexfit (const std::vector< const xAOD::TrackParticle * > &, const xAOD::TrackParticleContainer *, const xAOD::TrackParticleContainer *) const
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 ()
static double getInvariantMass (const std::vector< const xAOD::TrackParticle * > &trk, double mass1, double mass2, double mass3)

Protected Member Functions

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.

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

bool m_piMassHyp {}
bool m_kMassHyp {}
double m_trkThresholdPt {}
double m_trkMaxEta {}
double m_BThresholdPt {}
double m_BMassUpper {}
double m_BMassLower {}
SG::ReadHandleKey< xAOD::VertexContainerm_jpsiCollectionKey
double m_jpsiMassUpper {}
double m_jpsiMassLower {}
SG::ReadHandleKey< xAOD::TrackParticleContainerm_TrkParticleCollection
SG::ReadHandleKey< xAOD::MuonContainerm_MuonsUsedInJpsi
bool m_excludeJpsiMuonsOnly {}
bool m_excludeCrossJpsiTracks {}
ToolHandle< Trk::IVertexFitterm_iVertexFitter
ToolHandle< Trk::ITrackSelectorToolm_trkSelector
Trk::TrkVKalVrtFitterm_VKVFitter {}
bool m_useMassConst {}
double m_altMassConst {}
double m_chi2cut {}
double m_trkTrippletMassUpper {}
double m_trkTrippletMassLower {}
double m_trkTrippletPt {}
double m_trkDeltaZ {}
int m_requiredNMuons =0
int m_requiredNElectrons =0
std::vector< double > m_muonMasses
std::vector< int > m_useGSFTrackIndices
SG::ReadHandleKey< xAOD::TrackParticleContainerm_TrkParticleGSFCollection
std::bitset< 3 > m_useGSFTrack
SG::ReadHandleKey< xAOD::ElectronContainerm_electronCollectionKey
bool m_skipNoElectron = false
size_t m_candidateLimit {}
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 40 of file JpsiPlus1Track.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

◆ JpsiPlus1Track()

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

Definition at line 73 of file JpsiPlus1Track.cxx.

73 : AthAlgTool(t,n,p),
74 m_piMassHyp(false),
75 m_kMassHyp(true),
77 m_trkMaxEta(102.5),
78 m_BThresholdPt(0.0),
79 m_BMassUpper(0.0),
80 m_BMassLower(0.0),
81 m_jpsiCollectionKey("JpsiCandidates"),
82 m_jpsiMassUpper(0.0),
83 m_jpsiMassLower(0.0),
84 m_TrkParticleCollection("TrackParticleCandidate"),
88 m_iVertexFitter("Trk::TrkVKalVrtFitter"),
89 m_trkSelector("InDet::TrackSelectorTool"),
90 m_useMassConst(true),
91 m_altMassConst(-1.0),
92 m_chi2cut(-1.0),
95 m_trkTrippletPt(-1.0),
96 m_trkDeltaZ(-1.0),
99 m_candidateLimit(std::numeric_limits<size_t>::max())
100 {
101 declareInterface<JpsiPlus1Track>(this);
102 declareProperty("pionHypothesis",m_piMassHyp);
103 declareProperty("kaonHypothesis",m_kMassHyp);
104 declareProperty("trkThresholdPt",m_trkThresholdPt);
105 declareProperty("trkMaxEta",m_trkMaxEta);
106 declareProperty("BThresholdPt",m_BThresholdPt);
107 declareProperty("BMassUpper",m_BMassUpper);
108 declareProperty("BMassLower",m_BMassLower);
109 declareProperty("JpsiContainerKey",m_jpsiCollectionKey);
110 declareProperty("JpsiMassUpper",m_jpsiMassUpper);
111 declareProperty("JpsiMassLower",m_jpsiMassLower);
112 declareProperty("TrackParticleCollection",m_TrkParticleCollection);
113 declareProperty("MuonsUsedInJpsi",m_MuonsUsedInJpsi);
114 declareProperty("ExcludeJpsiMuonsOnly",m_excludeJpsiMuonsOnly);
115 declareProperty("ExcludeCrossJpsiTracks",m_excludeCrossJpsiTracks); //Essential when trying to make vertices out of multiple muons (set to false)
116 declareProperty("TrkVertexFitterTool",m_iVertexFitter);
117 declareProperty("TrackSelectorTool", m_trkSelector);
118 declareProperty("UseMassConstraint", m_useMassConst);
119 declareProperty("AlternativeMassConstraint",m_altMassConst);
120
121 // additional cuts by Daniel Scheirich copied from 2Tracks by Adam Barton
122 declareProperty("Chi2Cut",m_chi2cut);
123 declareProperty("TrkTrippletMassUpper" ,m_trkTrippletMassUpper);
124 declareProperty("TrkTrippletMassLower" ,m_trkTrippletMassLower);
125 declareProperty("TrkQuadrupletPt" ,m_trkTrippletPt );
126 declareProperty("TrkDeltaZ" ,m_trkDeltaZ );
127 declareProperty("RequireNMuonTracks" ,m_requiredNMuons );
128 declareProperty("RequireNElectronTracks", m_requiredNElectrons );
129 declareProperty("AlternativeMassConstraintTrack", m_muonMasses );
130 declareProperty("UseGSFTrackIndices", m_useGSFTrackIndices );
132 declareProperty("ElectronCollection", m_electronCollectionKey);
133 declareProperty("SkipNoElectron", m_skipNoElectron);
134 declareProperty("CandidateLimit", m_candidateLimit);
135 }
std::vector< double > m_muonMasses
SG::ReadHandleKey< xAOD::MuonContainer > m_MuonsUsedInJpsi
SG::ReadHandleKey< xAOD::VertexContainer > m_jpsiCollectionKey
std::vector< int > m_useGSFTrackIndices
ToolHandle< Trk::ITrackSelectorTool > m_trkSelector
SG::ReadHandleKey< xAOD::ElectronContainer > m_electronCollectionKey
ToolHandle< Trk::IVertexFitter > m_iVertexFitter
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_TrkParticleCollection
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_TrkParticleGSFCollection
AthAlgTool()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ ~JpsiPlus1Track()

Analysis::JpsiPlus1Track::~JpsiPlus1Track ( )

Definition at line 137 of file JpsiPlus1Track.cxx.

137{}

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 }

◆ 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

◆ fit()

xAOD::Vertex * Analysis::JpsiPlus1Track::fit ( const std::vector< const xAOD::TrackParticle * > & inputTracks,
const xAOD::TrackParticleContainer * importedTrackCollection,
const xAOD::TrackParticleContainer * gsfCollection ) const

Definition at line 429 of file JpsiPlus1Track.cxx.

429 {
430
431 std::unique_ptr<Trk::IVKalState> state = m_VKVFitter->makeState();
432
433 // Set the mass constraint if requested by user (default=true)
434 // Can be set by user (m_altMassConstraint) - default is -1.0.
435 // If < 0.0, uses J/psi (default)
436 // If > 0.0, uses the value provided
437 constexpr double jpsiTableMass = ParticleConstants::JpsiMassInMeV;
438
439 if (m_useMassConst) {
440 m_VKVFitter->setMassInputParticles(m_muonMasses, *state);
441 std::array<int,2> indices = {1, 2};
442 if (m_altMassConst<0.0) m_VKVFitter->setMassForConstraint(jpsiTableMass,indices, *state);
443 if (m_altMassConst>0.0) m_VKVFitter->setMassForConstraint(m_altMassConst,indices, *state);
444 }
445
446 // Do the fit itself.......
447 Amg::Vector3D startingPoint(0,0,0);
448 StatusCode sc=m_VKVFitter->VKalVrtFitFast(inputTracks, startingPoint, *state);
449 if(sc.isFailure()){
450 startingPoint = Amg::Vector3D(0,0,0);
451 }
452 xAOD::Vertex* theResult = m_VKVFitter->fit(inputTracks, startingPoint, *state);
453
454 // Added by ASC
455 if(theResult != 0){
456 std::vector<ElementLink<DataVector<xAOD::TrackParticle> > > newLinkVector;
457 for(unsigned int i=0; i< theResult->trackParticleLinks().size(); i++)
458 {
459 ElementLink<DataVector<xAOD::TrackParticle> > mylink=theResult->trackParticleLinks()[i]; //makes a copy (non-const)
460 mylink.setStorableObject( m_useGSFTrack[i] ? *gsfCollection : *importedTrackCollection, true);
461 newLinkVector.push_back( mylink );
462 }
463 theResult->clearTracks();
464 theResult->setTrackParticleLinks( newLinkVector );
465 }
466
467
468 return theResult;
469
470 }
static Double_t sc
std::bitset< 3 > m_useGSFTrack
Trk::TrkVKalVrtFitter * m_VKVFitter
void setTrackParticleLinks(const TrackParticleLinks_t &trackParticles)
Set all track particle links at once.
void clearTracks()
Remove all tracks from the vertex.
const TrackParticleLinks_t & trackParticleLinks() const
Get all the particles associated with the vertex.
Eigen::Matrix< double, 3, 1 > Vector3D
::StatusCode StatusCode
StatusCode definition for legacy code.
std::pair< long int, long int > indices
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ getInvariantMass()

double Analysis::JpsiPlus1Track::getInvariantMass ( const std::vector< const xAOD::TrackParticle * > & trk,
double mass1,
double mass2,
double mass3 )
static

Definition at line 474 of file JpsiPlus1Track.cxx.

476 {
477 const auto trk1V = trk[0]->p4();
478 double px1 = trk1V.Px();
479 double py1 = trk1V.Py();
480 double pz1 = trk1V.Pz();
481 double e1 = sqrt(px1*px1+py1*py1+pz1*pz1+mass1*mass1);
482
483 const auto trk2V = trk[1]->p4();
484 double px2 = trk2V.Px();
485 double py2 = trk2V.Py();
486 double pz2 = trk2V.Pz();
487 double e2 = sqrt(px2*px2+py2*py2+pz2*pz2+mass2*mass2);
488
489 const auto trk3V = trk[2]->p4();
490 double px3 = trk3V.Px();
491 double py3 = trk3V.Py();
492 double pz3 = trk3V.Pz();
493 double e3 = sqrt(px3*px3+py3*py3+pz3*pz3+mass3*mass3);
494
495
496 double pxSum=px1+px2+px3;
497 double pySum=py1+py2+py3;
498 double pzSum=pz1+pz2+pz3;
499 double eSum=e1+e2+e3;
500
501 double M=sqrt((eSum*eSum)-(pxSum*pxSum)-(pySum*pySum)-(pzSum*pzSum));
502
503 return M;
504
505 }
double e2(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 2nd sampling
double e1(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 1st sampling

◆ initialize()

StatusCode Analysis::JpsiPlus1Track::initialize ( )
overridevirtual

Definition at line 38 of file JpsiPlus1Track.cxx.

38 {
39
40 // retrieving vertex Fitter
41 ATH_CHECK(m_iVertexFitter.retrieve());
42 m_VKVFitter = dynamic_cast<Trk::TrkVKalVrtFitter*>(&(*m_iVertexFitter));
43
44 // Get the track selector tool from ToolSvc
45 ATH_CHECK( m_trkSelector.retrieve());
46
47 // Get the vertex point estimator tool from ToolSvc
48 ATH_CHECK(m_jpsiCollectionKey.initialize());
50 if(m_MuonsUsedInJpsi.key() == "NONE") m_MuonsUsedInJpsi = "";//for backwards compatability
51 ATH_CHECK(m_MuonsUsedInJpsi.initialize(!m_MuonsUsedInJpsi.key().empty()));
54
56 ATH_MSG_FATAL("Invalid configuration");
57 return StatusCode::FAILURE;
58 }
59
60 if(m_muonMasses.empty()){
61 m_muonMasses.assign(2, muMass);
62 }
63
64 m_useGSFTrack.reset();
65 for(int i : m_useGSFTrackIndices) m_useGSFTrack.set(i, true);
66 ATH_MSG_DEBUG("Initialize successful");
67
68 return StatusCode::SUCCESS;
69
70 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_DEBUG(x)
constexpr double muMass

◆ 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 & Analysis::JpsiPlus1Track::interfaceID ( )
inlinestatic

Definition at line 47 of file JpsiPlus1Track.h.

47{ return IID_JpsiPlus1Track;};
static const InterfaceID IID_JpsiPlus1Track("JpsiPlus1Track", 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.

◆ performSearch()

StatusCode Analysis::JpsiPlus1Track::performSearch ( const EventContext & ctx,
xAOD::VertexContainer & bContainer ) const
overridevirtual

Implements Analysis::ICandidateSearch.

Definition at line 142 of file JpsiPlus1Track.cxx.

143 {
144 ATH_MSG_DEBUG( "JpsiPlus1Track::performSearch" );
145
146 // Get the J/psis from StoreGate
147 const xAOD::VertexContainer* importedJpsiCollection(0);
148 SG::ReadHandle<xAOD::VertexContainer> jpsiCollectionhandle(m_jpsiCollectionKey,ctx);
149 if(!jpsiCollectionhandle.isValid()){
150 ATH_MSG_ERROR("No VertexContainer with key " << m_jpsiCollectionKey.key() << " found in StoreGate. BCandidates will be EMPTY!");
151 return StatusCode::FAILURE;
152 }else{
153 importedJpsiCollection = jpsiCollectionhandle.cptr();
154 ATH_MSG_DEBUG("Found VxCandidate container with key "<<m_jpsiCollectionKey.key());
155 }
156 ATH_MSG_DEBUG("VxCandidate container size " << importedJpsiCollection->size());
157
158 // Get tracks
159 const xAOD::TrackParticleContainer* importedTrackCollection{nullptr};
160 SG::ReadHandle<xAOD::TrackParticleContainer> TrkParticleHandle(m_TrkParticleCollection,ctx);
161 if(!TrkParticleHandle.isValid()){
162 ATH_MSG_ERROR("No track particle collection with name " << m_TrkParticleCollection << " found in StoreGate!");
163 return StatusCode::FAILURE;
164 } else {
165 importedTrackCollection = TrkParticleHandle.cptr();
166 ATH_MSG_DEBUG("Found track particle collection " << m_TrkParticleCollection << " in StoreGate!");
167 }
168 ATH_MSG_DEBUG("Track container size "<< importedTrackCollection->size());
169
170 const xAOD::TrackParticleContainer* importedGSFTrackCollection(nullptr);
171 if(m_useGSFTrack.any() && !m_TrkParticleGSFCollection.key().empty()){
172 SG::ReadHandle<xAOD::TrackParticleContainer> h(m_TrkParticleGSFCollection,ctx);
173 ATH_CHECK(h.isValid());
174 importedGSFTrackCollection = h.cptr();
175 }
176
177
178 // Get the muon collection used to build the J/psis
179 const xAOD::MuonContainer* importedMuonCollection(0);
180 if (!m_MuonsUsedInJpsi.key().empty()) {
181 SG::ReadHandle<xAOD::MuonContainer> h(m_MuonsUsedInJpsi,ctx);
182 if (!h.isValid()){
183 ATH_MSG_ERROR("No muon collection with name " << m_MuonsUsedInJpsi.key() << " found in StoreGate!");
184 return StatusCode::FAILURE;
185 } else {
186 importedMuonCollection = h.cptr();
187 ATH_MSG_DEBUG("Found muon collection " << m_MuonsUsedInJpsi.key() << " in StoreGate!");
188 }
189 ATH_MSG_DEBUG("Muon container size "<< importedMuonCollection->size());
190 }
191
192 // Get the electrons from StoreGate
193 const xAOD::ElectronContainer* importedElectronCollection = nullptr;
194 if (!m_electronCollectionKey.empty()){
195 SG::ReadHandle<xAOD::ElectronContainer> h(m_electronCollectionKey,ctx);
196 if (!h.isValid()){
197 ATH_MSG_ERROR("No Electron collection with name " << m_electronCollectionKey.key() << " found in StoreGate!");
198 return StatusCode::FAILURE;
199 } else {
200 importedElectronCollection = h.cptr();
201 ATH_MSG_DEBUG("Found Electron collection " << m_electronCollectionKey.key() << " in StoreGate!");
202 }
203 ATH_MSG_DEBUG("Electron container size "<< importedElectronCollection->size());
204 }
205
206
207 // Typedef for vectors of tracks and VxCandidates
208 typedef std::vector<const xAOD::TrackParticle*> TrackBag;
209 typedef std::vector<const xAOD::Electron*> ElectronBag;
210
211 // Select the inner detector tracks
212 const xAOD::Vertex* vx = nullptr;
213 TrackBag theIDTracksAfterSelection;
214 for (auto trkPBItr=importedTrackCollection->cbegin(); trkPBItr!=importedTrackCollection->cend(); ++trkPBItr) {
215 const xAOD::TrackParticle* tp (*trkPBItr);
216 if ( tp->pt()<m_trkThresholdPt ) continue;
217 if ( fabs(tp->eta())>m_trkMaxEta ) continue;
218 if (importedMuonCollection!=NULL && !m_excludeJpsiMuonsOnly) {
219 if (JpsiUpsilonCommon::isContainedIn(tp,importedMuonCollection)) continue;
220 }
221 if ( m_trkSelector->decision(*tp, vx) ) theIDTracksAfterSelection.push_back(tp);
222 }
223
224 // Add GSF tracks
225 if(m_useGSFTrack.any()){
226 ATH_MSG_DEBUG("importedTrackCollection contains GSF tracks " << importedGSFTrackCollection->size());
227 for (auto trkPBItr=importedGSFTrackCollection->cbegin(); trkPBItr!=importedGSFTrackCollection->cend(); ++trkPBItr) {
228 const xAOD::TrackParticle* tp (*trkPBItr);
229 if ( tp->pt()<m_trkThresholdPt ) continue;
230 if ( fabs(tp->eta())>m_trkMaxEta ) continue;
231 if (!importedMuonCollection->empty() && !m_excludeJpsiMuonsOnly) {
232 if (JpsiUpsilonCommon::isContainedIn(tp,importedMuonCollection)) continue;
233 }
234 if (m_trkSelector->decision(*tp, vx)) theIDTracksAfterSelection.push_back(tp);
235 }
236 }
237 if (theIDTracksAfterSelection.empty()) return StatusCode::SUCCESS;
238 ATH_MSG_DEBUG("Number of tracks after ID trkSelector: " << theIDTracksAfterSelection.size());
239
240 // Loop over J/psi candidates, select, collect up tracks used to build a J/psi
241 std::vector<const xAOD::Vertex*> selectedJpsiCandidates;
242 std::vector<const xAOD::TrackParticle*> jpsiTracks;
243 for(auto vxcItr=importedJpsiCollection->cbegin(); vxcItr!=importedJpsiCollection->cend(); ++vxcItr) {
244 // Check J/psi candidate invariant mass and skip if need be
245 if (m_jpsiMassUpper>0.0 || m_jpsiMassLower >0.0) {
246 xAOD::BPhysHelper jpsiCandidate(*vxcItr);
247 double jpsiMass = jpsiCandidate.totalP(m_muonMasses).M();
249 if (!pass) continue;
250 }
251 selectedJpsiCandidates.push_back(*vxcItr);
252
253 // Collect up tracks
255 // Extract tracks from J/psi
256 const xAOD::TrackParticle* jpsiTP1 = (*vxcItr)->trackParticle(0);
257 const xAOD::TrackParticle* jpsiTP2 = (*vxcItr)->trackParticle(1);
258 jpsiTracks.push_back(jpsiTP1);
259 jpsiTracks.push_back(jpsiTP2);
260 }
261 }
262
263
264 std::vector<double> massHypotheses;
265 if (m_kMassHyp) massHypotheses.push_back(kMass);
266 if (m_piMassHyp) massHypotheses.push_back(piMass);
267 if (!m_kMassHyp && !m_piMassHyp && m_BMassUpper>0.0) {
268 massHypotheses.push_back(kMass); massHypotheses.push_back(piMass);
269 }
270 std::vector<double> tripletMasses(m_muonMasses);
271 // Attempt to fit each track with the two tracks from the J/psi candidates
272 // Loop over J/psis
273
274 std::vector<double> massCuts;
275
276 TrackBag muonTracks;
277 if (importedMuonCollection && m_excludeJpsiMuonsOnly) {
278 for(auto muon : *importedMuonCollection){
279 if(!muon->inDetTrackParticleLink().isValid()) continue;
280 auto track = muon->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
281 if(track==nullptr) continue;
282 if(!JpsiUpsilonCommon::isContainedIn(track, theIDTracksAfterSelection)) continue;
283 muonTracks.push_back(track);
284 }
285 }
286
287 TrackBag electronTracks;
288 ElectronBag theElectronsAfterSelection;
289 if (importedElectronCollection && !importedElectronCollection->empty()) {
290 for(auto electron : *importedElectronCollection) {
291 if (!electron->trackParticleLink().isValid()) continue; // No electrons without ID tracks
292 const xAOD::TrackParticle* elTrk(0);
293 elTrk = electron->trackParticleLink().cachedElement();
294 if (!elTrk) continue;
295 theElectronsAfterSelection.push_back(electron);
296 electronTracks.push_back(elTrk);
297 }
298 if (m_skipNoElectron && theElectronsAfterSelection.size() == 0) return StatusCode::SUCCESS;
299 ATH_MSG_DEBUG("Number of electrons after selection: " << theElectronsAfterSelection.size());
300 }
301
302 std::vector<const xAOD::TrackParticle*> tracks(3, nullptr);
303
304 for(auto jpsiItr=selectedJpsiCandidates.cbegin(); jpsiItr!=selectedJpsiCandidates.cend(); ++jpsiItr) {
305
306 // Extract tracks from J/psi
307
308 const xAOD::TrackParticle* jpsiTP1 = tracks[0] = (*jpsiItr)->trackParticle(0);
309 const xAOD::TrackParticle* jpsiTP2 = tracks[1] = (*jpsiItr)->trackParticle(1);
310
311 //If requested, only exclude duplicates in the same triplet
313 jpsiTracks.resize(2);
314 jpsiTracks[0] = jpsiTP1;
315 jpsiTracks[1] = jpsiTP2;
316 }
317
318 // Loop over ID tracks, call vertexing
319 for (auto trkItr=theIDTracksAfterSelection.cbegin(); trkItr!=theIDTracksAfterSelection.cend(); ++trkItr) {
320 if(bContainer.size() >= m_candidateLimit ){
321 ATH_MSG_WARNING("Hard Limit exceeded N=" << bContainer.size());
322 return StatusCode::SUCCESS;
323 }
324 if (!m_excludeJpsiMuonsOnly && JpsiUpsilonCommon::isContainedIn(*trkItr,jpsiTracks)) continue; // remove tracks which were used to build J/psi
325 int linkedMuonTrk = 0;
327 linkedMuonTrk = JpsiUpsilonCommon::isContainedIn(*trkItr, muonTracks);
328 if (linkedMuonTrk) ATH_MSG_DEBUG("This id track is a muon track!");
329 if (JpsiUpsilonCommon::isContainedIn(*trkItr,jpsiTracks)) {
330 if (linkedMuonTrk) ATH_MSG_DEBUG("ID track removed: id track is slected to build Jpsi!");
331 continue;
332
333 }
334 }
335 if(!m_electronCollectionKey.empty()){
336 int linkedElectronTrk = 0;
337 linkedElectronTrk = JpsiUpsilonCommon::isContainedIn(*trkItr, electronTracks);
338 if (linkedElectronTrk) ATH_MSG_DEBUG("This id track is an electron track!");
339 if (JpsiUpsilonCommon::isContainedIn(*trkItr,jpsiTracks)) {
340 if (linkedElectronTrk) ATH_MSG_DEBUG("ID track removed: id track is selected to build Jpsi!");
341 continue;
342 }
343 if (linkedElectronTrk < m_requiredNElectrons) {
344 ATH_MSG_DEBUG("Skipping Tracks with Electron " << linkedElectronTrk << " Limited to " << m_requiredNElectrons);
345 continue;
346 }
347 }
348
349 // Convert to TrackParticleBase
350 const xAOD::TrackParticle* theThirdTP = tracks[2] = *trkItr;
351
352 if (m_trkTrippletPt>0 && JpsiUpsilonCommon::getPt(jpsiTP1, jpsiTP2, theThirdTP) < m_trkTrippletPt ) continue; // track tripplet pT cut (daniel Scheirich)
353 if(m_trkDeltaZ>0 &&
354 fabs(theThirdTP->z0() + theThirdTP->vz() - (*jpsiItr)->z()) > m_trkDeltaZ )
355 continue;
356 if (linkedMuonTrk < m_requiredNMuons) {
357 ATH_MSG_DEBUG("Skipping Tracks with Muons " << linkedMuonTrk << " Limited to " << m_requiredNMuons);
358 continue;
359 }
360
361 // apply mass cut on track tripplet if requested
362 bool passRoughMassCuts(true);
363
365 massCuts.clear();
366 if(m_kMassHyp) massCuts.push_back(getInvariantMass(tracks, m_muonMasses[0], m_muonMasses[1], kMass));
367 if(m_piMassHyp) massCuts.push_back(getInvariantMass(tracks, m_muonMasses[0], m_muonMasses[1], piMass));
369 }
370 if (!passRoughMassCuts) continue;
371
372
373 //Managed pointer, "release" if you don't want it deleted. Automatically "deleted" otherwise
374 std::unique_ptr<xAOD::Vertex> bVertex( fit(tracks, importedTrackCollection, importedGSFTrackCollection));
375 if (bVertex) {
376
377 // Chi2/DOF cut
378 double bChi2DOF = bVertex->chiSquared()/bVertex->numberDoF();
379 ATH_MSG_DEBUG("Candidate chi2/DOF is " << bChi2DOF);
380
381 bool chi2CutPassed = (m_chi2cut <= 0.0 || bChi2DOF < m_chi2cut);
382 if(!chi2CutPassed) { ATH_MSG_DEBUG("Chi Cut failed!"); continue; }
383
384 // create the helper class
385 xAOD::BPhysHelper bHelper(bVertex.get());
386 // setRefTrks needs to be called after BPhysHelper is created if you want to access refitted track parameters
387 bHelper.setRefTrks();
388 // Decide whether to keep the candidate
389 bool masspTpassed = true;
390 if (m_BMassUpper>0.0 || m_BThresholdPt >0.0 || m_BMassLower > 0.0) {
391 masspTpassed = false;
392 for (double masshypo3rd : massHypotheses) {
393 tripletMasses.push_back(masshypo3rd); //Add third mass
394 TLorentzVector bMomentum = bHelper.totalP(tripletMasses);//Calulcate result
395 tripletMasses.pop_back(); //Remove 3rd mass - now same as beginning
396 double bpt = bMomentum.Pt();
397 bool PtPassed = m_BThresholdPt <= 0.0 || (bpt >= m_BThresholdPt);
398 double bMass = bMomentum.M();
399 ATH_MSG_DEBUG("candidate pt/mass under track mass hypothesis of " << masshypo3rd << " is " << bpt << " / " << bMass);
400 bool masscut = JpsiUpsilonCommon::cutRange(bMass, m_BMassLower, m_BMassUpper); //( bMass >= m_BMassLower && bMass <= m_BMassUpper );
401 if(masscut && PtPassed) { masspTpassed = true; break; }
402 }
403 }
404 if (masspTpassed) {
405 // Set links to J/psi
406 std::vector<const xAOD::Vertex*> theJpsiPreceding;
407 theJpsiPreceding.push_back(*jpsiItr);
408 bHelper.setPrecedingVertices(theJpsiPreceding, importedJpsiCollection);
409 bContainer.push_back(bVertex.release());
410 } else {
411 ATH_MSG_DEBUG("Cuts failed!");
412 }
413 } else {ATH_MSG_DEBUG("Fitter failed!");}
414
415 } // End of loop over tracks
416
417 } // End of loop over J/psis
418 ATH_MSG_DEBUG("bContainer size " << bContainer.size());
419 return StatusCode::SUCCESS;
420
421 }
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
xAOD::Vertex * fit(const std::vector< const xAOD::TrackParticle * > &, const xAOD::TrackParticleContainer *, const xAOD::TrackParticleContainer *) const
static double getInvariantMass(const std::vector< const xAOD::TrackParticle * > &trk, double mass1, double mass2, double mass3)
static bool cutRange(double value, double min, double max)
static bool isContainedIn(const xAOD::TrackParticle *, const std::vector< const xAOD::TrackParticle * > &)
static bool cutRangeOR(const std::vector< double > &values, double min, double max)
static double getPt(const xAOD::TrackParticle *, const xAOD::TrackParticle *)
value_type push_back(value_type pElem)
Add an element to the end of the collection.
const_iterator cbegin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
const_iterator cend() const noexcept
Return a const_iterator pointing past the end of the collection.
bool empty() const noexcept
Returns true if the collection is empty.
float z0() const
Returns the parameter.
float vz() const
The z origin for the parameters.
constexpr double piMass
constexpr double kMass
std::vector< const xAOD::TrackParticle * > TrackBag
ElectronContainer_v1 ElectronContainer
Definition of the current "electron container version".
TrackParticle_v1 TrackParticle
Reference the current persistent version:
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".

◆ 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_altMassConst

double Analysis::JpsiPlus1Track::m_altMassConst {}
private

Definition at line 78 of file JpsiPlus1Track.h.

78{};

◆ m_BMassLower

double Analysis::JpsiPlus1Track::m_BMassLower {}
private

Definition at line 66 of file JpsiPlus1Track.h.

66{};

◆ m_BMassUpper

double Analysis::JpsiPlus1Track::m_BMassUpper {}
private

Definition at line 65 of file JpsiPlus1Track.h.

65{};

◆ m_BThresholdPt

double Analysis::JpsiPlus1Track::m_BThresholdPt {}
private

Definition at line 64 of file JpsiPlus1Track.h.

64{};

◆ m_candidateLimit

size_t Analysis::JpsiPlus1Track::m_candidateLimit {}
private

Definition at line 92 of file JpsiPlus1Track.h.

92{};

◆ m_chi2cut

double Analysis::JpsiPlus1Track::m_chi2cut {}
private

Definition at line 79 of file JpsiPlus1Track.h.

79{};

◆ 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_electronCollectionKey

SG::ReadHandleKey<xAOD::ElectronContainer> Analysis::JpsiPlus1Track::m_electronCollectionKey
private

Definition at line 90 of file JpsiPlus1Track.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_excludeCrossJpsiTracks

bool Analysis::JpsiPlus1Track::m_excludeCrossJpsiTracks {}
private

Definition at line 73 of file JpsiPlus1Track.h.

73{}; //Added by Matteo Bedognetti

◆ m_excludeJpsiMuonsOnly

bool Analysis::JpsiPlus1Track::m_excludeJpsiMuonsOnly {}
private

Definition at line 72 of file JpsiPlus1Track.h.

72{};

◆ m_iVertexFitter

ToolHandle< Trk::IVertexFitter > Analysis::JpsiPlus1Track::m_iVertexFitter
private

Definition at line 74 of file JpsiPlus1Track.h.

◆ m_jpsiCollectionKey

SG::ReadHandleKey<xAOD::VertexContainer> Analysis::JpsiPlus1Track::m_jpsiCollectionKey
private

Definition at line 67 of file JpsiPlus1Track.h.

◆ m_jpsiMassLower

double Analysis::JpsiPlus1Track::m_jpsiMassLower {}
private

Definition at line 69 of file JpsiPlus1Track.h.

69{};

◆ m_jpsiMassUpper

double Analysis::JpsiPlus1Track::m_jpsiMassUpper {}
private

Definition at line 68 of file JpsiPlus1Track.h.

68{};

◆ m_kMassHyp

bool Analysis::JpsiPlus1Track::m_kMassHyp {}
private

Definition at line 61 of file JpsiPlus1Track.h.

61{};

◆ m_muonMasses

std::vector<double> Analysis::JpsiPlus1Track::m_muonMasses
private

Definition at line 86 of file JpsiPlus1Track.h.

◆ m_MuonsUsedInJpsi

SG::ReadHandleKey<xAOD::MuonContainer> Analysis::JpsiPlus1Track::m_MuonsUsedInJpsi
private

Definition at line 71 of file JpsiPlus1Track.h.

◆ m_piMassHyp

bool Analysis::JpsiPlus1Track::m_piMassHyp {}
private

Definition at line 60 of file JpsiPlus1Track.h.

60{};

◆ m_requiredNElectrons

int Analysis::JpsiPlus1Track::m_requiredNElectrons =0
private

Definition at line 85 of file JpsiPlus1Track.h.

◆ m_requiredNMuons

int Analysis::JpsiPlus1Track::m_requiredNMuons =0
private

Definition at line 84 of file JpsiPlus1Track.h.

◆ m_skipNoElectron

bool Analysis::JpsiPlus1Track::m_skipNoElectron = false
private

Definition at line 91 of file JpsiPlus1Track.h.

◆ m_trkDeltaZ

double Analysis::JpsiPlus1Track::m_trkDeltaZ {}
private

Definition at line 83 of file JpsiPlus1Track.h.

83{};

◆ m_trkMaxEta

double Analysis::JpsiPlus1Track::m_trkMaxEta {}
private

Definition at line 63 of file JpsiPlus1Track.h.

63{};

◆ m_TrkParticleCollection

SG::ReadHandleKey<xAOD::TrackParticleContainer> Analysis::JpsiPlus1Track::m_TrkParticleCollection
private

Definition at line 70 of file JpsiPlus1Track.h.

◆ m_TrkParticleGSFCollection

SG::ReadHandleKey<xAOD::TrackParticleContainer> Analysis::JpsiPlus1Track::m_TrkParticleGSFCollection
private

Definition at line 88 of file JpsiPlus1Track.h.

◆ m_trkSelector

ToolHandle< Trk::ITrackSelectorTool > Analysis::JpsiPlus1Track::m_trkSelector
private

Definition at line 75 of file JpsiPlus1Track.h.

◆ m_trkThresholdPt

double Analysis::JpsiPlus1Track::m_trkThresholdPt {}
private

Definition at line 62 of file JpsiPlus1Track.h.

62{};

◆ m_trkTrippletMassLower

double Analysis::JpsiPlus1Track::m_trkTrippletMassLower {}
private

Definition at line 81 of file JpsiPlus1Track.h.

81{};

◆ m_trkTrippletMassUpper

double Analysis::JpsiPlus1Track::m_trkTrippletMassUpper {}
private

Definition at line 80 of file JpsiPlus1Track.h.

80{};

◆ m_trkTrippletPt

double Analysis::JpsiPlus1Track::m_trkTrippletPt {}
private

Definition at line 82 of file JpsiPlus1Track.h.

82{};

◆ m_useGSFTrack

std::bitset<3> Analysis::JpsiPlus1Track::m_useGSFTrack
private

Definition at line 89 of file JpsiPlus1Track.h.

◆ m_useGSFTrackIndices

std::vector<int> Analysis::JpsiPlus1Track::m_useGSFTrackIndices
private

Definition at line 87 of file JpsiPlus1Track.h.

◆ m_useMassConst

bool Analysis::JpsiPlus1Track::m_useMassConst {}
private

Definition at line 77 of file JpsiPlus1Track.h.

77{};

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.

◆ m_VKVFitter

Trk::TrkVKalVrtFitter* Analysis::JpsiPlus1Track::m_VKVFitter {}
private

Definition at line 76 of file JpsiPlus1Track.h.

76{};

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