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

#include <JpsiFinder.h>

Inheritance diagram for Analysis::JpsiFinder:
Collaboration diagram for Analysis::JpsiFinder:

Public Member Functions

 JpsiFinder (const std::string &t, const std::string &n, const IInterface *p)
 ~JpsiFinder ()
virtual StatusCode initialize () override
virtual StatusCode performSearch (const EventContext &ctx, xAOD::VertexContainer &vxContainer) const override
std::vector< JpsiCandidategetPairs (const std::vector< const xAOD::TrackParticle * > &) const
std::vector< JpsiCandidategetPairs (const std::vector< const xAOD::Muon * > &) const
std::vector< JpsiCandidategetPairs2Colls (const std::vector< const xAOD::TrackParticle * > &, const std::vector< const xAOD::Muon * > &, bool) const
double getInvariantMass (const JpsiCandidate &, std::span< const double >) const
std::vector< JpsiCandidateselectCharges (const std::vector< JpsiCandidate > &) const
xAOD::Vertexfit (const std::vector< const xAOD::TrackParticle * > &, const xAOD::TrackParticleContainer *importedTrackCollection) const
bool passesMCPCuts (const xAOD::Muon *) const
bool isContainedIn (const xAOD::TrackParticle *, const xAOD::TrackParticleContainer *) const
TVector3 trackMomentum (const xAOD::Vertex *vxCandidate, int trkIndex) const

Private Attributes

bool m_mumu
bool m_mutrk
bool m_trktrk
bool m_allMuons
bool m_combOnly
bool m_atLeastOneComb
bool m_useCombMeasurement
bool m_useV0Fitter
bool m_diMuons
double m_trk1M
double m_trk2M
double m_thresholdPt
double m_higherPt
double m_trkThresholdPt
double m_invMassUpper
double m_invMassLower
double m_collAngleTheta
double m_collAnglePhi
double m_Chi2Cut
bool m_oppChOnly
bool m_sameChOnly
bool m_allChCombs
SG::ReadHandleKey< xAOD::MuonContainerm_muonCollectionKey {this, "muonCollectionKey", "Muons"}
SG::ReadHandleKey< xAOD::TrackParticleContainerm_TrkParticleCollection {this, "TrackParticleCollection", "InDetTrackParticles" }
SG::ReadHandleKeyArray< xAOD::TrackParticleContainerm_MuonTrackKeys {this, "MuonTrackKeys", {}}
PublicToolHandle< Trk::IVertexFitterm_iVertexFitter {this, "TrkVertexFitterTool", "Trk::TrkVKalVrtFitter"}
PublicToolHandle< Trk::IVertexFitterm_iV0VertexFitter {this, "V0VertexFitterTool", "Trk::V0VertexFitter"}
PublicToolHandle< Trk::ITrackSelectorToolm_trkSelector {this, "TrackSelectorTool", "InDet::TrackSelectorTool"}
PublicToolHandle< InDet::VertexPointEstimatorm_vertexEstimator {this, "VertexPointEstimator", "InDet::VertexPointEstimator"}
ServiceHandle< IPartPropSvc > m_partPropSvc {this, "PartPropSvc", "PartPropSvc"}
bool m_mcpCuts
bool m_doTagAndProbe
bool m_forceTagAndProbe

Detailed Description

Definition at line 47 of file JpsiFinder.h.

Constructor & Destructor Documentation

◆ JpsiFinder()

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

Definition at line 102 of file JpsiFinder.cxx.

102 : base_class(t,n,p),
103 m_mumu(true),
104 m_mutrk(false),
105 m_trktrk(false),
106 m_allMuons(false),
107 m_combOnly(false),
108 m_atLeastOneComb(true),
110 m_useV0Fitter(false),
111 m_diMuons(true),
114 m_thresholdPt(0.0),
115 m_higherPt(0.0),
116 m_trkThresholdPt(0.0),
117 m_invMassUpper(100000.0),
118 m_invMassLower(0.0),
119 m_collAngleTheta(0.0),
120 m_collAnglePhi(0.0),
121 m_Chi2Cut(50.),
122 m_oppChOnly(true),
123 m_sameChOnly(false),
124 m_allChCombs(false),
125 m_mcpCuts(true),
126 m_doTagAndProbe(false),
127 m_forceTagAndProbe(false) //forcing T&P method for any charge combinations
128
129 {
130 declareProperty("muAndMu",m_mumu);
131 declareProperty("muAndTrack",m_mutrk);
132 declareProperty("TrackAndTrack",m_trktrk);
133 declareProperty("allMuons",m_allMuons);
134 declareProperty("combOnly",m_combOnly);
135 declareProperty("atLeastOneComb",m_atLeastOneComb);
136 declareProperty("useCombinedMeasurement",m_useCombMeasurement);
137 declareProperty("useV0Fitter",m_useV0Fitter);
138 declareProperty("assumeDiMuons",m_diMuons);
139 declareProperty("track1Mass",m_trk1M);
140 declareProperty("track2Mass",m_trk2M);
141 declareProperty("muonThresholdPt",m_thresholdPt);
142 declareProperty("higherPt",m_higherPt);
143 declareProperty("trackThresholdPt",m_trkThresholdPt);
144 declareProperty("invMassUpper",m_invMassUpper);
145 declareProperty("invMassLower",m_invMassLower);
146 declareProperty("collAngleTheta",m_collAngleTheta);
147 declareProperty("collAnglePhi",m_collAnglePhi);
148 declareProperty("Chi2Cut",m_Chi2Cut);
149 declareProperty("oppChargesOnly",m_oppChOnly);
150 declareProperty("sameChargesOnly",m_sameChOnly);
151 declareProperty("allChargeCombinations",m_allChCombs);
152 declareProperty("useMCPCuts",m_mcpCuts);
153 declareProperty("doTagAndProbe",m_doTagAndProbe);
154 declareProperty("forceTagAndProbe",m_forceTagAndProbe);
155 }
constexpr double muonMassInMeV
the mass of the muon (in MeV)

