ATLAS Offline Software
Loading...
Searching...
No Matches
FourMuonTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5// ****************************************************************************
6// ----------------------------------------------------------------------------
7// FourMuonTool
8// James Catmore <James.Catmore@cern.ch>
9// Evelina Bouhova-Thacker <e.bouhova@cern.ch>
10// ----------------------------------------------------------------------------
11// ****************************************************************************
12
13#include "FourMuonTool.h"
14#include "BPhysPVTools.h"
19#include "AthLinks/ElementLink.h"
20
21#include "xAODTracking/Vertex.h"
27
28#include <algorithm>
29
30namespace DerivationFramework {
31
33
34 // retrieving vertex Fitter
35 ATH_CHECK( m_iVertexFitter.retrieve() );
36 ATH_MSG_DEBUG("Retrieved tool " << m_iVertexFitter);
37
38 // retrieving V0 Fitter
39 ATH_CHECK( m_iV0VertexFitter.retrieve(DisableTool{!m_useV0Fitter}));
40
41 // Get the track selector tool from ToolSvc
42 ATH_CHECK ( m_trkSelector.retrieve() );
43 ATH_MSG_DEBUG("Retrieved tool " << m_trkSelector);
44
45 // Get the beam spot service
46 ATH_CHECK( m_eventInfo_key.initialize() );
47 ATH_CHECK(m_muonCollectionKey.initialize());
49 ATH_CHECK(m_muonIndex.initialize());
50 ATH_MSG_DEBUG("Initialize successful");
51
52 return StatusCode::SUCCESS;
53
54 }
55
56 FourMuonTool::FourMuonTool(const std::string& t, const std::string& n, const IInterface* p) : AthAlgTool(t,n,p)
57 {
58 }
59
61
62 //-------------------------------------------------------------------------------------
63 // Find the candidates
64 //-------------------------------------------------------------------------------------
65 StatusCode FourMuonTool::performSearch(xAOD::VertexContainer* pairVxContainer, xAOD::VertexContainer* quadVxContainer, bool &selectEvent, const EventContext& ctx) const
66 {
67 ATH_MSG_DEBUG( "FourMuonTool::performSearch" );
68 selectEvent = false;
69
70 // Get the muons from StoreGate
72 ATH_CHECK(importedMuonCollection.isValid());
73 ATH_MSG_DEBUG("Muon container size "<<importedMuonCollection->size());
74
75 // Get ID tracks
77 ATH_CHECK(importedTrackCollection.isValid());
78 ATH_MSG_DEBUG("ID TrackParticle container size "<< importedTrackCollection->size());
79
80 // Select the muons
81 std::vector<const xAOD::Muon*> theMuonsAfterSelection;
83 unsigned int nCombMuons = 0;
84 unsigned int nSegmentTaggedMuons = 0;
85
86 for (const auto * muon : *importedMuonCollection) {
87 if ( !muon ) continue;
88 muonDecorator(*muon) = -1; // all muons must be decorated
89 if ( (muon->muonType() != xAOD::Muon::MuonType::Combined ) && (muon->muonType() != xAOD::Muon::MuonType::SegmentTagged ) ) continue;
90 const xAOD::TrackParticle* muonTrk = muon->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle);
91 if ( !muonTrk ) continue;
92 const xAOD::Vertex* vx{};
93 if ( !m_trkSelector->decision(*muonTrk, vx) ) continue; // all ID tracks must pass basic tracking cuts
94 if ( fabs(muonTrk->pt())<m_ptCut ) continue; // pt cut
95 if ( fabs(muonTrk->eta())>m_etaCut ) continue; // eta cut
96 if ( muon->muonType() == xAOD::Muon::MuonType::Combined ) ++nCombMuons;
97 if ( muon->muonType() == xAOD::Muon::MuonType::SegmentTagged ) ++nSegmentTaggedMuons;
98 theMuonsAfterSelection.push_back(muon);
99 }
100 unsigned int nSelectedMuons = theMuonsAfterSelection.size();
101 ATH_MSG_DEBUG("Number of muons after selection: " << nSelectedMuons);
102 ATH_MSG_DEBUG("of which " << nCombMuons << " are combined");
103 ATH_MSG_DEBUG("and " << nSegmentTaggedMuons << " are segment tagged");
104 if ( (nSelectedMuons < 4) || (nCombMuons < 1) ) {
105 ATH_MSG_DEBUG("Muon criteria not met. Skipping event.");
106 return StatusCode::SUCCESS;
107 }
108 selectEvent = true; // if we got this far we should definitively accept the event
109
110 // Decorators
111 SG::AuxElement::Decorator< std::string > indexDecorator("CombinationCode");
112 SG::AuxElement::Decorator< std::string > chargeDecorator("ChargeCode");
113
114 // Order by pT
115 std::sort(theMuonsAfterSelection.begin(), theMuonsAfterSelection.end(), [](const xAOD::Muon *a, const xAOD::Muon *b) {
116 return b->pt() < a->pt();
117 });
118
119 // Decorate the selected muons (now pT ordered) with their index
120 unsigned int muonIndex(0);
121 for (auto selMuon : theMuonsAfterSelection) {
122 muonDecorator(*selMuon) = muonIndex;
123 ++muonIndex;
124 }
125
126 // Quadruplet combinatorics
127 std::vector<Combination> quadruplets;
128 std::vector<Combination> pairs;
129 buildCombinations(theMuonsAfterSelection,pairs,quadruplets,nSelectedMuons);
130 if (quadruplets.size()==0) {
131 ATH_MSG_DEBUG("No acceptable quadruplets");
132 return StatusCode::SUCCESS;
133 }
134
135 // Get the beam spot (for the vertexing starting point)
137 ATH_CHECK(evt.isValid());
138 const Amg::Vector3D beamSpot(evt->beamPosX(), evt->beamPosY(), evt->beamPosZ());
139
140 // fit pairs
141 ATH_MSG_DEBUG("Successful pairs.....");
142 for (std::vector<Combination>::iterator pairItr = pairs.begin(); pairItr!=pairs.end(); ++pairItr) {
143 std::vector<const xAOD::TrackParticle*> theTracks = (*pairItr).trackParticles("pair1");
144 std::unique_ptr<xAOD::Vertex> pairVxCandidate = fit(ctx,theTracks,importedTrackCollection.get(),beamSpot); // This line actually does the fitting and object making
145 if (pairVxCandidate) {
146 // decorate the candidate with its codes
147 indexDecorator(*pairVxCandidate) = (*pairItr).combinationIndices();
148 chargeDecorator(*pairVxCandidate) = (*pairItr).combinationCharges();
149 // decorate the candidate with refitted tracks and muons via the BPhysHelper
150 xAOD::BPhysHelper helper(pairVxCandidate.get());
151 helper.setRefTrks();
152 std::vector<const xAOD::Muon*> theStoredMuons;
153 theStoredMuons = (*pairItr).muons;
154 helper.setMuons(theStoredMuons,importedMuonCollection.get());
155 ATH_MSG_DEBUG("..... indices: " << (*pairItr).combinationIndices() <<
156 " charges: " << (*pairItr).combinationCharges() <<
157 " chi2: " << pairVxCandidate->chiSquared());
158 // Retain the vertex
159 pairVxContainer->push_back(std::move(pairVxCandidate));
160 } else { // fit failed
161 ATH_MSG_DEBUG("Fitter failed!");
162 }
163 }
164 ATH_MSG_DEBUG("pairContainer size " << pairVxContainer->size());
165
166 // fit quadruplets
167 ATH_MSG_DEBUG("Successful quadruplets.....");
168 for (std::vector<Combination>::iterator quadItr = quadruplets.begin(); quadItr!=quadruplets.end(); ++quadItr) {
169 std::vector<const xAOD::TrackParticle*> theDCTracks; theDCTracks.clear();
170 theDCTracks = (*quadItr).trackParticles("DC");
171 std::unique_ptr<xAOD::Vertex> dcVxCandidate = fit(ctx,theDCTracks,importedTrackCollection.get(), beamSpot);
172 if (dcVxCandidate != 0) {
173 // decorate the candidate with its codes
174 indexDecorator(*dcVxCandidate) = (*quadItr).combinationIndices();
175 chargeDecorator(*dcVxCandidate) = (*quadItr).combinationCharges();
176 // Decorate the DC candidate with the differences between its chi2 and the other
177 double dcChi2 = dcVxCandidate->chiSquared();
178 // decorate the candidate with refitted tracks and muons via the BPhysHelper
179 xAOD::BPhysHelper helper(dcVxCandidate.get());
180 helper.setRefTrks();
181 const std::vector<const xAOD::Muon*> &theStoredMuons = (*quadItr).muons;
182 helper.setMuons(theStoredMuons,importedMuonCollection.get());
183 // Retain the vertex
184 quadVxContainer->push_back(std::move(dcVxCandidate));
185 ATH_MSG_DEBUG("..... indices: " << (*quadItr).combinationIndices() <<
186 " charges: " << (*quadItr).combinationCharges() <<
187 " chi2(DC): " << dcChi2);
188 } else { // fit failed
189 ATH_MSG_DEBUG("Fitter failed!");
190 }
191 }
192 ATH_MSG_DEBUG("quadruplet container size " << quadVxContainer->size());
193
194 return StatusCode::SUCCESS;;
195 }
196
197 // *********************************************************************************
198
199 // ---------------------------------------------------------------------------------
200 // fit - does the fit
201 // ---------------------------------------------------------------------------------
202
203 std::unique_ptr<xAOD::Vertex> FourMuonTool::fit(const EventContext& ctx,
204 const std::vector<const xAOD::TrackParticle*> &inputTracks,
205 const xAOD::TrackParticleContainer* importedTrackCollection,
206 const Amg::Vector3D &beamSpot) const {
207
208 const Trk::TrkV0VertexFitter* concreteVertexFitter=0;
209 if (m_useV0Fitter) {
210 // making a concrete fitter for the V0Fitter
211 concreteVertexFitter = dynamic_cast<const Trk::TrkV0VertexFitter * >(&(*m_iV0VertexFitter));
212 if(concreteVertexFitter == 0) {
213 ATH_MSG_FATAL("The vertex fitter passed is not a V0 Vertex Fitter");
214 return nullptr;
215 }
216 }
217
218 std::unique_ptr<xAOD::Vertex> myVxCandidate;
219 if (m_useV0Fitter) {
220 myVxCandidate = concreteVertexFitter->fit(ctx, inputTracks, beamSpot /*vertex startingPoint*/ );
221 } else {
222 myVxCandidate = m_iVertexFitter->fit(ctx, inputTracks, beamSpot /*vertex startingPoint*/ );
223 }
224
225 if(myVxCandidate) BPhysPVTools::PrepareVertexLinks(myVxCandidate.get(), importedTrackCollection);
226
227 return myVxCandidate;
228
229 } // End of fit method
230
231
232 // *********************************************************************************
233
234 // ---------------------------------------------------------------------------------
235 // getQuadIndices: forms up index lists
236 // ---------------------------------------------------------------------------------
237
238 std::vector<std::vector<unsigned int> > FourMuonTool::getQuadIndices(unsigned int length) {
239
240 std::vector<std::vector<unsigned int> > quadIndices = mFromN(4,length);
241 return(quadIndices);
242 }
243
244
245 // *********************************************************************************
246
247 // ---------------------------------------------------------------------------------
248 // mFromN and combinatorics
249 // ---------------------------------------------------------------------------------
250 std::vector<std::vector<unsigned int> > FourMuonTool::mFromN(unsigned int m, unsigned int N) {
251
252 std::vector<std::vector<unsigned int> > allCombinations;
253 std::vector<unsigned int> mainList;
254 std::vector<unsigned int> combination;
255 for (unsigned int i=0; i<N; ++i) mainList.push_back(i);
256 combinatorics(0,m,combination,mainList,allCombinations);
257 return allCombinations;
258 }
259
260 void FourMuonTool::combinatorics(unsigned int offset,
261 unsigned int k,
262 std::vector<unsigned int> &combination,
263 std::vector<unsigned int> &mainList,
264 std::vector<std::vector<unsigned int> > &allCombinations) {
265 if (k==0) {
266 allCombinations.push_back(combination);
267 return;
268 }
269 if (k>0) {
270 for (unsigned int i=offset; i<=mainList.size()-k; ++i) {
271 combination.push_back(mainList[i]);
272 combinatorics(i+1,k-1,combination,mainList,allCombinations);
273 combination.pop_back();
274 }
275 }
276 }
277
278 // ---------------------------------------------------------------------------------
279 // getPairIndices
280 // ---------------------------------------------------------------------------------
281
282 std::vector<std::pair<unsigned int, unsigned int> > FourMuonTool::getPairIndices(unsigned int length){
283
284 std::vector<std::pair<unsigned int, unsigned int> > uniquePairs;
285 std::vector<std::vector<unsigned int> > doublets = mFromN(2,length);
286 for (std::vector<std::vector<unsigned int> >::iterator it=doublets.begin(); it!=doublets.end(); ++it) {
287 std::pair<unsigned int, unsigned int> tmpPair = std::make_pair((*it).at(0),(*it).at(1));
288 uniquePairs.push_back(tmpPair);
289 }
290
291 return(uniquePairs);
292 }
293
294
295
296 // *********************************************************************************
297
298 // ---------------------------------------------------------------------------------
299 // buildCombinations: forms up the quadruplet of muons/tracks
300 // ---------------------------------------------------------------------------------
301
302 void FourMuonTool::buildCombinations(const std::vector<const xAOD::Muon*> &muonsIn,
303 std::vector<Combination> &pairs,
304 std::vector<Combination> &quadruplets,
305 unsigned int nSelectedMuons) {
306
307 std::vector<std::vector<unsigned int> > quadrupletIndices = getQuadIndices(nSelectedMuons);
308 std::vector<std::pair<unsigned int, unsigned int> > pairIndices = getPairIndices(nSelectedMuons);
309
310 // Quadruplets
311 std::vector<std::vector<unsigned int> >::iterator quadItr;
312 for (quadItr=quadrupletIndices.begin(); quadItr!=quadrupletIndices.end(); ++quadItr) {
313 const std::vector<unsigned int> &quad = (*quadItr);
314 std::vector<const xAOD::Muon*> theMuons = {muonsIn[quad[0]],muonsIn[quad[1]],muonsIn[quad[2]],muonsIn[quad[3]]};
315 if (!passesQuadSelection(theMuons)) continue;
316 Combination tmpQuad;
317 tmpQuad.muons = std::move(theMuons);
318 tmpQuad.quadIndices = quad;
319 quadruplets.emplace_back(std::move(tmpQuad));
320 }
321 if (quadruplets.size() == 0) return;
322
323 // pairs
324 std::vector<std::pair<unsigned int, unsigned int> >::iterator pairItr;
325 for (pairItr=pairIndices.begin(); pairItr!=pairIndices.end(); ++pairItr) {
326 std::pair<unsigned int, unsigned int> pair = (*pairItr);
327 Combination tmpPair;
328 std::vector<const xAOD::Muon*> theMuons = {muonsIn[pair.first],muonsIn[pair.second]};
329 tmpPair.muons = std::move(theMuons);
330 tmpPair.pairIndices = pair;
331 pairs.emplace_back(std::move(tmpPair));
332 }
333
334 return;
335
336 }
337
338 // *********************************************************************************
339
340 // ---------------------------------------------------------------------------------
341 // passesQuadSelection: 4-muon selection
342 // ---------------------------------------------------------------------------------
343
344 bool FourMuonTool::passesQuadSelection(const std::vector<const xAOD::Muon*> &muons) {
345 bool accept(false);
346 bool charges(true);
347 bool quality(false);
348 if (( muons.at(0)->muonType() == xAOD::Muon::MuonType::Combined ) ||
349 ( muons.at(1)->muonType() == xAOD::Muon::MuonType::Combined ) ||
350 ( muons.at(2)->muonType() == xAOD::Muon::MuonType::Combined ) ||
351 ( muons.at(3)->muonType() == xAOD::Muon::MuonType::Combined )
352 ) quality = true;
353 if (charges && quality) accept = true;
354 return accept;
355 }
356
357}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_DEBUG(x)
: B-physics xAOD helpers.
double length(const pvec &v)
static Double_t a
Handle class for adding a decoration to an object.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
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.
static void PrepareVertexLinks(xAOD::Vertex *theResult, const xAOD::TrackParticleContainer *importedTrackCollection)
FourMuonTool(const std::string &t, const std::string &n, const IInterface *p)
static void combinatorics(unsigned int offset, unsigned int k, std::vector< unsigned int > &combination, std::vector< unsigned int > &mainList, std::vector< std::vector< unsigned int > > &allCombinations)
static bool passesQuadSelection(const std::vector< const xAOD::Muon * > &muonsIn)
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfo_key
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_TrkParticleCollection
PublicToolHandle< Trk::ITrackSelectorTool > m_trkSelector
PublicToolHandle< Trk::IVertexFitter > m_iVertexFitter
static std::vector< std::vector< unsigned int > > getQuadIndices(unsigned int length)
SG::ReadHandleKey< xAOD::MuonContainer > m_muonCollectionKey
Gaudi::Property< bool > m_useV0Fitter
static void buildCombinations(const std::vector< const xAOD::Muon * > &muonsIn, std::vector< Combination > &pairs, std::vector< Combination > &quadruplets, unsigned int nSelectedMuons)
StatusCode performSearch(xAOD::VertexContainer *pairVxContainer, xAOD::VertexContainer *quadVxContainer, bool &acceptEvent, const EventContext &ctx) const
Gaudi::Property< double > m_etaCut
static std::vector< std::pair< unsigned int, unsigned int > > getPairIndices(unsigned int length)
SG::WriteDecorHandleKey< xAOD::MuonContainer > m_muonIndex
Gaudi::Property< double > m_ptCut
static std::vector< std::vector< unsigned int > > mFromN(unsigned int m, unsigned int n)
std::unique_ptr< xAOD::Vertex > fit(const EventContext &ctx, const std::vector< const xAOD::TrackParticle * > &, const xAOD::TrackParticleContainer *importedTrackCollection, const Amg::Vector3D &beamSpot) const
PublicToolHandle< Trk::IVertexFitter > m_iV0VertexFitter
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type get() const
Dereference the pointer, but don't cache anything.
Handle class for adding a decoration to an object.
This class implements a vertex fitting algorithm optimised for V0 finding.
virtual std::unique_ptr< xAOD::Vertex > fit(const EventContext &ctx, const std::vector< const xAOD::TrackParticle * > &vectorTrk, const Amg::Vector3D &startingPoint) const override
Interface for xAOD::TrackParticle with Amg::Vector3D starting point.
STL class.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Eigen::Matrix< double, 3, 1 > Vector3D
THE reconstruction tool.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
Muon_v1 Muon
Reference the current persistent version:
std::vector< const xAOD::Muon * > muons
std::vector< unsigned int > quadIndices
std::pair< unsigned int, unsigned int > pairIndices