ATLAS Offline Software
Loading...
Searching...
No Matches
TrigBmumuxComboHypo.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#include <algorithm>
6#include <numeric>
7#include <iterator>
8
10
11#include "xAODMuon/Muon.h"
13#include "xAODTracking/Vertex.h"
20
24
25#include "AthViews/View.h"
26#include "AthViews/ViewHelper.h"
28
29#include "Math/GenVector/VectorUtil.h"
30#include "Math/Vector2D.h"
31
32
37using ROOT::Math::XYVector;
38
39
40
41const std::vector<std::vector<double>> TrigBmumuxComboHypo::s_trkMass{
42 {PDG::mMuon, PDG::mMuon}, // {Psi.mu1, Psi.mu2}
43 {PDG::mMuon, PDG::mMuon, PDG::mKaon}, // {Psi.mu1, Psi.mu2, trk1}
44 {PDG::mMuon, PDG::mMuon, PDG::mKaon, PDG::mKaon}, // {Psi.mu1, Psi.mu2, trk1, trk2}
45 {PDG::mKaon, PDG::mKaon, PDG::mPion}, // {D_s+.K+, D_s+.K-, D_s+.pi+}
46 {PDG::mPion, PDG::mPion, PDG::mKaon}, // {D+.pi+, D+.pi+, D+.K-}
47 {PDG::mKaon, PDG::mPion}, // {D0.K-, D0.pi+}
48 {PDG::mMuon, PDG::mMuon, PDG::mPion}, // {Psi.mu1, Psi.mu2, D*+.pi+}
49 {PDG::mPion, PDG::mPion}, // {trk1, trk2}
50 {PDG::mMuon, PDG::mMuon, PDG::mPion} // {mu1, mu2, trk1}
51};
52
53TrigBmumuxComboHypo::TrigBmumuxComboHypo(const std::string& name, ISvcLocator* pSvcLocator)
54 : ::ComboHypo(name, pSvcLocator) {}
55
56
58 ATH_MSG_DEBUG( "TrigBmumuxComboHypo::initialize()" );
59
61
62 ATH_CHECK( m_muonContainerKey.initialize() );
66 ATH_CHECK( m_trigBphysContainerKey.initialize() );
67 ATH_CHECK( m_beamSpotKey.initialize() );
68
69 ATH_CHECK( m_vertexFitter.retrieve() );
71 ATH_CHECK( m_trackToVertexTool.retrieve() );
72
73 // allowed IDs to filter out incoming decisions at L2 level
74 for (const auto& item : triggerMultiplicityMap()) {
75 const HLT::Identifier id = HLT::Identifier(item.first);
76 m_allowedIDs.insert(id.numeric());
77 if (item.second.size() > 1) {
78 for (size_t i = 0; i < item.second.size(); i++) {
79 m_allowedIDs.insert(TrigCompositeUtils::createLegName(id, i).numeric());
80 }
81 }
82 }
83 if (msgLvl(MSG::DEBUG)) {
84 ATH_MSG_DEBUG( "Allowed decisions:" );
85 for (const DecisionID& id : m_allowedIDs) {
86 ATH_MSG_DEBUG( " +++ " << HLT::Identifier(id) );
87 }
88 }
89
90 if (!m_monTool.empty()) {
91 ATH_CHECK( m_monTool.retrieve() );
92 ATH_MSG_DEBUG( "GenericMonitoringTool name:" << m_monTool );
93 }
94 else {
95 ATH_MSG_DEBUG( "No GenericMonitoringTool configured: no monitoring histograms will be available" );
96 }
97
98 return StatusCode::SUCCESS;
99}
100
101
102StatusCode TrigBmumuxComboHypo::execute(const EventContext& context) const {
103
104 ATH_MSG_DEBUG( "TrigBmumuxComboHypo::execute() starts" );
105
106 ATH_MSG_DEBUG( "decision input key: " << decisionsInput().at(0).key() );
107 auto previousDecisionsHandle = SG::makeHandle(decisionsInput().at(0), context);
108 ATH_CHECK( previousDecisionsHandle.isValid() );
109 ATH_MSG_DEBUG( "Running with " << previousDecisionsHandle->size() << " previous decisions" );
110
112
113 auto trigBphysHandle = SG::makeHandle(m_trigBphysContainerKey, context);
114 ATH_CHECK( trigBphysHandle.record(std::make_unique<xAOD::TrigBphysContainer>(),
115 std::make_unique<xAOD::TrigBphysAuxContainer>()) );
116
118 ATH_CHECK( beamSpotHandle.isValid() );
119
120 auto state = std::make_unique<TrigBmumuxState>(context, *previousDecisionsHandle, *outputDecisionsHandle, trigBphysHandle.ptr(), *beamSpotHandle);
121
124
125 if (!state->dimuons.empty()) {
129 }
130
131 ATH_MSG_DEBUG( "TrigBmumuxComboHypo::execute() terminates with StatusCode::SUCCESS" );
132 return StatusCode::SUCCESS;
133}
134
135
137
138 auto& muons = state.muons;
139 muons.clear();
140
141 // all muons from views are already connected with previous decisions by TrigMuonEFHypoAlg
142 for (const Decision* decision : state.previousDecisions()) {
144 auto muonEL = decision->objectLink<xAOD::MuonContainer>(TrigCompositeUtils::featureString());
145 const xAOD::Muon* muon = *muonEL;
146 if (!muon->trackParticle(xAOD::Muon::TrackParticleType::CombinedTrackParticle)) continue;
147
148 auto decisionEL = TrigCompositeUtils::decisionToElementLink(decision, state.context());
149 auto itr = std::find_if(muons.begin(), muons.end(), [this, muon](const auto& x){ return isIdenticalTracks(muon, *x.link); });
150 if (itr == muons.end()) {
151 muons.push_back({muonEL, std::vector<ElementLink<DecisionContainer>>(1, decisionEL), DecisionIDContainer()});
152 }
153 else {
154 (*itr).decisionLinks.push_back(decisionEL);
155 }
156 }
157
158 // muon->pt() is equal to muon->trackParticle(xAOD::Muon::TrackParticleType::CombinedTrackParticle)->pt()
159 // and the later is used by TrigMuonEFHypoTool for classification of muEFCB candidates
160 std::sort(muons.begin(), muons.end(), [](const auto& lhs, const auto& rhs){ return ((*lhs.link)->pt() > (*rhs.link)->pt()); });
161
162 // for each muon we extract DecisionIDs stored in the associated Decision objects and copy them at muon.decisionIDs
163 for (auto& item : muons) {
164 for (const ElementLink<xAOD::TrigCompositeContainer>& decisionEL : item.decisionLinks) {
165 TrigCompositeUtils::decisionIDs(*decisionEL, item.decisionIDs);
166 }
167 }
168
169 if (msgLvl(MSG::DEBUG)) {
170 ATH_MSG_DEBUG( "Dump found muons before vertex fit: " << muons.size() << " candidates" );
171 for (const auto& item : muons) {
172 const xAOD::Muon* muon = *item.link;
173 const xAOD::TrackParticle* track = muon->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle);
174 ATH_MSG_DEBUG( " -- muon InDetTrackParticle pt/eta/phi/q: " << track->pt() << " / " << track->eta() << " / " << track->phi() << " / " << track->charge() );
175 ATH_MSG_DEBUG( " muon CombinedTrackParticle pt: " << muon->pt() );
176 ATH_MSG_DEBUG( " allowed decisions:" );
177 for (const DecisionID& id : item.decisionIDs) {
178 ATH_MSG_DEBUG( " +++ " << HLT::Identifier(id) );
179 }
180 }
181 }
182
183 return StatusCode::SUCCESS;
184}
185
186
188
189 auto& tracks = state.tracks;
190 tracks.clear();
191
192 size_t viewCounter = 0;
193 for (const Decision* decision : state.previousDecisions()) {
195 ATH_CHECK( viewLinkInfo.isValid() );
196 auto view = *viewLinkInfo.link;
197
199 ATH_CHECK( roiLinkInfo.isValid() );
200 const auto roi = *roiLinkInfo.link;
201
202 auto tracksHandle = ViewHelper::makeHandle(view, m_trackParticleContainerKey, state.context());
203 ATH_CHECK( tracksHandle.isValid() );
204 ATH_MSG_DEBUG( "tracks handle " << m_trackParticleContainerKey << " size: " << tracksHandle->size() );
205
206 std::vector<ElementLink<xAOD::TrackParticleContainer>> tracksFromView;
207 tracksFromView.reserve(tracksHandle->size());
208 for (size_t idx = 0; idx < tracksHandle->size(); ++idx) {
209 tracksFromView.emplace_back(ViewHelper::makeLink<xAOD::TrackParticleContainer>(view, tracksHandle, idx));
210 }
211
212 for (const auto& trackEL : tracksFromView) {
213 const xAOD::TrackParticle* track = *trackEL;
214 if (track->definingParametersCovMatrixVec().empty()) continue;
215
216 if (viewCounter == 0 ||
217 std::find_if(tracks.begin(), tracks.end(),
218 [this, track](const auto& x){ return isIdenticalTracks(track, *x); }) == tracks.end()) {
219 tracks.emplace_back(trackEL);
220 }
221 }
222 viewCounter++;
223 if (roi->composite()) {
224 state.isCompositeRoI = true;
225 break;
226 }
227 }
228 std::sort(tracks.begin(), tracks.end(), [](const auto& lhs, const auto& rhs){ return ((*lhs)->pt() > (*rhs)->pt()); });
229
230 if (msgLvl(MSG::DEBUG)) {
231 ATH_MSG_DEBUG( "Dump found tracks before vertex fit: " << tracks.size() << " candidates" );
232 for (const auto& trackEL : tracks) {
233 const xAOD::TrackParticle* track = *trackEL;
234 ATH_MSG_DEBUG( " -- track pt/eta/phi/q: " << track->pt() << " / " << track->eta() << " / " << track->phi() << " / " << track->charge() );
235 }
236 }
237 return StatusCode::SUCCESS;
238}
239
240
242
243 auto mon_nDimuon = Monitored::Scalar<int>("nDimuon", 0);
244 auto group = Monitored::Group(m_monTool, mon_nDimuon);
245
246 const auto& muons = state.muons;
247 std::vector<ElementLink<xAOD::TrackParticleContainer>> trackParticleLinks(2);
248 std::vector<const DecisionIDContainer*> previousDecisionIDs(2, nullptr);
249
250 for (size_t itrk1 = 0; itrk1 < muons.size(); ++itrk1) {
251 const xAOD::Muon* mu1 = *muons[itrk1].link;
252 trackParticleLinks[0] = ::linkTrack(mu1->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle));
253 previousDecisionIDs[0] = &muons[itrk1].decisionIDs;
254 const xAOD::TrackParticle* trk1 = *trackParticleLinks[0];
255 auto p1 = trk1->genvecP4();
256 p1.SetM(PDG::mMuon);
257 auto charge1 = trk1->charge();
258
259 for (size_t itrk2 = itrk1 + 1; itrk2 < muons.size(); ++itrk2) {
260 const xAOD::Muon* mu2 = *muons[itrk2].link;
261 trackParticleLinks[1] = ::linkTrack(mu2->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle));
262 previousDecisionIDs[1] = &muons[itrk2].decisionIDs;
263 const xAOD::TrackParticle* trk2 = *trackParticleLinks[1];
264 auto p2 = trk2->genvecP4();
265 p2.SetM(PDG::mMuon);
266 auto charge2 = trk2->charge();
267
268 double mass = (p1 + p2).M();
269
270 ATH_MSG_DEBUG( "muon 1: " << p1.Pt()<< " / " << p1.Eta() << " / " << p1.Phi() << " / " << trk1->charge() );
271 ATH_MSG_DEBUG( "muon 2: " << p2.Pt()<< " / " << p2.Eta() << " / " << p2.Phi() << " / " << trk2->charge() );
272 ATH_MSG_DEBUG( "track pair mass: " << mass );
273
274 if (m_dimuon_rejectSameChargeTracks && charge1 * charge2 > 0.) {
275 ATH_MSG_DEBUG( "muon pair is rejected by opposite charge check" );
276 continue;
277 }
278
279 if (!passDimuonTrigger(previousDecisionIDs)) {
280 ATH_MSG_DEBUG( "muon pair did not pass passDimuonTrigger() check" );
281 continue;
282 }
283
284 if (!isInMassRange(mass, m_dimuon_massRange)) {
285 ATH_MSG_DEBUG( "muon pair is out of the requested mass range" );
286 continue;
287 }
288
289 // fit muons to the common vertex and pass this vertex to the dimuons collection which also takes the ownership of the created object
290 auto vertex = fit(state.context(), trackParticleLinks);
291 if (!vertex) continue;
292
293 // convert vertex to trigger object and add it to the output xAOD::TrigBphysContainer
294 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vertex);
295 if (!trigBphys) {
296 ATH_MSG_ERROR( "xAOD::Vertex could not be converted to xAOD::TrigBphys object: please enable MakeExtendedVertex option in vertex fitter " << m_vertexFitter->name() );
297 return StatusCode::FAILURE;
298 }
299
300 // the dimuon vertex in the xAOD::TrigBphysContainer is used as a reference for the upcoming B-candidates and to fire bBmumux_idperf chain
301 state.dimuons.push_back(vertex.release());
302 state.trigBphysMuonIndices.emplace_back(std::array<size_t, 2>{itrk1, itrk2});
303 state.trigBphysCollection().push_back(trigBphys);
304 }
305 }
306 mon_nDimuon = state.dimuons.size();
307 ATH_MSG_DEBUG( "Found " << state.dimuons.size() << " dimuon candidates" );
308
309 return StatusCode::SUCCESS;
310}
311
312
314
315 std::vector<int> nSelectedTrk;
316 auto mon_nTrk = Monitored::Scalar<int>("nTrk", 0);
317 auto mon_nSelectedTrk = Monitored::Collection("nSelectedTrk", nSelectedTrk);
318 auto mon_nBPhysObject = Monitored::Scalar<int>("nBPhysObject", 0);
319 auto group = Monitored::Group(m_monTool, mon_nTrk, mon_nSelectedTrk, mon_nBPhysObject);
320
321 for (size_t dimuonIndex = 0; dimuonIndex < state.dimuons.size(); ++dimuonIndex) {
322
323 ATH_CHECK( findBmumuxCandidates_selectTracks(state, dimuonIndex) );
324 if (state.selectedTracks.empty()) continue;
325 nSelectedTrk.push_back(state.selectedTracks.size());
326
327 ATH_CHECK( findBmumuxCandidates_fit(state, dimuonIndex, true) );
328 ATH_CHECK( findBmumuxCandidates_fastFit(state, dimuonIndex) );
329 ATH_CHECK( findBmumuxCandidates_fit(state, dimuonIndex) );
330 }
331
332 mon_nTrk = state.tracks.size();
333 mon_nBPhysObject = state.trigBphysCollection().size() - state.dimuons.size();
334
335 return StatusCode::SUCCESS;
336}
337
338
340
341 auto& selectedTracks = state.selectedTracks;
342 auto& selectedTrackZ0 = state.selectedTrackZ0;
343
344 selectedTracks.clear();
345 selectedTrackZ0.clear();
346
347 const xAOD::Vertex* dimuon = state.dimuons.get(dimuonIndex);
348
349 std::vector<const xAOD::Muon*> muons(2, nullptr);
350 const auto& muonIndices = state.trigBphysMuonIndices.at(dimuonIndex);
351 for (size_t i = 0; i < 2; ++i) {
352 const auto& muon = state.muons.at(muonIndices[i]);
353 muons[i] = *muon.link;
354 }
355
356 const xAOD::TrackParticle* mu1 = dimuon->trackParticle(0);
357 const xAOD::TrackParticle* mu2 = dimuon->trackParticle(1);
358
359 // check impact parameter of the track with respect to the fitted dimuon vertex
360 // we can safely omit tracks with large z0
361 state.selectedTracks.reserve(state.tracks.size());
362 for (const auto& trackEL : state.tracks) {
363 if (state.isCompositeRoI && !isInSameRoI(muons[0], *trackEL) && !isInSameRoI(muons[1], *trackEL)) continue;
364 if (m_trkZ0 > 0.) {
365 std::unique_ptr<const Trk::Perigee> perigee(m_trackToVertexTool->perigeeAtVertex(state.context(), **trackEL, dimuon->position()));
366 if (perigee && std::abs(perigee->parameters()[Trk::z0]) < m_trkZ0) {
367 selectedTracks.push_back(trackEL);
368 selectedTrackZ0[*trackEL] = perigee->parameters()[Trk::z0];
369 }
370 }
371 else {
372 selectedTracks.push_back(trackEL);
373 selectedTrackZ0[*trackEL] = -1000.;
374 }
375 }
376
377 // remove muon duplicates
378 if (selectedTracks.size() < 2) {
379 ATH_MSG_DEBUG( "Found no tracks consistent with dimuon vertex " << dimuonIndex );
380 selectedTracks.clear();
381 selectedTrackZ0.clear();
382 return StatusCode::SUCCESS;
383 }
384 std::sort(selectedTracks.begin(), selectedTracks.end(), [p_mu=mu1->genvecP4()](const auto& lhs, const auto& rhs){ return ROOT::Math::VectorUtil::DeltaR(p_mu, (*lhs)->genvecP4()) > ROOT::Math::VectorUtil::DeltaR(p_mu, (*rhs)->genvecP4()); });
385 if (isIdenticalTracks(mu1, *selectedTracks.back())) selectedTracks.pop_back();
386 std::sort(selectedTracks.begin(), selectedTracks.end(), [p_mu=mu2->genvecP4()](const auto& lhs, const auto& rhs){ return ROOT::Math::VectorUtil::DeltaR(p_mu, (*lhs)->genvecP4()) > ROOT::Math::VectorUtil::DeltaR(p_mu, (*rhs)->genvecP4()); });
387 if (isIdenticalTracks(mu2, *selectedTracks.back())) selectedTracks.pop_back();
388 std::sort(selectedTracks.begin(), selectedTracks.end(), [](const auto& lhs, const auto& rhs){ return ((*lhs)->pt() > (*rhs)->pt()); });
389
390 ATH_MSG_DEBUG( "Found " << selectedTracks.size() << " tracks consistent with dimuon vertex " << dimuonIndex );
391
392 return StatusCode::SUCCESS;
393}
394
395
396StatusCode TrigBmumuxComboHypo::findBmumuxCandidates_fit(TrigBmumuxState& state, size_t dimuonIndex, bool makeCombinations) const {
397
398 const auto& selectedTracks = state.selectedTracks;
399 if (makeCombinations) state.trackCombinations.clear();
400
401 const xAOD::Vertex* dimuon = state.dimuons.get(dimuonIndex);
402 auto dimuonTriggerObjectEL = ElementLink<xAOD::TrigBphysContainer>(state.trigBphysCollection(), dimuonIndex);
403 ATH_CHECK( dimuonTriggerObjectEL.isValid() );
404
405 // vtx1 = {mu1, mu2, trk1}
406 std::vector<ElementLink<xAOD::TrackParticleContainer>> trackParticleLinks_vtx1(dimuon->trackParticleLinks());
407 trackParticleLinks_vtx1.emplace_back();
408
409 // vtx2 = {mu1, mu2, trk1, trk2}
410 std::vector<ElementLink<xAOD::TrackParticleContainer>> trackParticleLinks_vtx2(trackParticleLinks_vtx1);
411 trackParticleLinks_vtx2.emplace_back();
412
413 // vtx3 = {trk1, trk2, trk3}
414 std::vector<ElementLink<xAOD::TrackParticleContainer>> trackParticleLinks_vtx3(3);
415
416 const xAOD::TrackParticle* mu1 = *trackParticleLinks_vtx1[0];
417 const xAOD::TrackParticle* mu2 = *trackParticleLinks_vtx1[1];
418 auto p_dimuon = mu1->genvecP4().SetM(PDG::mMuon) + mu2->genvecP4().SetM(PDG::mMuon);
419
420 size_t iterations = 0;
421 size_t nTrigBphysObjects = state.trigBphysCollection().size();
422 bool isOverWarningThreshold = false;
423 // dimuon + 1 track
424 for (size_t itrk1 = 0; itrk1 < selectedTracks.size(); ++itrk1) {
425 const xAOD::TrackParticle* trk1 = *selectedTracks[itrk1];
426
427 trackParticleLinks_vtx1[2] = selectedTracks[itrk1];
428 auto p_trk1 = trk1->genvecP4();
429 auto charge1 = trk1->charge();
430
431 std::unique_ptr<xAOD::Vertex> vtx1;
432 bool makeFit_vtx1 = !makeCombinations;
433 bool passFastFit_vtx1 = (!makeCombinations && !state.isBadCombination(itrk1));
434
435 // B+ -> mu+ mu- K+
436 if (m_BplusToMuMuKaon &&
437 p_trk1.Pt() > m_BplusToMuMuKaon_minKaonPt &&
438 isInMassRange((p_dimuon + p_trk1.SetM(PDG::mKaon)).M(), m_BplusToMuMuKaon_massRange)) {
439
440 if (makeCombinations && m_BplusToMuMuKaon_useFastFit) state.addTrackCombination(itrk1);
441 if (!vtx1 && makeFit_vtx1 && (passFastFit_vtx1 || !m_BplusToMuMuKaon_useFastFit)) {
442 vtx1 = fit(state.context(), trackParticleLinks_vtx1, kB_2mu1trk, dimuon);
443 makeFit_vtx1 = false;
444 ++iterations;
445 }
446 if (vtx1 && vtx1->chiSquared() < m_BplusToMuMuKaon_chi2) {
447 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx1, xAOD::TrigBphys::BKMUMU, {PDG::mMuon, PDG::mMuon, PDG::mKaon}, dimuonTriggerObjectEL);
448 ATH_CHECK( state.addTriggerObject(trigBphys) );
449 }
450 }
451
452 // B_c+ -> J/psi(-> mu+ mu-) pi+
453 if (m_BcToMuMuPion &&
454 p_trk1.Pt() > m_BcToMuMuPion_minPionPt &&
456 isInMassRange((p_dimuon + p_trk1.SetM(PDG::mPion)).M() - p_dimuon.M() + PDG::mJpsi, m_BcToMuMuPion_massRange)) {
457
458 if (makeCombinations && m_BcToMuMuPion_useFastFit) state.addTrackCombination(itrk1);
459 if (!vtx1 && makeFit_vtx1 && (passFastFit_vtx1 || !m_BcToMuMuPion_useFastFit)) {
460 vtx1 = fit(state.context(), trackParticleLinks_vtx1, kB_2mu1trk, dimuon);
461 makeFit_vtx1 = false;
462 ++iterations;
463 }
464 if (vtx1 && vtx1->chiSquared() < m_BcToMuMuPion_chi2) {
465 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx1, xAOD::TrigBphys::BCPIMUMU, {PDG::mMuon, PDG::mMuon, PDG::mPion}, dimuonTriggerObjectEL);
466 ATH_CHECK( state.addTriggerObject(trigBphys) );
467 }
468 }
469 vtx1.reset();
470
471 // dimuon + 2 tracks
472 for (size_t itrk2 = itrk1 + 1; itrk2 < selectedTracks.size(); ++itrk2) {
473 const xAOD::TrackParticle* trk2 = *selectedTracks[itrk2];
474
475 trackParticleLinks_vtx2[2] = selectedTracks[itrk1];
476 trackParticleLinks_vtx2[3] = selectedTracks[itrk2];
477 auto p_trk2 = trk2->genvecP4();
478 auto charge2 = trk2->charge();
479
480 std::unique_ptr<xAOD::Vertex> vtx2;
481 bool makeFit_vtx2 = !makeCombinations;
482 bool passFastFit_vtx2 = (!makeCombinations && !state.isBadCombination(itrk1) && !state.isBadCombination(itrk2));
483
484 // B_s0 -> mu+ mu- phi(-> K+ K-)
485 if (m_BsToMuMuPhi1020 &&
486 (!m_BsToMuMuPhi1020_rejectSameChargeTracks || charge1 * charge2 < 0.) &&
487 p_trk1.Pt() > m_BsToMuMuPhi1020_minKaonPt &&
488 p_trk2.Pt() > m_BsToMuMuPhi1020_minKaonPt &&
489 isInMassRange((p_trk1.SetM(PDG::mKaon) + p_trk2.SetM(PDG::mKaon)).M(), m_BsToMuMuPhi1020_phiMassRange) &&
490 isInMassRange((p_dimuon + p_trk1.SetM(PDG::mKaon) + p_trk2.SetM(PDG::mKaon)).M(), m_BsToMuMuPhi1020_massRange)) {
491
492 if (makeCombinations && m_BsToMuMuPhi1020_useFastFit) {
493 state.addTrackCombination(itrk1);
494 state.addTrackCombination(itrk2);
495 }
496 if (!vtx2 && makeFit_vtx2 && (passFastFit_vtx2 || !m_BsToMuMuPhi1020_useFastFit)) {
497 vtx2 = fit(state.context(), trackParticleLinks_vtx2, kB_2mu2trk, dimuon);
498 makeFit_vtx2 = false;
499 ++iterations;
500 }
501 if (vtx2 && vtx2->chiSquared() < m_BsToMuMuPhi1020_chi2) {
502 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx2, xAOD::TrigBphys::BSPHIMUMU, {PDG::mMuon, PDG::mMuon, PDG::mKaon, PDG::mKaon}, dimuonTriggerObjectEL);
503 ATH_CHECK( state.addTriggerObject(trigBphys) );
504 }
505 }
506
507 // B0 -> mu+ mu- K*0(-> K+ pi-)
508 if (m_BdToMuMuKstar0 &&
509 (!m_BdToMuMuKstar0_rejectSameChargeTracks || charge1 * charge2 < 0.) &&
510 p_trk1.Pt() > m_BdToMuMuKstar0_minKaonPt &&
511 p_trk2.Pt() > m_BdToMuMuKstar0_minPionPt &&
512 isInMassRange((p_trk1.SetM(PDG::mKaon) + p_trk2.SetM(PDG::mPion)).M(), m_BdToMuMuKstar0_KstarMassRange) &&
513 isInMassRange((p_dimuon + p_trk1.SetM(PDG::mKaon) + p_trk2.SetM(PDG::mPion)).M(), m_BdToMuMuKstar0_massRange)) {
514
515 if (makeCombinations && m_BdToMuMuKstar0_useFastFit) {
516 state.addTrackCombination(itrk1);
517 state.addTrackCombination(itrk2);
518 }
519 if (!vtx2 && makeFit_vtx2 && (passFastFit_vtx2 || !m_BdToMuMuKstar0_useFastFit)) {
520 vtx2 = fit(state.context(), trackParticleLinks_vtx2, kB_2mu2trk, dimuon);
521 makeFit_vtx2 = false;
522 ++iterations;
523 }
524 if (vtx2 && vtx2->chiSquared() < m_BdToMuMuKstar0_chi2) {
525 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx2, xAOD::TrigBphys::BDKSTMUMU, {PDG::mMuon, PDG::mMuon, PDG::mKaon, PDG::mPion}, dimuonTriggerObjectEL);
526 ATH_CHECK( state.addTriggerObject(trigBphys) );
527 }
528 }
529 // anti-B0 -> mu+ mu- anti-K*0(-> K- pi+)
530 if (m_BdToMuMuKstar0 &&
531 (!m_BdToMuMuKstar0_rejectSameChargeTracks || charge1 * charge2 < 0.) &&
532 p_trk1.Pt() > m_BdToMuMuKstar0_minPionPt &&
533 p_trk2.Pt() > m_BdToMuMuKstar0_minKaonPt &&
534 isInMassRange((p_trk1.SetM(PDG::mPion) + p_trk2.SetM(PDG::mKaon)).M(), m_BdToMuMuKstar0_KstarMassRange) &&
535 isInMassRange((p_dimuon + p_trk1.SetM(PDG::mPion) + p_trk2.SetM(PDG::mKaon)).M(), m_BdToMuMuKstar0_massRange)) {
536
537 if (makeCombinations && m_BdToMuMuKstar0_useFastFit) {
538 state.addTrackCombination(itrk1);
539 state.addTrackCombination(itrk2);
540 }
541 if (!vtx2 && makeFit_vtx2 && (passFastFit_vtx2 || !m_BdToMuMuKstar0_useFastFit)) {
542 vtx2 = fit(state.context(), trackParticleLinks_vtx2, kB_2mu2trk, dimuon);
543 makeFit_vtx2 = false;
544 ++iterations;
545 }
546 if (vtx2 && vtx2->chiSquared() < m_BdToMuMuKstar0_chi2) {
547 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx2, xAOD::TrigBphys::BDKSTMUMU, {PDG::mMuon, PDG::mMuon, PDG::mPion, PDG::mKaon}, dimuonTriggerObjectEL);
548 ATH_CHECK( state.addTriggerObject(trigBphys) );
549 }
550 }
551
552 // Lambda_b0 -> J/psi(-> mu+ mu-) p K-
556 (p_trk1.SetM(PDG::mKaon) + p_trk2.SetM(PDG::mPion)).M() > m_LambdaBToMuMuProtonKaon_minKstarMass &&
557 (p_trk1.SetM(PDG::mPion) + p_trk2.SetM(PDG::mKaon)).M() > m_LambdaBToMuMuProtonKaon_minKstarMass &&
559 isInMassRange((p_dimuon + p_trk1.SetM(PDG::mProton) + p_trk2.SetM(PDG::mKaon)).M() - p_dimuon.M() + PDG::mJpsi, m_LambdaBToMuMuProtonKaon_massRange)) {
560
561 if (makeCombinations && m_LambdaBToMuMuProtonKaon_useFastFit) {
562 state.addTrackCombination(itrk1);
563 state.addTrackCombination(itrk2);
564 }
565 if (!vtx2 && makeFit_vtx2 && (passFastFit_vtx2 || !m_LambdaBToMuMuProtonKaon_useFastFit)) {
566 vtx2 = fit(state.context(), trackParticleLinks_vtx2, kB_2mu2trk, dimuon);
567 makeFit_vtx2 = false;
568 ++iterations;
569 }
570 if (vtx2 && vtx2->chiSquared() < m_LambdaBToMuMuProtonKaon_chi2 && Lxy(state.beamSpotPosition(), *vtx2) > 0.) {
571 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx2, xAOD::TrigBphys::LBPQMUMU, {PDG::mMuon, PDG::mMuon, PDG::mProton, PDG::mKaon}, dimuonTriggerObjectEL);
572 ATH_CHECK( state.addTriggerObject(trigBphys) );
573 }
574 }
575 // anti-Lambda_b0 -> J/psi(-> mu+ mu-) anti-p K+
579 (p_trk1.SetM(PDG::mKaon) + p_trk2.SetM(PDG::mPion)).M() > m_LambdaBToMuMuProtonKaon_minKstarMass &&
580 (p_trk1.SetM(PDG::mPion) + p_trk2.SetM(PDG::mKaon)).M() > m_LambdaBToMuMuProtonKaon_minKstarMass &&
582 isInMassRange((p_dimuon + p_trk1.SetM(PDG::mKaon) + p_trk2.SetM(PDG::mProton)).M() - p_dimuon.M() + PDG::mJpsi, m_LambdaBToMuMuProtonKaon_massRange)) {
583
584 if (makeCombinations && m_LambdaBToMuMuProtonKaon_useFastFit) {
585 state.addTrackCombination(itrk1);
586 state.addTrackCombination(itrk2);
587 }
588 if (!vtx2 && makeFit_vtx2 && (passFastFit_vtx2 || !m_LambdaBToMuMuProtonKaon_useFastFit)) {
589 vtx2 = fit(state.context(), trackParticleLinks_vtx2, kB_2mu2trk, dimuon);
590 makeFit_vtx2 = false;
591 ++iterations;
592 }
593 if (vtx2 && vtx2->chiSquared() < m_LambdaBToMuMuProtonKaon_chi2 && Lxy(state.beamSpotPosition(), *vtx2) > 0.) {
594 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx2, xAOD::TrigBphys::LBPQMUMU, {PDG::mMuon, PDG::mMuon, PDG::mKaon, PDG::mProton}, dimuonTriggerObjectEL);
595 ATH_CHECK( state.addTriggerObject(trigBphys) );
596 }
597 }
598 vtx2.reset();
599
600 for (size_t itrk3 = 0; itrk3 < selectedTracks.size(); ++itrk3) {
601 const xAOD::TrackParticle* trk3 = *selectedTracks[itrk3];
602 if (itrk3 == itrk1 || itrk3 == itrk2) continue;
603
604 trackParticleLinks_vtx3[0] = selectedTracks[itrk1];
605 trackParticleLinks_vtx3[1] = selectedTracks[itrk2];
606 trackParticleLinks_vtx3[2] = selectedTracks[itrk3];
607 auto p_trk3 = trk3->genvecP4();
608 auto charge3 = trk3->charge();
609
610 std::unique_ptr<xAOD::Vertex> vtx3;
611 bool makeFit_vtx3 = !makeCombinations;
612 bool passFastFit_vtx3 = (!makeCombinations && !state.isBadCombination(itrk1, itrk2, itrk3));
613
614 // B_c+ -> J/psi(-> mu+ mu-) D_s+(->phi(-> K+ K-) pi+)
615 p_trk1.SetM(PDG::mKaon); // D_s+.phi.K+
616 p_trk2.SetM(PDG::mKaon); // D_s+.phi.K-
617 p_trk3.SetM(PDG::mPion); // D_s+.pi+
618 if (m_BcToDsMuMu &&
619 charge1 * charge2 < 0. &&
620 p_trk1.Pt() > m_BcToDsMuMu_minKaonPt &&
621 p_trk2.Pt() > m_BcToDsMuMu_minKaonPt &&
622 p_trk3.Pt() > m_BcToDsMuMu_minPionPt &&
624 isInMassRange((p_trk1 + p_trk2).M(), m_BcToDsMuMu_phiMassRange) &&
625 isInMassRange((p_trk1 + p_trk2 + p_trk3).M(), m_BcToDsMuMu_DsMassRange) &&
626 isInMassRange((p_dimuon + p_trk1 + p_trk2 + p_trk3).M() - p_dimuon.M() + PDG::mJpsi, m_BcToDsMuMu_massRange)) {
627
628 if (makeCombinations && m_BcToDsMuMu_useFastFit) state.addTrackCombination(itrk1, itrk2, itrk3);
629 if (!vtx3 && makeFit_vtx3 && (passFastFit_vtx3 || !m_BcToDsMuMu_useFastFit)) {
630 vtx3 = fit(state.context(), trackParticleLinks_vtx3, kDs, dimuon);
631 makeFit_vtx3 = false;
632 ++iterations;
633 }
634 if (vtx3 && vtx3->chiSquared() < m_BcToDsMuMu_chi2) {
635 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx3, xAOD::TrigBphys::BCDSMUMU, {PDG::mKaon, PDG::mKaon, PDG::mPion}, dimuonTriggerObjectEL);
636 ATH_CHECK( state.addTriggerObject(trigBphys) );
637 }
638 }
639
640 // B_c+ -> J/psi(-> mu+ mu-) D+(-> K- pi+ pi+)
641 p_trk1.SetM(PDG::mPion); // D+.pi+
642 p_trk2.SetM(PDG::mPion); // D+.pi+
643 p_trk3.SetM(PDG::mKaon); // D+.K-
644 if (m_BcToDplusMuMu &&
645 charge1 * charge2 > 0. && charge1 * charge3 < 0. &&
646 p_trk1.Pt() > m_BcToDplusMuMu_minPionPt &&
647 p_trk2.Pt() > m_BcToDplusMuMu_minPionPt &&
648 p_trk3.Pt() > m_BcToDplusMuMu_minKaonPt &&
650 isInMassRange((p_trk1 + p_trk2 + p_trk3).M(), m_BcToDplusMuMu_DplusMassRange) &&
651 isInMassRange((p_dimuon + p_trk1 + p_trk2 + p_trk3).M() - p_dimuon.M() + PDG::mJpsi, m_BcToDplusMuMu_massRange)) {
652
653 if (makeCombinations && m_BcToDplusMuMu_useFastFit) state.addTrackCombination(itrk1, itrk2, itrk3);
654 if (!vtx3 && makeFit_vtx3 && (passFastFit_vtx3 || !m_BcToDplusMuMu_useFastFit)) {
655 vtx3 = fit(state.context(), trackParticleLinks_vtx3, kDplus, dimuon);
656 makeFit_vtx3 = false;
657 ++iterations;
658 }
659 if (vtx3 && vtx3->chiSquared() < m_BcToDplusMuMu_chi2 && Lxy(dimuon->position(), *vtx3) > 0.) {
660 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *vtx3, xAOD::TrigBphys::BCDPMUMU, {PDG::mPion, PDG::mPion, PDG::mKaon}, dimuonTriggerObjectEL);
661 ATH_CHECK( state.addTriggerObject(trigBphys) );
662 }
663 }
664 vtx3.reset();
665
666 }
667 }
668
669 if (iterations > m_fitAttemptsWarningThreshold && !isOverWarningThreshold) {
670 ATH_MSG_WARNING( "Dimuon + tracks: " << state.trigBphysCollection().size() - nTrigBphysObjects << " vertices created after " << iterations << " vertex fitter calls" );
671 isOverWarningThreshold = true;
672 }
673 if (iterations > m_fitAttemptsBreakThreshold) {
674 ATH_MSG_WARNING( "Dimuon + tracks: the number of fit attempts has exceeded the limit; breaking the loop at this point" );
675 break;
676 }
677 }
678 ATH_MSG_DEBUG( "Dimuon + tracks: " << state.trigBphysCollection().size() - nTrigBphysObjects << " vertices created after " << iterations << " vertex fitter calls" );
679
680 iterations = 0;
681 nTrigBphysObjects = state.trigBphysCollection().size();
682 isOverWarningThreshold = false;
683 // B_c+ -> J/psi(-> mu+ mu-) D*+(-> D0(-> K- pi+) pi+)
684 if (!makeCombinations && m_BcToDstarMuMu && isInMassRange(p_dimuon.M(), m_BcToDstarMuMu_dimuonMassRange)) {
685 std::vector<ElementLink<xAOD::TrackParticleContainer>> trackParticleLinks_D0(2); // {K-, pi+}
686
687 for (size_t itrk1 = 0; itrk1 < selectedTracks.size(); ++itrk1) {
688 const xAOD::TrackParticle* trk1 = *selectedTracks[itrk1];
689
690 trackParticleLinks_D0[0] = selectedTracks[itrk1];
691 auto p_trk1 = trk1->genvecP4();
692 p_trk1.SetM(PDG::mKaon);
693 auto charge1 = trk1->charge();
694
695 for (size_t itrk2 = 0; itrk2 < selectedTracks.size(); ++itrk2) {
696 if (itrk2 == itrk1) continue;
697 const xAOD::TrackParticle* trk2 = *selectedTracks[itrk2];
698
699 trackParticleLinks_D0[1] = selectedTracks[itrk2];
700 auto p_trk2 = trk2->genvecP4();
701 p_trk2.SetM(PDG::mPion);
702 auto charge2 = trk2->charge();
703
704 std::unique_ptr<xAOD::Vertex> D0;
705 if (charge1 * charge2 < 0. &&
706 p_trk1.Pt() > m_BcToDstarMuMu_minD0KaonPt &&
707 p_trk2.Pt() > m_BcToDstarMuMu_minD0PionPt &&
708 isInMassRange((p_trk1 + p_trk2).M(), m_BcToDstarMuMu_D0MassRange) &&
709 isInMassRange((p_dimuon + p_trk1 + p_trk2).M() - p_dimuon.M() + PDG::mJpsi, m_BcToDstarMuMu_massRange)) {
710 D0 = fit(state.context(), trackParticleLinks_D0, kD0, dimuon);
711 ++iterations;
712 }
713 bool isValidD0 = false;
714 if (D0 && D0->chiSquared() < m_BcToDstarMuMu_chi2 && Lxy(dimuon->position(), *D0) > 0.) {
715 isValidD0 = true;
716 ATH_MSG_DEBUG( "Partially reconstructed B_c+(-> mu+ mu- D0 X) candidate has been created" );
717 xAOD::TrigBphys* trigBphys = makeTriggerObject(state, *D0, xAOD::TrigBphys::DZKPI, s_trkMass[kD0], dimuonTriggerObjectEL);
718 ATH_CHECK( state.addTriggerObject(trigBphys) );
719 }
720
721 if (m_BcToDstarMuMu_makeDstar && isValidD0) { // full B_c+ reconstruction
723
724 for (size_t itrk3 = 0; itrk3 < selectedTracks.size(); ++itrk3) {
725 const xAOD::TrackParticle* trk3 = *selectedTracks[itrk3];
726 if (itrk3 == itrk1 || itrk3 == itrk2) continue;
727
728 // J/psi + pion from D*+
729 trackParticleLinks_vtx1[2] = selectedTracks[itrk3];
730 auto p_trk3 = trk3->genvecP4();
731 p_trk3.SetM(PDG::mPion);
732
733 if (p_trk3.Pt() > m_BcToDstarMuMu_minDstarPionPt &&
735 isInMassRange((p_D0 + p_trk3).M() - p_D0.M() + PDG::mD0, m_BcToDstarMuMu_DstarMassRange)) {
736 auto Bc_vtx1 = fit(state.context(), trackParticleLinks_vtx1, kB_PsiPi);
737 ++iterations;
738
739 if (Bc_vtx1 && Bc_vtx1->chiSquared() < m_BcToDstarMuMu_chi2) {
740 ATH_MSG_DEBUG( "Decay vertex(mu+ mu- D*+.pi+) for B_c+ candidate has been created" );
741 xAOD::TrigBphys* triggerObject_vtx1 = makeTriggerObject(state, *Bc_vtx1, xAOD::TrigBphys::DSTDZPI, s_trkMass[kB_PsiPi], dimuonTriggerObjectEL);
742 auto triggerObjectEL_vtx1 = ElementLink<xAOD::TrigBphysContainer>(state.trigBphysCollection(), state.trigBphysCollection().size());
743 ATH_CHECK( state.addTriggerObject(triggerObject_vtx1) );
744 ATH_CHECK( triggerObjectEL_vtx1.isValid() );
745
746 // refit D0 vertex
747 auto Bc_vtx2 = fit(state.context(), trackParticleLinks_D0, kD0, Bc_vtx1.get());
748 ++iterations;
749 if (Bc_vtx2 && Bc_vtx2->chiSquared() < m_BcToDstarMuMu_chi2) {
750 ATH_MSG_DEBUG( "Fully reconstructed B_c+(-> mu+ mu- D*+) candidate has been created" );
751 xAOD::TrigBphys* triggerObject_vtx2 = makeTriggerObject(state, *Bc_vtx2, xAOD::TrigBphys::BCDSTMUMU, s_trkMass[kD0], triggerObjectEL_vtx1);
752 ATH_CHECK( state.addTriggerObject(triggerObject_vtx2) );
753 }
754 }
755 }
756 }
757 } // end of full B_c+ reconstruction
758
759 }
760
761 if (iterations > m_fitAttemptsWarningThreshold && !isOverWarningThreshold) {
762 ATH_MSG_WARNING( "B_c+ -> mu+ mu- D*+: " << state.trigBphysCollection().size() - nTrigBphysObjects << " vertices created after " << iterations << " vertex fitter calls" );
763 isOverWarningThreshold = true;
764 }
765 if (iterations > m_fitAttemptsBreakThreshold) {
766 ATH_MSG_WARNING( "B_c+ -> mu+ mu- D*+: the number of fit attempts has exceeded the limit; breaking the loop at this point" );
767 break;
768 }
769 }
770 ATH_MSG_DEBUG( "B_c+ -> mu+ mu- D*+: " << state.trigBphysCollection().size() - nTrigBphysObjects << " vertices created after " << iterations << " vertex fitter calls" );
771
772 } // end of B_c+ -> J/psi D*+ topology
773
774 return StatusCode::SUCCESS;
775}
776
777
778StatusCode TrigBmumuxComboHypo::findBmumuxCandidates_fastFit(TrigBmumuxState& state, size_t dimuonIndex) const {
779
780 state.badTrackCombinations.clear();
781
782 const xAOD::Vertex* dimuon = state.dimuons.get(dimuonIndex);
783
784 // {mu1, mu2, trk1}
785 std::vector<ElementLink<xAOD::TrackParticleContainer>> trackParticleLinks_2mu1trk(dimuon->trackParticleLinks());
786 trackParticleLinks_2mu1trk.emplace_back();
787
788 // {trk1, trk2}
789 std::vector<ElementLink<xAOD::TrackParticleContainer>> trackParticleLinks_2trk(2);
790
791 size_t n = state.selectedTracks.size();
792 size_t iterations = 0;
793 for (const auto& item : state.trackCombinations) {
794 if (item.second < 5) continue;
795
796 size_t key = item.first;
797 if (key < n) { // dimuon + track
798 trackParticleLinks_2mu1trk[2] = state.selectedTracks[key];
799 auto vertex = fit(state.context(), trackParticleLinks_2mu1trk, kFastFit_2mu1trk, dimuon);
800 iterations++;
801 if (!vertex || vertex->chiSquared() > m_fastFit_2mu1trk_chi2) state.badTrackCombinations.push_back(key);
802 }
803 else { // track + track
804 trackParticleLinks_2trk[0] = state.selectedTracks[key % n];
805 trackParticleLinks_2trk[1] = state.selectedTracks[key / n];
806 auto vertex = fit(state.context(), trackParticleLinks_2trk, kFastFit_2trk);
807 iterations++;
808 if (!vertex || vertex->chiSquared() > m_fastFit_2trk_chi2) state.badTrackCombinations.push_back(key);
809 }
810 }
811 std::sort(state.badTrackCombinations.begin(), state.badTrackCombinations.end());
812 ATH_MSG_DEBUG( "Fast fit found " << state.badTrackCombinations.size() << " bad combinations after " << iterations << " iterations" );
813
814 return StatusCode::SUCCESS;
815}
816
817
819
820 for (const xAOD::TrigBphys* triggerObject : state.trigBphysCollection()) {
821 // skip the first vertex of the cascade decay B_c+(-> mu+ mu- D*+), it is already linked to the xAOD::TrigBphys::BCDSTMUMU trigger object via lowerChainLink()
822 if (triggerObject->particleType() == xAOD::TrigBphys::DSTDZPI) continue;
823
824 ATH_MSG_DEBUG( "Found xAOD::TrigBphys object: mass = " << triggerObject->mass() );
825
826 auto triggerObjectEL = ElementLink<xAOD::TrigBphysContainer>(state.trigBphysCollection(), triggerObject->index());
827 ATH_CHECK( triggerObjectEL.isValid() );
828
829 const xAOD::TrigBphys* dimuonTriggerObject = (triggerObject->particleType() == xAOD::TrigBphys::MULTIMU ? triggerObject : triggerObject->lowerChain());
830 if (triggerObject->particleType() == xAOD::TrigBphys::BCDSTMUMU && dimuonTriggerObject) dimuonTriggerObject = dimuonTriggerObject->lowerChain();
831 if (!dimuonTriggerObject) {
832 ATH_MSG_ERROR( "Failed to found a valid link for preceding dimuon trigger object" );
833 return StatusCode::FAILURE;
834 }
835
836 // need to get the references to the original muon objects used to build the dimuon vertex
837 // the position of this vertex in state.dimuons container is the same as for dimuonTriggerObject in trigBphysCollection
838 // dimuon vertex has already been decorated with muon indices
839 auto dimuonIndex = dimuonTriggerObject->index();
840 const xAOD::Vertex* dimuon = state.dimuons.get(dimuonIndex);
841 if ( !dimuon || dimuonIndex >= state.trigBphysMuonIndices.size() ) {
842 ATH_MSG_ERROR( "Failed to find original muons the dimuon vertex had been built from" );
843 return StatusCode::FAILURE;
844 }
845
846 // create a new output Decision object, backed by the 'decisions' container.
848
849 std::vector<const DecisionIDContainer*> previousDecisionIDs;
850 for (const size_t& i : state.trigBphysMuonIndices.at(dimuonIndex)) {
851 const auto& muon = state.muons.at(i);
852 // attach all previous decisions: if the same previous decision is called twice, that's fine - internally takes care of that
853 // we already have an array of links to the previous decisions, so there is no need to use TrigCompositeUtils::linkToPrevious()
854 decision->addObjectCollectionLinks(TrigCompositeUtils::seedString(), muon.decisionLinks);
855 previousDecisionIDs.push_back(&muon.decisionIDs);
856 }
857
858 // set mandatory link to the trigger object
860
861 for (const auto& tool : hypoTools()) {
862 ATH_MSG_DEBUG( "Go to " << tool );
863 ATH_CHECK( tool->decideOnSingleObject(decision, previousDecisionIDs) );
864 }
865 }
866
867 return StatusCode::SUCCESS;
868}
869
870
871std::unique_ptr<xAOD::Vertex> TrigBmumuxComboHypo::fit(
872 const EventContext& context,
873 const std::vector<ElementLink<xAOD::TrackParticleContainer>>& trackParticleLinks,
874 Decay decay,
875 const xAOD::Vertex* dimuon) const {
876
877 ATH_MSG_VERBOSE( "Perform vertex fit" );
878
879 if (trackParticleLinks.size() < 2) {
880 ATH_MSG_WARNING( "At least two tracks should be given to the vertex fitter" );
881 return nullptr;
882 }
883
884 std::vector<const xAOD::TrackParticle*> tracklist(trackParticleLinks.size(), nullptr);
885 std::transform(trackParticleLinks.begin(), trackParticleLinks.end(), tracklist.begin(),
886 [](const ElementLink<xAOD::TrackParticleContainer>& link){ return *link; });
887
888 Amg::Vector3D startingPoint = Amg::Vector3D::Zero(3);
889 if (dimuon) {
890 startingPoint = Amg::Vector3D(dimuon->x(), dimuon->y(), dimuon->z());
891 }
892 else {
893 if (decay != Decay::kPsi_2mu && decay != Decay::kB_PsiPi && decay != Decay::kFastFit_2trk) {
894 ATH_MSG_WARNING( "Already fitted dimuon vertex should be provided for B -> mu1 mu2 trk1 .. trkN decay as a starting point for fitter" );
895 }
896 int flag = 0;
897 int errorcode = 0;
898 const Trk::Perigee& perigee1 = tracklist[0]->perigeeParameters();
899 const Trk::Perigee& perigee2 = tracklist[1]->perigeeParameters();
900 startingPoint = m_vertexPointEstimator->getCirclesIntersectionPoint(&perigee1, &perigee2, flag, errorcode);
901 if (errorcode != 0) startingPoint = Amg::Vector3D::Zero(3);
902 }
903 ATH_MSG_VERBOSE( "Starting point: (" << startingPoint(0) << ", " << startingPoint(1) << ", " << startingPoint(2) << ")" );
904
905 auto fitterState = m_vertexFitter->makeState(context);
906 m_vertexFitter->setMassInputParticles(s_trkMass[static_cast<size_t>(decay)], *fitterState);
907
908 // the combined momentum of D+/D_s+ candidate is constrained to point to the dimuon vertex
909 if (decay == Decay::kDs || decay == Decay::kDplus || decay == Decay::kD0) {
910 m_vertexFitter->setVertexForConstraint(*dimuon, *fitterState);
911 m_vertexFitter->setCnstType(8, *fitterState);
912 }
913
914 std::unique_ptr<xAOD::Vertex> vertex(m_vertexFitter->fit(tracklist, startingPoint, *fitterState));
915 if (!vertex) {
916 ATH_MSG_VERBOSE( "Vertex fit fails" );
917 return vertex;
918 }
919 if (vertex->chiSquared() > 150. || (decay == Decay::kPsi_2mu && vertex->chiSquared() > m_dimuon_chi2)) {
920 ATH_MSG_VERBOSE( "Fit is successful, but vertex chi2 is too high, we are not going to save it (chi2 = " << vertex->chiSquared() << ")" );
921 vertex.reset();
922 return vertex;
923 }
924 ATH_MSG_VERBOSE( "Fit is successful" );
925
926 // update trackParticleLinks()
927 vertex->clearTracks();
928 vertex->setTrackParticleLinks(trackParticleLinks);
929
930 return vertex;
931}
932
933
935 TrigBmumuxState& state,
936 const xAOD::Vertex& vertex,
938 const std::vector<double>& trkMass,
939 const ElementLink<xAOD::TrigBphysContainer>& dimuonLink) const {
940
941 const xAOD::TrigBphys* dimuon = (type != xAOD::TrigBphys::MULTIMU ? *dimuonLink : nullptr);
943
944 // refitted track momentum as a 4-vector for mass hypothesis defined by the given decay value
946 std::vector<xAOD::TrackParticle::GenVecFourMom_t> momenta;
947 if (!vertex.vxTrackAtVertexAvailable()) return nullptr;
948 for (size_t i = 0; i < vertex.vxTrackAtVertex().size(); ++i) {
949 const Trk::TrackParameters* perigee = vertex.vxTrackAtVertex()[i].perigeeAtVertex();
950 if (!perigee) return nullptr;
951 const Amg::Vector3D& p = perigee->momentum();
952 momenta.emplace_back(p.x(), p.y(), p.z(), trkMass[i]);
953 momentum += momenta.back();
954 }
955 if (isCascadeDecay) {
956 momentum += ROOT::Math::PtEtaPhiMVector(dimuon->pt(), dimuon->eta(), dimuon->phi(), dimuon->mass());
957 }
958
959 auto result = new xAOD::TrigBphys();
960 result->makePrivateStore();
961
962 float mass = momentum.M();
964 mass += PDG::mD0 - (momenta[0] + momenta[1]).M();
965 }
966 else if (type != xAOD::TrigBphys::MULTIMU) {
967 mass += PDG::mJpsi - (isCascadeDecay ? dimuon->mass() : (momenta[0] + momenta[1]).M());
968 }
969
970 result->initialise(0, momentum.Eta(), momentum.Phi(), momentum.Pt(), type, mass, xAOD::TrigBphys::EF);
971
972 result->setFitmass(momentum.M());
973 result->setFitx(vertex.x());
974 result->setFity(vertex.y());
975 result->setFitz(vertex.z());
976 result->setFitchi2(vertex.chiSquared());
977 result->setFitndof(vertex.numberDoF());
978
979 Amg::Vector3D productionVertex = (isCascadeDecay ? Amg::Vector3D(dimuon->fitx(), dimuon->fity(), dimuon->fitz()) : state.beamSpotPosition());
980 result->setLxy(Lxy(productionVertex, vertex));
981
982 // set all the particles associated with the decay
983 result->setTrackParticleLinks(vertex.trackParticleLinks());
984
985 // use lowerChainLink() as a link to the preceding dimuon trigger object
987 result->setLowerChainLink(dimuonLink);
988 }
989
991 "TrigBphys object:\n\t " <<
992 "roiId: " << result->roiId() << "\n\t " <<
993 "particleType: " << result->particleType() << "\n\t " <<
994 "level: " << result->level() << "\n\t " <<
995 "eta: " << result->eta() << "\n\t " <<
996 "phi: " << result->phi() << "\n\t " <<
997 "mass: " << result->mass() << "\n\t " <<
998 "fitmass: " << result->fitmass() << "\n\t " <<
999 "chi2/NDF: " << result->fitchi2() << " / " << result->fitndof() << "\n\t " <<
1000 "vertex: (" << result->fitx() << ", " << result->fity() << ", " << result->fitz() << ")" << "\n\t " <<
1001 "Lxy: " << result->lxy() );
1002
1003 return result;
1004}
1005
1006
1008
1009 if (lhs->charge() * rhs->charge() < 0.) return false;
1010 return (ROOT::Math::VectorUtil::DeltaR(lhs->genvecP4(), rhs->genvecP4()) < m_deltaR);
1011}
1012
1013
1014bool TrigBmumuxComboHypo::isIdenticalTracks(const xAOD::Muon* lhs, const xAOD::Muon* rhs) const {
1015
1016 return isIdenticalTracks(lhs->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle),
1017 rhs->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle));
1018}
1019
1020
1022
1023 auto p_mu = muon->genvecP4();
1024 auto p_trk = track->genvecP4();
1025 return (ROOT::Math::VectorUtil::DeltaPhi(p_mu, p_trk) < m_roiPhiWidth && std::abs(p_mu.eta() - p_trk.eta()) < m_roiEtaWidth);
1026}
1027
1028
1029double TrigBmumuxComboHypo::Lxy(const Amg::Vector3D& productionVertex, const xAOD::Vertex& decayVertex) const {
1030
1031 XYVector R(decayVertex.x() - productionVertex.x(), decayVertex.y() - productionVertex.y());
1032 XYVector pT;
1033
1034 if (!decayVertex.vxTrackAtVertexAvailable()) return -100.;
1035 const auto& tracks = decayVertex.vxTrackAtVertex();
1036 for (const auto& track : tracks) {
1037 const Trk::TrackParameters* perigee = track.perigeeAtVertex();
1038 if (!perigee) return -100.;
1039 const Amg::Vector3D& momentum = perigee->momentum();
1040 pT += XYVector(momentum.x(), momentum.y());
1041 }
1042 return R.Dot(pT.unit());
1043}
1044
1045
1046xAOD::TrackParticle::GenVecFourMom_t TrigBmumuxComboHypo::momentum(const xAOD::Vertex& vertex, const std::vector<double>& trkMass) const {
1047
1049 for (size_t i = 0; i < vertex.vxTrackAtVertex().size(); ++i) {
1050 const Trk::TrackParameters* perigee = vertex.vxTrackAtVertex()[i].perigeeAtVertex();
1051 const Amg::Vector3D& p = perigee->momentum();
1052 momentum += xAOD::TrackParticle::GenVecFourMom_t(p.x(), p.y(), p.z(), trkMass[i]);
1053 }
1054 return momentum;
1055}
1056
1057
1058bool TrigBmumuxComboHypo::passDimuonTrigger(const std::vector<const DecisionIDContainer*>& previousDecisionIDs) const {
1059
1060 if (previousDecisionIDs.size() != 2) {
1061 ATH_MSG_WARNING( "TrigBmumuxComboHypo::passDimuonTrigger() expects exactly two containers with previous decision IDs" );
1062 return false;
1063 }
1064
1065 for (const auto& tool : hypoTools()) {
1066 const std::vector<HLT::Identifier>& legDecisionIDs = tool->legDecisionIds();
1067 if (legDecisionIDs.size() == 1 && tool->legMultiplicity().at(0) >= 2) {
1068 // trigger with symmetric legs like HLT_2mu4_bBmumux_BsmumuPhi_L12MU4
1069 const DecisionID id = legDecisionIDs[0].numeric();
1070 if (TrigCompositeUtils::passed(id, *previousDecisionIDs[0]) && TrigCompositeUtils::passed(id, *previousDecisionIDs[1])) return true;
1071 }
1072 else if (legDecisionIDs.size() == 2) {
1073 // trigger with asymmetric legs like HLT_mu6_mu4_bBmumux_BsmumuPhi_L1MU6_2MU4
1074 bool direct = true;
1075 bool inverse = true;
1076 for (size_t i = 0; i < 2; ++i) {
1077 if (direct && !TrigCompositeUtils::passed(legDecisionIDs[i].numeric(), *previousDecisionIDs[i])) direct = false;
1078 if (inverse && !TrigCompositeUtils::passed(legDecisionIDs[i].numeric(), *previousDecisionIDs[1-i])) inverse = false;
1079 }
1080 if (direct || inverse) return true;
1081 }
1082 else {
1083 ATH_MSG_WARNING( "TrigBmumuxComboHypo can not check decisions for " << tool->name() );
1084 }
1085 }
1086 return false;
1087}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
int charge3(const T &p)
Definition AtlasPID.h:995
Base class for elements of a container that can have aux data.
ElementLink< xAOD::TrackParticleContainer > linkTrack(const xAOD::TrackParticle *trk)
size_t size() const
Number of registered mappings.
#define x
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)
bool msgLvl(const MSG::Level lvl) const
ComboHypo(const std::string &name, ISvcLocator *pSvcLocator)
Definition ComboHypo.cxx:13
const SG::WriteHandleKeyArray< TrigCompositeUtils::DecisionContainer > & decisionsOutput() const
Definition ComboHypo.h:42
const Combo::MultiplicityReqMap & triggerMultiplicityMap() const
Definition ComboHypo.h:43
ToolHandleArray< ComboHypoToolBase > & hypoTools()
Definition ComboHypo.h:45
const SG::ReadHandleKeyArray< TrigCompositeUtils::DecisionContainer > & decisionsInput() const
Definition ComboHypo.h:41
virtual StatusCode initialize() override
Definition ComboHypo.cxx:22
const T * get(size_type n) const
Access an element, as an rvalue.
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.
TrigCompositeUtils::DecisionContainer & decisions()
const TrigCompositeUtils::DecisionContainer & previousDecisions() const
xAOD::TrigBphysContainer & trigBphysCollection()
Amg::Vector3D beamSpotPosition() const
const EventContext & context() const
Group of local monitoring quantities and retain correlation when filling histograms
Declare a monitored scalar variable.
pointer_type ptr()
Dereference the pointer.
Gaudi::Property< bool > m_dimuon_rejectSameChargeTracks
Gaudi::Property< bool > m_BcToDplusMuMu_useFastFit
StatusCode findBmumuxCandidates_fit(TrigBmumuxState &, size_t dimuonIndex, bool makeCombinations=false) const
Perform fit of B decays for the topologies described above if makeCombinations = false.
Gaudi::Property< std::pair< double, double > > m_BcToMuMuPion_massRange
Gaudi::Property< float > m_BdToMuMuKstar0_chi2
ToolHandle< Trk::TrkVKalVrtFitter > m_vertexFitter
Gaudi::Property< bool > m_BdToMuMuKstar0_rejectSameChargeTracks
Gaudi::Property< bool > m_BdToMuMuKstar0_useFastFit
Gaudi::Property< double > m_BdToMuMuKstar0_minPionPt
Gaudi::Property< std::pair< double, double > > m_BcToDplusMuMu_massRange
std::unique_ptr< xAOD::Vertex > fit(const EventContext &context, const std::vector< ElementLink< xAOD::TrackParticleContainer > > &trackParticleLinks, Decay decay=kPsi_2mu, const xAOD::Vertex *dimuon=nullptr) const
Perform a vertex fit on selected tracks.
Gaudi::Property< std::pair< double, double > > m_BcToDsMuMu_phiMassRange
Gaudi::Property< std::pair< double, double > > m_BcToDstarMuMu_D0MassRange
Gaudi::Property< std::pair< double, double > > m_BcToDplusMuMu_DplusMassRange
Gaudi::Property< bool > m_BcToDsMuMu
Gaudi::Property< double > m_fastFit_2mu1trk_chi2
Gaudi::Property< bool > m_BcToMuMuPion
Gaudi::Property< double > m_BplusToMuMuKaon_minKaonPt
Gaudi::Property< std::pair< double, double > > m_BsToMuMuPhi1020_massRange
Gaudi::Property< double > m_BcToDsMuMu_minPionPt
Gaudi::Property< std::pair< double, double > > m_BdToMuMuKstar0_massRange
TrigCompositeUtils::DecisionIDContainer m_allowedIDs
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
Gaudi::Property< bool > m_LambdaBToMuMuProtonKaon
Gaudi::Property< double > m_roiPhiWidth
Gaudi::Property< double > m_BcToMuMuPion_minPionPt
SG::WriteHandleKey< xAOD::TrigBphysContainer > m_trigBphysContainerKey
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trackParticleContainerKey
Gaudi::Property< double > m_dimuon_chi2
Gaudi::Property< float > m_BsToMuMuPhi1020_chi2
Gaudi::Property< bool > m_BcToDsMuMu_useFastFit
Gaudi::Property< double > m_BsToMuMuPhi1020_minKaonPt
ToolHandle< InDet::VertexPointEstimator > m_vertexPointEstimator
bool passDimuonTrigger(const std::vector< const TrigCompositeUtils::DecisionIDContainer * > &previousDecisionIDs) const
Gaudi::Property< std::pair< double, double > > m_BsToMuMuPhi1020_phiMassRange
Gaudi::Property< float > m_BcToDsMuMu_chi2
Gaudi::Property< std::pair< double, double > > m_BplusToMuMuKaon_massRange
Gaudi::Property< std::pair< double, double > > m_BcToDsMuMu_DsMassRange
Gaudi::Property< float > m_BcToMuMuPion_chi2
Gaudi::Property< double > m_BcToDstarMuMu_minD0KaonPt
Gaudi::Property< double > m_BdToMuMuKstar0_minKaonPt
double Lxy(const Amg::Vector3D &productionVertex, const xAOD::Vertex &decayVertex) const
Calculate the Lxy (~distance between vertices) It is defined as the transverse distance between the p...
Gaudi::Property< bool > m_LambdaBToMuMuProtonKaon_useFastFit
ToolHandle< Reco::ITrackToVertex > m_trackToVertexTool
Gaudi::Property< std::pair< double, double > > m_LambdaBToMuMuProtonKaon_massRange
StatusCode findBmumuxCandidates(TrigBmumuxState &) const
Find B decays by appling next three subprocedures to each found dimuon candidate.
Gaudi::Property< std::pair< double, double > > m_dimuon_massRange
StatusCode mergeMuonsFromDecisions(TrigBmumuxState &) const
Go through state.previousDecisions(), fetch xAOD::Muons objects attached to decisions and save links ...
Gaudi::Property< float > m_BplusToMuMuKaon_chi2
Gaudi::Property< bool > m_BsToMuMuPhi1020
Gaudi::Property< double > m_LambdaBToMuMuProtonKaon_minKaonPt
Gaudi::Property< bool > m_BplusToMuMuKaon
Gaudi::Property< bool > m_BdToMuMuKstar0
Gaudi::Property< double > m_BcToDstarMuMu_minD0PionPt
StatusCode findBmumuxCandidates_selectTracks(TrigBmumuxState &, size_t dimuonIndex) const
Select tracks in vicinity of given dimuon vertex.
StatusCode findBmumuxCandidates_fastFit(TrigBmumuxState &, size_t dimuonIndex) const
Go through (dimuon+track) and (track+track) combinations found by findBmumuxCandidates_fit(makeCombin...
Gaudi::Property< std::pair< double, double > > m_BcToDstarMuMu_DstarMassRange
xAOD::TrigBphys * makeTriggerObject(TrigBmumuxState &state, const xAOD::Vertex &vertex, xAOD::TrigBphys::pType type=xAOD::TrigBphys::MULTIMU, const std::vector< double > &trkMass={PDG::mMuon, PDG::mMuon}, const ElementLink< xAOD::TrigBphysContainer > &dimuonLink=ElementLink< xAOD::TrigBphysContainer >()) const
Construct the trigger object that may be stored for debugging or matching.
Gaudi::Property< double > m_BcToDstarMuMu_minDstarPionPt
Gaudi::Property< std::pair< double, double > > m_LambdaBToMuMuProtonKaon_dimuonMassRange
Gaudi::Property< bool > m_BplusToMuMuKaon_useFastFit
Gaudi::Property< bool > m_BsToMuMuPhi1020_useFastFit
bool isInSameRoI(const xAOD::Muon *, const xAOD::TrackParticle *) const
Attempts to identify if the track is in the same RoI as the muon by comparing the angle with the RoI ...
Gaudi::Property< std::pair< double, double > > m_BcToDplusMuMu_dimuonMassRange
TrigBmumuxComboHypo()=delete
virtual StatusCode initialize() override
Gaudi::Property< double > m_roiEtaWidth
Gaudi::Property< std::pair< double, double > > m_BcToDsMuMu_massRange
Gaudi::Property< double > m_BcToDplusMuMu_minPionPt
Gaudi::Property< bool > m_BcToDstarMuMu
Gaudi::Property< std::pair< double, double > > m_BdToMuMuKstar0_KstarMassRange
static const std::vector< std::vector< double > > s_trkMass
Gaudi::Property< bool > m_BcToDstarMuMu_makeDstar
Gaudi::Property< std::pair< double, double > > m_BcToDstarMuMu_massRange
SG::ReadHandleKey< xAOD::MuonContainer > m_muonContainerKey
ToolHandle< GenericMonitoringTool > m_monTool
Gaudi::Property< double > m_BcToDstarMuMu_maxDstarPionZ0
Gaudi::Property< bool > m_BcToDplusMuMu
Gaudi::Property< double > m_LambdaBToMuMuProtonKaon_minProtonPt
Gaudi::Property< bool > m_BcToMuMuPion_useFastFit
Gaudi::Property< size_t > m_fitAttemptsBreakThreshold
Gaudi::Property< std::pair< double, double > > m_BcToDsMuMu_dimuonMassRange
Gaudi::Property< float > m_BcToDstarMuMu_chi2
Gaudi::Property< double > m_deltaR
Gaudi::Property< size_t > m_fitAttemptsWarningThreshold
Gaudi::Property< float > m_LambdaBToMuMuProtonKaon_chi2
Gaudi::Property< std::pair< double, double > > m_BcToMuMuPion_dimuonMassRange
Gaudi::Property< double > m_fastFit_2trk_chi2
StatusCode mergeTracksFromViews(TrigBmumuxState &) const
Go through state.previousDecisions() and fetch xAOD::TrackParticle objects associated with the neares...
virtual StatusCode execute(const EventContext &context) const override
Gaudi::Property< bool > m_BsToMuMuPhi1020_rejectSameChargeTracks
StatusCode createDecisionObjects(TrigBmumuxState &) const
Create a decision for each xAOD::TrigBphys object from state.trigBphysCollection() and save it to sta...
Gaudi::Property< float > m_BcToDplusMuMu_chi2
StatusCode findDimuonCandidates(TrigBmumuxState &) const
Make all possible dimuon combinations from state.muons(), fit muon InDet tracks to the common vertex,...
bool isInMassRange(double mass, const std::pair< double, double > &range) const
Checks that the given mass value falls into the specified range.
bool isIdenticalTracks(const xAOD::TrackParticle *lhs, const xAOD::TrackParticle *rhs) const
Attempts to identify identical tracks by selection on DeltaR.
Gaudi::Property< double > m_trkZ0
Gaudi::Property< double > m_BcToDplusMuMu_minKaonPt
Gaudi::Property< double > m_LambdaBToMuMuProtonKaon_minKstarMass
Gaudi::Property< std::pair< double, double > > m_BcToDstarMuMu_dimuonMassRange
Gaudi::Property< double > m_BcToDsMuMu_minKaonPt
xAOD::TrackParticle::GenVecFourMom_t momentum(const xAOD::Vertex &vertex, const std::vector< double > &trkMass) const
Calculate 4-momentum of the fitted vertex particle assuming the given masses.
State class for TrigBmumuxComboHypo algorithm.
std::vector< Muon > muons
std::map< const xAOD::TrackParticle *, double > selectedTrackZ0
void addTrackCombination(size_t i1)
bool isBadCombination(size_t i1) const
std::vector< ElementLink< xAOD::TrackParticleContainer > > selectedTracks
std::map< size_t, size_t > trackCombinations
StatusCode addTriggerObject(xAOD::TrigBphys *triggerObject)
std::vector< size_t > badTrackCombinations
xAOD::VertexContainer dimuons
std::vector< ElementLink< xAOD::TrackParticleContainer > > tracks
std::vector< std::array< size_t, 2 > > trigBphysMuonIndices
const Amg::Vector3D & momentum() const
Access method for the momentum.
const TrackParticle * trackParticle(TrackParticleType type) const
Returns a pointer (which can be NULL) to the TrackParticle used in identification of this muon.
Definition Muon_v1.cxx:482
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzM4D< double > > GenVecFourMom_t
Base 4 Momentum type for TrackParticle.
GenVecFourMom_t genvecP4() const
The full 4-momentum of the particle : GenVector form.
float charge() const
Returns the charge.
float pt() const
accessor method: pt
float eta() const
accessor method: eta
float phi() const
accessor method: phi
const TrigBphys_v1 * lowerChain() const
accessor method: lowerChain decay particle
float fitx() const
accessor method: x position of vertex
float fitz() const
accessor method: z position of vertex
float mass() const
accessor method: mass
pType
enum for different particle types
float fity() const
accessor method: y position of vertex
bool setObjectLink(const std::string &name, const ElementLink< CONTAINER > &link)
Set the link to an object.
bool addObjectCollectionLinks(const std::string &collectionName, const std::vector< ElementLink< CONTAINER > > &links)
Add links to multiple objects within a collection. Performs de-duplication.
float z() const
Returns the z position.
const TrackParticleLinks_t & trackParticleLinks() const
Get all the particles associated with the vertex.
const TrackParticle * trackParticle(size_t i) const
Get the pointer to a given track that was used in vertex reco.
float y() const
Returns the y position.
bool vxTrackAtVertexAvailable() const
Check if VxTrackAtVertices are attached to the object.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
const Amg::Vector3D & position() const
Returns the 3-pos.
float x() const
Returns the x position.
Eigen::Matrix< double, 3, 1 > Vector3D
ValuesCollection< T > Collection(std::string name, const T &collection)
Declare a monitored (double-convertible) collection.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
HLT::Identifier createLegName(const HLT::Identifier &chainIdentifier, size_t counter)
Generate the HLT::Identifier which corresponds to a specific leg of a given chain.
unsigned int DecisionID
const std::string & viewString()
const std::string & roiString()
Decision * newDecisionIn(DecisionContainer *dc, const std::string &name)
Helper method to create a Decision object, place it in the container and return a pointer to it.
const std::string & featureString()
bool passed(DecisionID id, const DecisionIDContainer &idSet)
checks if required decision ID is in the set of IDs in the container
const std::string & comboHypoAlgNodeName()
std::set< DecisionID > DecisionIDContainer
SG::WriteHandle< DecisionContainer > createAndStore(const SG::WriteHandleKey< DecisionContainer > &key, const EventContext &ctx)
Creates and right away records the DecisionContainer with the key.
LinkInfo< T > findLink(const Decision *start, const std::string &linkName, const bool suppressMultipleLinksWarning=false)
Perform a recursive search for ElementLinks of type T and name 'linkName', starting from Decision obj...
const std::string & seedString()
void decisionIDs(const Decision *d, DecisionIDContainer &destination)
Extracts DecisionIDs stored in the Decision object.
ElementLink< DecisionContainer > decisionToElementLink(const Decision *d, const EventContext &ctx)
Takes a raw pointer to a Decision and returns an ElementLink to the Decision.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ z0
Definition ParamDefs.h:64
ParametersBase< TrackParametersDim, Charged > TrackParameters
ElementLink< T > makeLink(const SG::View *view, const SG::ReadHandle< T > &handle, size_t index)
Create EL to a collection in view.
Definition ViewHelper.h:309
auto makeHandle(const SG::View *view, const KEY &key, const EventContext &ctx)
Create a view handle from a handle key.
Definition ViewHelper.h:273
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
TrigBphys_v1 TrigBphys
Definition TrigBphys.h:18
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.
Muon_v1 Muon
Reference the current persistent version:
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".
TrigBphysContainer_v1 TrigBphysContainer