◆ ~JpsiFinder()

Analysis::JpsiFinder::~JpsiFinder ( )

Definition at line 157 of file JpsiFinder.cxx.

157{ }

Member Function Documentation

◆ fit()

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

Definition at line 410 of file JpsiFinder.cxx.

410 {
411
412 const Trk::TrkV0VertexFitter* concreteVertexFitter=0;
413 if (m_useV0Fitter) {
414 // making a concrete fitter for the V0Fitter
415 concreteVertexFitter = dynamic_cast<const Trk::TrkV0VertexFitter * >(m_iV0VertexFitter.get());
416 if(concreteVertexFitter == 0) {
417 ATH_MSG_FATAL("The vertex fitter passed is not a V0 Vertex Fitter");
418 return NULL;
419 }
420 }
421
422 const Trk::Perigee& aPerigee1 = inputTracks[0]->perigeeParameters();
423 const Trk::Perigee& aPerigee2 = inputTracks[1]->perigeeParameters();
424 int sflag = 0;
425 int errorcode = 0;
426 Amg::Vector3D startingPoint = m_vertexEstimator->getCirclesIntersectionPoint(&aPerigee1,&aPerigee2,sflag,errorcode);
427 if (errorcode != 0) {startingPoint(0) = 0.0; startingPoint(1) = 0.0; startingPoint(2) = 0.0;}
428 if (m_useV0Fitter) {
429 xAOD::Vertex* myVxCandidate = concreteVertexFitter->fit(inputTracks, startingPoint);
430
431 // Added by ASC
432 if(myVxCandidate != 0){
433 std::vector<ElementLink<DataVector<xAOD::TrackParticle> > > newLinkVector;
434 for(unsigned int i=0; i< myVxCandidate->trackParticleLinks().size(); i++)
435 { ElementLink<DataVector<xAOD::TrackParticle> > mylink=myVxCandidate->trackParticleLinks()[i]; //makes a copy (non-const)
436 mylink.setStorableObject(*importedTrackCollection, true);
437 newLinkVector.push_back( mylink ); }
438
439 myVxCandidate->clearTracks();
440 myVxCandidate->setTrackParticleLinks( newLinkVector );
441 }
442
443
444
445 return myVxCandidate;
446 } else {
447 xAOD::Vertex* myVxCandidate = m_iVertexFitter->fit(inputTracks, startingPoint);
448
449 // Added by ASC
450 if(myVxCandidate != 0){
451 std::vector<ElementLink<DataVector<xAOD::TrackParticle> > > newLinkVector;
452 for(unsigned int i=0; i< myVxCandidate->trackParticleLinks().size(); i++)
453 { ElementLink<DataVector<xAOD::TrackParticle> > mylink=myVxCandidate->trackParticleLinks()[i]; //makes a copy (non-const)
454 mylink.setStorableObject(*importedTrackCollection, true);
455 newLinkVector.push_back( mylink ); }
456
457 myVxCandidate->clearTracks();
458 myVxCandidate->setTrackParticleLinks( newLinkVector );
459 }
460
461
462 return myVxCandidate;
463 }
464
465
466
467 return NULL;
468
469 } // End of fit method
#define ATH_MSG_FATAL(x)
PublicToolHandle< InDet::VertexPointEstimator > m_vertexEstimator
Definition JpsiFinder.h:98
PublicToolHandle< Trk::IVertexFitter > m_iV0VertexFitter
Definition JpsiFinder.h:96
PublicToolHandle< Trk::IVertexFitter > m_iVertexFitter
Definition JpsiFinder.h:95
virtual xAOD::Vertex * fit(const std::vector< const xAOD::TrackParticle * > &vectorTrk, const Amg::Vector3D &startingPoint) const override
Interface for xAOD::TrackParticle with Amg::Vector3D starting point.
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
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ getInvariantMass()

double Analysis::JpsiFinder::getInvariantMass ( const JpsiCandidate & jpsiIn,
std::span< const double > massHypotheses ) const

Definition at line 582 of file JpsiFinder.cxx.

582 {
583
584 // construct 4-vectors from track perigee parameters using given mass hypotheses.
585 // NOTE: in new data model (xAOD) the defining parameters are expressed as perigee parameters w.r.t. the beamspot
586 // NOTE2: TrackParticle::p4() method already returns TLorentzVector, however, we want to enforce our own mass hypothesis
587 TLorentzVector mu1;
588 TLorentzVector mu2;
589 mu1.SetVectM(jpsiIn.trackParticle1->p4().Vect(), massHypotheses[0]);
590 mu2.SetVectM(jpsiIn.trackParticle2->p4().Vect(), massHypotheses[1]);
591
592 return (mu1+mu2).M();
593
594 }

◆ getPairs() [1/2]

std::vector< JpsiCandidate > Analysis::JpsiFinder::getPairs ( const std::vector< const xAOD::Muon * > & ) const

◆ getPairs() [2/2]

std::vector< JpsiCandidate > Analysis::JpsiFinder::getPairs ( const std::vector< const xAOD::TrackParticle * > & TracksIn) const

Definition at line 478 of file JpsiFinder.cxx.

478 {
479
480 std::vector<JpsiCandidate> myPairs;
481 JpsiCandidate pair;
482 std::vector<const xAOD::TrackParticle*>::const_iterator outerItr;
483 std::vector<const xAOD::TrackParticle*>::const_iterator innerItr;
484
485 if(TracksIn.size()>=2){
486 for(outerItr=TracksIn.begin();outerItr<TracksIn.end();++outerItr){
487 for(innerItr=(outerItr+1);innerItr!=TracksIn.end();++innerItr){
488 pair.trackParticle1 = *innerItr;
489 pair.trackParticle2 = *outerItr;
490 pair.pairType = TRKTRK;
491 myPairs.push_back(pair);
492 }
493 }
494 }
495
496 return(myPairs);
497 }

◆ getPairs2Colls()

std::vector< JpsiCandidate > Analysis::JpsiFinder::getPairs2Colls ( const std::vector< const xAOD::TrackParticle * > & tracks,
const std::vector< const xAOD::Muon * > & muons,
bool tagAndProbe ) const

Definition at line 535 of file JpsiFinder.cxx.

535 {
536
537 std::vector<JpsiCandidate> myPairs;
538 JpsiCandidate pair;
539
540 // Unless user is running in tag and probe mode, remove tracks which are also identified as muons
541 std::vector<const xAOD::TrackParticle*> tracksToKeep;
542 if (!tagAndProbe) {
543 if(tracks.size()>=1 && muons.size()>=1){
544 for (const xAOD::TrackParticle* trk : tracks) {
545 bool trackIsMuon(false);
546 for (const xAOD::Muon* mu : muons) {
547 auto& link = mu->inDetTrackParticleLink();
548 if ( link.isValid() && *link == trk ) {
549 trackIsMuon=true;
550 break;
551 }
552 }
553 if (!trackIsMuon) tracksToKeep.push_back(trk);
554 }
555 }
556 } else {tracksToKeep = tracks;}
557
558 if(tracksToKeep.size()>=1 && muons.size()>=1){
559 for (const xAOD::TrackParticle* trk : tracksToKeep) {
560 for (const xAOD::Muon* mu : muons) {
561 pair.muon1 = mu;
562 // Muon track 1st
563 pair.trackParticle1 = mu->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
564 pair.trackParticle2 = trk;
565 pair.pairType = MUTRK;
566 myPairs.push_back(pair);
567 }
568 }
569 }
570
571 return(myPairs);
572 }
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Muon_v1 Muon
Reference the current persistent version:

◆ initialize()

StatusCode Analysis::JpsiFinder::initialize ( )
overridevirtual

Definition at line 29 of file JpsiFinder.cxx.

29 {
30
31 // retrieving vertex Fitter
32 ATH_CHECK(m_iVertexFitter.retrieve());
33
34 // retrieving V0 Fitter
35 ATH_CHECK(m_iV0VertexFitter.retrieve(DisableTool{!m_useV0Fitter }));
36
37 // Get the track selector tool from ToolSvc
38 ATH_CHECK(m_trkSelector.retrieve());
39
40 // Get the vertex point estimator tool from ToolSvc
41 ATH_CHECK(m_vertexEstimator.retrieve());
42
43
44 ATH_CHECK(m_muonCollectionKey.initialize());
46 ATH_CHECK(m_MuonTrackKeys.initialize(m_MuonTrackKeys.size() != 0));
47
48 if (m_diMuons) {
49 // Get the Particle Properties Service
50 ATH_CHECK(m_partPropSvc.retrieve());
51 auto particleDataTable = m_partPropSvc->PDT();
52 const HepPDT::ParticleData* pd_mu = particleDataTable->particle(MC::MUON);
53 m_trk1M = pd_mu->mass();
54 m_trk2M = pd_mu->mass();
55 }
56
57 if (m_doTagAndProbe) ATH_MSG_WARNING("You have requested tag and probe mode. Duplicate mu+trk pairs WILL be allowed, charge ordering WILL NOT be done. Tag track will be first in each candidate");
58
59
60// // Check that the user's settings are sensible
61 bool illogicalOptions(false);
62 if ( (m_mumu && m_mutrk) || (m_mumu && m_trktrk) || (m_mutrk && m_trktrk) ) {
63 ATH_MSG_WARNING("You are requesting incompatible combinations of muons and tracks in the pairs. JpsiCandidates will be EMPTY!");
64 illogicalOptions=true;
65 };
67 ATH_MSG_WARNING("You are requesting Tag and Probe analysis but have not requested mu+trk mode. This is impossible. JpsiCandidates will be EMPTY!");
68 illogicalOptions=true;
69 };
70 if ( (m_mutrk || m_trktrk ) && m_useCombMeasurement ) {
71 ATH_MSG_WARNING("You are requesting a combined muon measurement to be used with pairs involving non-muon tracks. Not permitted. JpsiCandidates will be EMPTY!");
72 illogicalOptions=true;
73 };
75 ATH_MSG_WARNING("You are requesting incompatible combinations of combined muons non-combined muons in the pairs. JpsiCandidates will be EMPTY!");
76 illogicalOptions=true;
77 };
79 ATH_MSG_WARNING("You are requesting at least one muon to be combined in a tag and probe analysis. This doesn't make sense. JpsiCandidates will be EMPTY!");
80 illogicalOptions=true;
81 };
83 ATH_MSG_WARNING("You are requesting incompatible combinations of charges in the pairs. JpsiCandidates will be EMPTY!");
84 illogicalOptions=true;
85 };
87 if (!m_forceTagAndProbe){ //if m_forceTagAndProbe=TRUE then T&P will work with any charge combinations
88 ATH_MSG_WARNING("You are requesting same-sign or all-sign combinations in a tag and probe analysis. This doesn't make sense. JpsiCandidates will be EMPTY!");
89 illogicalOptions=true;
90 }
91 }
92 if (illogicalOptions) return StatusCode::FAILURE;;
93
94
95 ATH_MSG_DEBUG("Initialize successful");
96
97 return StatusCode::SUCCESS;
98
99 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
SG::ReadHandleKeyArray< xAOD::TrackParticleContainer > m_MuonTrackKeys
Definition JpsiFinder.h:94
PublicToolHandle< Trk::ITrackSelectorTool > m_trkSelector
Definition JpsiFinder.h:97
SG::ReadHandleKey< xAOD::MuonContainer > m_muonCollectionKey
Definition JpsiFinder.h:92
ServiceHandle< IPartPropSvc > m_partPropSvc
Definition JpsiFinder.h:99
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_TrkParticleCollection
Definition JpsiFinder.h:93
static const int MUON

◆ isContainedIn()

bool Analysis::JpsiFinder::isContainedIn ( const xAOD::TrackParticle * theTrack,
const xAOD::TrackParticleContainer * theCollection ) const

Definition at line 652 of file JpsiFinder.cxx.

652 {
653 return std::find(theCollection->begin(), theCollection->end(), theTrack) != theCollection->end();
654 }
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.

◆ passesMCPCuts()

bool Analysis::JpsiFinder::passesMCPCuts ( const xAOD::Muon * muon) const

Definition at line 642 of file JpsiFinder.cxx.

642 {
643
644 return muon->passesIDCuts();
645
646 }

◆ performSearch()

StatusCode Analysis::JpsiFinder::performSearch ( const EventContext & ctx,
xAOD::VertexContainer & vxContainer ) const
overridevirtual

Definition at line 162 of file JpsiFinder.cxx.

163 {
164 ATH_MSG_DEBUG( "JpsiFinder::performSearch" );
165
166 SG::ReadHandle<xAOD::MuonContainer> muonhandle(m_muonCollectionKey,ctx);
167 if(!muonhandle.isValid()){
168 ATH_MSG_ERROR("Could not retrieve " << m_muonCollectionKey.key());
169 return StatusCode::FAILURE;
170 }
171 // Get the muons from StoreGate
172 const xAOD::MuonContainer* importedMuonCollection = muonhandle.cptr();
173 ATH_MSG_DEBUG("Muon container size "<<importedMuonCollection->size());
174
175 // Get muon tracks if combined measurement requested
176 std::vector<const xAOD::TrackParticleContainer*> importedMuonTrackCollections;
178
179 for (SG::ReadHandle<xAOD::TrackParticleContainer>& handle: m_MuonTrackKeys.makeHandles(ctx)) {
180 if(!handle.isValid()){
181 ATH_MSG_WARNING("No muon TrackParticle collection with name " << handle.key() << " found in StoreGate!");
182 return StatusCode::FAILURE;
183 } else {
184 ATH_MSG_DEBUG("Found muon TrackParticle collection " << handle.key() << " in StoreGate!");
185 ATH_MSG_DEBUG("Muon TrackParticle container size "<< handle->size());
186 importedMuonTrackCollections.push_back(handle.cptr());
187 }
188 }
189 }
190
191 // Get ID tracks
192 SG::ReadHandle<xAOD::TrackParticleContainer> handle(m_TrkParticleCollection,ctx);
193 const xAOD::TrackParticleContainer* importedTrackCollection{nullptr};
194 if(!handle.isValid()){
195 ATH_MSG_WARNING("No TrackParticle collection with name " << handle.key() << " found in StoreGate!");
196 return StatusCode::FAILURE;
197 } else {
198 importedTrackCollection = handle.cptr();
199 ATH_MSG_DEBUG("Found TrackParticle collection " << handle.key() << " in StoreGate!");
200 }
201 ATH_MSG_DEBUG("ID TrackParticle container size "<< handle->size());
202
203 // Typedef for vectors of tracks and muons
204 typedef std::vector<const xAOD::TrackParticle*> TrackBag;
205 typedef std::vector<const xAOD::Muon*> MuonBag;
206
207 // Select the inner detector tracks
208 const xAOD::Vertex* vx = nullptr;
209 TrackBag theIDTracksAfterSelection;
210 if (m_trktrk || m_mutrk) {
212 for (trkCItr=importedTrackCollection->begin(); trkCItr!=importedTrackCollection->end(); ++trkCItr) {
213 const xAOD::TrackParticle* TP = (*trkCItr);
214 if ( fabs(TP->pt())<m_trkThresholdPt ) continue;
215 if ( !m_trkSelector->decision(*TP, vx) ) continue;
216 theIDTracksAfterSelection.push_back(TP);
217 }
218 if (theIDTracksAfterSelection.size() == 0) return StatusCode::SUCCESS;;
219 ATH_MSG_DEBUG("Number of tracks after ID track selection: " << theIDTracksAfterSelection.size());
220 }
221
222 // Select the muons
223 MuonBag theMuonsAfterSelection;
224 if (m_mumu || m_mutrk) {
225 for (auto muItr=importedMuonCollection->begin(); muItr!=importedMuonCollection->end(); ++muItr) {
226 if ( *muItr == NULL ) continue;
227 if (!(*muItr)->inDetTrackParticleLink().isValid()) continue; // No muons without ID tracks
228 const xAOD::TrackParticle* muonTrk = *((*muItr)->inDetTrackParticleLink());
229 if ( muonTrk==NULL) continue;
230 if ( !m_trkSelector->decision(*muonTrk, vx) ) continue; // all ID tracks must pass basic tracking cuts
231 if ( fabs(muonTrk->pt())<m_thresholdPt ) continue; // higher pt cut if needed
232 if ( m_mcpCuts && !passesMCPCuts(*muItr)) continue; // MCP cuts
233 if ( m_combOnly && (*muItr)->muonType() != xAOD::Muon::Combined ) continue; // require combined muons
234 if ( (*muItr)->muonType() == xAOD::Muon::SiliconAssociatedForwardMuon && !m_useCombMeasurement) continue;
235 theMuonsAfterSelection.push_back(*muItr);
236 }
237 if (theMuonsAfterSelection.size() == 0) return StatusCode::SUCCESS;;
238 ATH_MSG_DEBUG("Number of muons after selection: " << theMuonsAfterSelection.size());
239 }
240
241 // Sort into pairs - end result will be a vector of JpsiCandidate structs
242 std::vector<JpsiCandidate> jpsiCandidates;
243 if (m_mumu) jpsiCandidates = getPairs(theMuonsAfterSelection);
244 if (m_trktrk) jpsiCandidates = getPairs(theIDTracksAfterSelection);
245 if (m_mutrk) jpsiCandidates = getPairs2Colls(theIDTracksAfterSelection,theMuonsAfterSelection,m_doTagAndProbe);
246
247
248 // (1) Enforce one combined muon
249 if (m_atLeastOneComb) {
250 std::vector<JpsiCandidate> selectCandidates;
251 for(auto &cand : jpsiCandidates) { if(cand.muonTypes!=TT) selectCandidates.push_back(cand); }
252 selectCandidates.swap(jpsiCandidates);
253 ATH_MSG_DEBUG("Number of candidates after requirement of at least 1 combined muon: " << jpsiCandidates.size() );
254 }
255
256 // (2) Establish track content for candidates
257 // and set the appropriate track collections for the combined muon tracks where appropriate (for saving to persistency later)
258
259 // mu+trk or trk+trk - always ID track collection
260 if (m_mutrk || m_trktrk) {
261 for (auto jpsiItr=jpsiCandidates.begin(); jpsiItr!=jpsiCandidates.end(); ++jpsiItr) {
262 (*jpsiItr).collection1 = importedTrackCollection;
263 (*jpsiItr).collection2 = importedTrackCollection;
264 }
265 }
266
267 if (m_mumu) {
268 for (auto jpsiItr=jpsiCandidates.begin(); jpsiItr!=jpsiCandidates.end(); ++jpsiItr) {
270 (*jpsiItr).trackParticle1 = (*jpsiItr).muon1->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
271 (*jpsiItr).trackParticle2 = (*jpsiItr).muon2->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
272 (*jpsiItr).collection1 = importedTrackCollection;
273 (*jpsiItr).collection2 = importedTrackCollection;
274 }
276 if (!(*jpsiItr).muon1->combinedTrackParticleLink().isValid()) {
277 (*jpsiItr).trackParticle1 = (*jpsiItr).muon1->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
278 (*jpsiItr).collection1 = importedTrackCollection;
279 }
280
281 if (!(*jpsiItr).muon2->combinedTrackParticleLink().isValid()) {
282 (*jpsiItr).trackParticle2 = (*jpsiItr).muon2->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
283 (*jpsiItr).collection2 = importedTrackCollection;
284 }
285
286 if ((*jpsiItr).muon1->combinedTrackParticleLink().isValid()) {
287 (*jpsiItr).trackParticle1 = (*jpsiItr).muon1->trackParticle( xAOD::Muon::CombinedTrackParticle );
288 bool foundCollection(false);
289 // Look for correct muon track container
290 for (auto muTrkCollItr=importedMuonTrackCollections.begin(); muTrkCollItr!=importedMuonTrackCollections.end(); ++muTrkCollItr) {
291 if (isContainedIn((*jpsiItr).trackParticle1,*muTrkCollItr)) { (*jpsiItr).collection1 = *muTrkCollItr; foundCollection=true; break;}
292 }
293 if (!foundCollection) { // didn't find the correct muon track container so go back to ID
294 (*jpsiItr).trackParticle1 = (*jpsiItr).muon1->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
295 (*jpsiItr).collection1 = importedTrackCollection;
296 ATH_MSG_WARNING("Muon track from muon of author " << (*jpsiItr).muon1->author() << " not found in muon track collections you have provided.");
297 ATH_MSG_WARNING("Defaulting to ID track collection - combined measurement will not be used");
298 }
299 }
300
301 if ((*jpsiItr).muon2->combinedTrackParticleLink().isValid()) {
302 (*jpsiItr).trackParticle2 = (*jpsiItr).muon2->trackParticle( xAOD::Muon::CombinedTrackParticle );
303 bool foundCollection(false);
304 // Look for correct muon track container
305 for (auto muTrkCollItr=importedMuonTrackCollections.begin(); muTrkCollItr!=importedMuonTrackCollections.end(); ++muTrkCollItr) {
306 if (isContainedIn((*jpsiItr).trackParticle2,*muTrkCollItr)) { (*jpsiItr).collection2 = *muTrkCollItr; foundCollection=true; break;}
307 }
308 if (!foundCollection) { // didn't find the correct muon track container so go back to ID
309 (*jpsiItr).trackParticle2 = (*jpsiItr).muon2->trackParticle( xAOD::Muon::InnerDetectorTrackParticle );
310 (*jpsiItr).collection2 = importedTrackCollection;
311 ATH_MSG_WARNING("Muon track from muon of author " << (*jpsiItr).muon2->author() << " not found in muon track collections you have provided.");
312 ATH_MSG_WARNING("Defaulting to ID track collection - combined measurement will not be used");
313 }
314 }
315 } // combined measurement
316 } // iteration over candidates
317 }
318
319
320 // (3) Enforce higher track pt if requested
321 if (m_higherPt>0.0) {
322 std::vector<JpsiCandidate> selectCandidates;
323 for(auto& cand : jpsiCandidates){
324 bool reject = (fabs(cand.trackParticle1->pt()) < m_higherPt) && (fabs(cand.trackParticle2->pt()) < m_higherPt);
325 if(!reject) selectCandidates.push_back(cand);
326 }
327 selectCandidates.swap(jpsiCandidates);
328 ATH_MSG_DEBUG("Number of candidates after higherPt cut: " << jpsiCandidates.size() );
329 }
330
331 // (4) Select all opp/same charged track pairs
332 std::vector<JpsiCandidate> sortedJpsiCandidates;
333 sortedJpsiCandidates = selectCharges(jpsiCandidates);
334
335 ATH_MSG_DEBUG("Number of candidates after charge selection: " << sortedJpsiCandidates.size() );
336
337 // (5) Select for decay angle, if requested
338 if (m_collAnglePhi>0.0 && m_collAngleTheta>0.0) {
339 std::vector<JpsiCandidate> selectCandidates;
340 for(auto& cand : sortedJpsiCandidates){
341 double deltatheta = fabs( cand.trackParticle1->theta() - cand.trackParticle2->theta() );
342 double deltaphi = std::abs(xAOD::P4Helpers::deltaPhi(cand.trackParticle1->phi0() , cand.trackParticle2->phi0()));
343 bool reject = (deltatheta > m_collAngleTheta) || (deltaphi > m_collAnglePhi);
344 if(!reject) selectCandidates.push_back(cand);
345 }
346 sortedJpsiCandidates.swap(selectCandidates);
347 ATH_MSG_DEBUG("Number of collimated candidates: " << sortedJpsiCandidates.size() );
348 }
349
350 // (6) Select for invariant mass, if requested
351 std::array<double,2> trkMasses{m_trk1M, m_trk2M};
352 if ( (m_invMassLower > 0.0) || (m_invMassUpper > 0.0) ) {
353 std::vector<JpsiCandidate> selectCandidates;
354 for(auto& cand : sortedJpsiCandidates){
355 double invMass = getInvariantMass(cand,trkMasses);
356 bool reject = invMass < m_invMassLower || invMass > m_invMassUpper;
357 if(!reject) selectCandidates.push_back(cand);
358 }
359 selectCandidates.swap(sortedJpsiCandidates);
360 ATH_MSG_DEBUG("Number of candidates passing invariant mass selection: " << sortedJpsiCandidates.size() );
361 }
362
363 ATH_MSG_DEBUG("Number of pairs passing all selections and going to vertexing: " << sortedJpsiCandidates.size() );
364 if (sortedJpsiCandidates.size() == 0) return StatusCode::SUCCESS;;
365 std::vector<const xAOD::TrackParticle*> theTracks;
366 std::vector<const xAOD::Muon*> theStoredMuons;
367 // Fit each pair of tracks to a vertex
368 for(auto jpsiItr=sortedJpsiCandidates.begin(); jpsiItr!=sortedJpsiCandidates.end(); ++jpsiItr) {
369 theTracks.clear();
370 theTracks.push_back((*jpsiItr).trackParticle1);
371 theTracks.push_back((*jpsiItr).trackParticle2);
372 std::unique_ptr<xAOD::Vertex> myVxCandidate {fit(theTracks,importedTrackCollection)}; // This line actually does the fitting and object making
373 if (myVxCandidate) {
374 // Chi2 cut if requested
375 double chi2 = myVxCandidate->chiSquared();
376 ATH_MSG_DEBUG("chi2 is: " << chi2);
377 if (m_Chi2Cut <= 0.0 || chi2 <= m_Chi2Cut) {
378 // decorate the candidate with refitted tracks and muons via the BPhysHelper
379 xAOD::BPhysHelper jpsiHelper(myVxCandidate.get());
380 bool validtrk = jpsiHelper.setRefTrks();
381 if(!validtrk) ATH_MSG_WARNING("Problem setting tracks " << __FILE__ << ':' << __LINE__);
382 if (m_mumu || m_mutrk) {
383 theStoredMuons.clear();
384 theStoredMuons.push_back((*jpsiItr).muon1);
385 if (m_mumu) theStoredMuons.push_back((*jpsiItr).muon2);
386 bool valid = jpsiHelper.setMuons(theStoredMuons,importedMuonCollection);
387 if(!valid) ATH_MSG_WARNING("Problem setting muons " << __FILE__ << ':' << __LINE__);
388 }
389 // Retain the vertex
390 vxContainer.push_back(std::move(myVxCandidate));
391 }
392 } else { // fit failed
393 ATH_MSG_DEBUG("Fitter failed!");
394 // Don't try to delete the object, since we arrived here,
395 // because this pointer is null...
396 //delete myVxCandidate;
397 }
398 }
399 ATH_MSG_DEBUG("vxContainer size " << vxContainer.size());
400
401 return StatusCode::SUCCESS;;
402 }
#define ATH_MSG_ERROR(x)
xAOD::Vertex * fit(const std::vector< const xAOD::TrackParticle * > &, const xAOD::TrackParticleContainer *importedTrackCollection) const
std::vector< JpsiCandidate > getPairs(const std::vector< const xAOD::TrackParticle * > &) const
std::vector< JpsiCandidate > getPairs2Colls(const std::vector< const xAOD::TrackParticle * > &, const std::vector< const xAOD::Muon * > &, bool) const
bool passesMCPCuts(const xAOD::Muon *) const
std::vector< JpsiCandidate > selectCharges(const std::vector< JpsiCandidate > &) const
bool isContainedIn(const xAOD::TrackParticle *, const xAOD::TrackParticleContainer *) const
double getInvariantMass(const JpsiCandidate &, std::span< const double >) const
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.
size_type size() const noexcept
Returns the number of elements in the collection.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
double chi2(TH1 *h0, TH1 *h1)
std::vector< const xAOD::TrackParticle * > TrackBag
double invMass(const I4Momentum &pA, const I4Momentum &pB)
invariant mass from two I4momentum references
Definition P4Helpers.h:252
list valid
Definition calibdata.py:44
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".

◆ selectCharges()

std::vector< JpsiCandidate > Analysis::JpsiFinder::selectCharges ( const std::vector< JpsiCandidate > & jpsisIn) const

Definition at line 602 of file JpsiFinder.cxx.

602 {
603
604 bool opposite(false),same(false),all(false);
605 opposite=m_oppChOnly;
608
609 JpsiCandidate tmpJpsi;
610 std::vector<JpsiCandidate> jpsis;
611 double qOverP1=0.;
612 double qOverP2=0.;
613 for(auto jpsiItr=jpsisIn.cbegin();jpsiItr!=jpsisIn.cend();jpsiItr++){
614 bool oppCh(false),sameCh(false);
615 tmpJpsi = *jpsiItr;
616 qOverP1=(*jpsiItr).trackParticle1->qOverP();
617 qOverP2=(*jpsiItr).trackParticle2->qOverP();
618 if(qOverP1*qOverP2<0.0) oppCh=true; // product charge < 0
619 if(qOverP1*qOverP2>0.0) sameCh=true; // product charge > 0
620 // +ve should be first so swap
621 // Don't do it for tag and probe analyses (because tag muon must not change position)
622 if (oppCh && qOverP1<0.0 && !m_doTagAndProbe && !m_mutrk) {
623 tmpJpsi.trackParticle1 = (*jpsiItr).trackParticle2;
624 tmpJpsi.trackParticle2 = (*jpsiItr).trackParticle1;
625 tmpJpsi.muon1 = (*jpsiItr).muon2;
626 tmpJpsi.muon2 = (*jpsiItr).muon1;
627 tmpJpsi.collection1 = (*jpsiItr).collection2;
628 tmpJpsi.collection2 = (*jpsiItr).collection1;
629 }
630 if (oppCh && (opposite || all) ) jpsis.push_back(tmpJpsi);
631 if (sameCh && (same || all) ) jpsis.push_back(tmpJpsi);
632
633 } // end of for loop
634
635 return(jpsis);
636 }

◆ trackMomentum()

TVector3 Analysis::JpsiFinder::trackMomentum ( const xAOD::Vertex * vxCandidate,
int trkIndex ) const

Definition at line 660 of file JpsiFinder.cxx.

661 {
662 double px = 0., py = 0., pz = 0.;
663 if (0 != vxCandidate) {
664 const Trk::TrackParameters* aPerigee = vxCandidate->vxTrackAtVertex()[trkIndex].perigeeAtVertex();
665 px = aPerigee->momentum()[Trk::px];
666 py = aPerigee->momentum()[Trk::py];
667 pz = aPerigee->momentum()[Trk::pz];
668 }
669 TVector3 mom(px,py,pz);
670 return mom;
671 }
const Amg::Vector3D & momentum() const
Access method for the momentum.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
@ pz
global momentum (cartesian)
Definition ParamDefs.h:61
@ px
Definition ParamDefs.h:59
@ py
Definition ParamDefs.h:60
ParametersBase< TrackParametersDim, Charged > TrackParameters

Member Data Documentation

◆ m_allChCombs

bool Analysis::JpsiFinder::m_allChCombs
private

Definition at line 91 of file JpsiFinder.h.

◆ m_allMuons

bool Analysis::JpsiFinder::m_allMuons
private

Definition at line 73 of file JpsiFinder.h.

◆ m_atLeastOneComb

bool Analysis::JpsiFinder::m_atLeastOneComb
private

Definition at line 75 of file JpsiFinder.h.

◆ m_Chi2Cut

double Analysis::JpsiFinder::m_Chi2Cut
private

Definition at line 88 of file JpsiFinder.h.

◆ m_collAnglePhi

double Analysis::JpsiFinder::m_collAnglePhi
private

Definition at line 87 of file JpsiFinder.h.

◆ m_collAngleTheta

double Analysis::JpsiFinder::m_collAngleTheta
private

Definition at line 86 of file JpsiFinder.h.

◆ m_combOnly

bool Analysis::JpsiFinder::m_combOnly
private

Definition at line 74 of file JpsiFinder.h.

◆ m_diMuons

bool Analysis::JpsiFinder::m_diMuons
private

Definition at line 78 of file JpsiFinder.h.

◆ m_doTagAndProbe

bool Analysis::JpsiFinder::m_doTagAndProbe
private

Definition at line 101 of file JpsiFinder.h.

◆ m_forceTagAndProbe

bool Analysis::JpsiFinder::m_forceTagAndProbe
private

Definition at line 102 of file JpsiFinder.h.

◆ m_higherPt

double Analysis::JpsiFinder::m_higherPt
private

Definition at line 82 of file JpsiFinder.h.

◆ m_invMassLower

double Analysis::JpsiFinder::m_invMassLower
private

Definition at line 85 of file JpsiFinder.h.

◆ m_invMassUpper

double Analysis::JpsiFinder::m_invMassUpper
private

Definition at line 84 of file JpsiFinder.h.

◆ m_iV0VertexFitter

PublicToolHandle< Trk::IVertexFitter > Analysis::JpsiFinder::m_iV0VertexFitter {this, "V0VertexFitterTool", "Trk::V0VertexFitter"}
private

Definition at line 96 of file JpsiFinder.h.

96{this, "V0VertexFitterTool", "Trk::V0VertexFitter"};

◆ m_iVertexFitter

PublicToolHandle< Trk::IVertexFitter > Analysis::JpsiFinder::m_iVertexFitter {this, "TrkVertexFitterTool", "Trk::TrkVKalVrtFitter"}
private

Definition at line 95 of file JpsiFinder.h.

95{this, "TrkVertexFitterTool", "Trk::TrkVKalVrtFitter"};

◆ m_mcpCuts

bool Analysis::JpsiFinder::m_mcpCuts
private

Definition at line 100 of file JpsiFinder.h.

◆ m_mumu

bool Analysis::JpsiFinder::m_mumu
private

Definition at line 70 of file JpsiFinder.h.

◆ m_muonCollectionKey

SG::ReadHandleKey<xAOD::MuonContainer> Analysis::JpsiFinder::m_muonCollectionKey {this, "muonCollectionKey", "Muons"}
private

Definition at line 92 of file JpsiFinder.h.

92{this, "muonCollectionKey", "Muons"};

◆ m_MuonTrackKeys

SG::ReadHandleKeyArray<xAOD::TrackParticleContainer> Analysis::JpsiFinder::m_MuonTrackKeys {this, "MuonTrackKeys", {}}
private

Definition at line 94 of file JpsiFinder.h.

94{this, "MuonTrackKeys", {}};

◆ m_mutrk

bool Analysis::JpsiFinder::m_mutrk
private

Definition at line 71 of file JpsiFinder.h.

◆ m_oppChOnly

bool Analysis::JpsiFinder::m_oppChOnly
private

Definition at line 89 of file JpsiFinder.h.

◆ m_partPropSvc

ServiceHandle<IPartPropSvc> Analysis::JpsiFinder::m_partPropSvc {this, "PartPropSvc", "PartPropSvc"}
private

Definition at line 99 of file JpsiFinder.h.

99{this, "PartPropSvc", "PartPropSvc"};

◆ m_sameChOnly

bool Analysis::JpsiFinder::m_sameChOnly
private

Definition at line 90 of file JpsiFinder.h.

◆ m_thresholdPt

double Analysis::JpsiFinder::m_thresholdPt
private

Definition at line 81 of file JpsiFinder.h.

◆ m_trk1M

double Analysis::JpsiFinder::m_trk1M
private

Definition at line 79 of file JpsiFinder.h.

◆ m_trk2M

double Analysis::JpsiFinder::m_trk2M
private

Definition at line 80 of file JpsiFinder.h.

◆ m_TrkParticleCollection

SG::ReadHandleKey<xAOD::TrackParticleContainer> Analysis::JpsiFinder::m_TrkParticleCollection {this, "TrackParticleCollection", "InDetTrackParticles" }
private

Definition at line 93 of file JpsiFinder.h.

93{this, "TrackParticleCollection", "InDetTrackParticles" };

◆ m_trkSelector

PublicToolHandle< Trk::ITrackSelectorTool > Analysis::JpsiFinder::m_trkSelector {this, "TrackSelectorTool", "InDet::TrackSelectorTool"}
private

Definition at line 97 of file JpsiFinder.h.

97{this, "TrackSelectorTool", "InDet::TrackSelectorTool"};

◆ m_trkThresholdPt

double Analysis::JpsiFinder::m_trkThresholdPt
private

Definition at line 83 of file JpsiFinder.h.

◆ m_trktrk

bool Analysis::JpsiFinder::m_trktrk
private

Definition at line 72 of file JpsiFinder.h.

◆ m_useCombMeasurement

bool Analysis::JpsiFinder::m_useCombMeasurement
private

Definition at line 76 of file JpsiFinder.h.

◆ m_useV0Fitter

bool Analysis::JpsiFinder::m_useV0Fitter
private

Definition at line 77 of file JpsiFinder.h.

◆ m_vertexEstimator

PublicToolHandle< InDet::VertexPointEstimator > Analysis::JpsiFinder::m_vertexEstimator {this, "VertexPointEstimator", "InDet::VertexPointEstimator"}
private

Definition at line 98 of file JpsiFinder.h.

98{this, "VertexPointEstimator", "InDet::VertexPointEstimator"};

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