ATLAS Offline Software
MuDstarFilter.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 // -------------------------------------------------------------
5 // File: GeneratorFilters/MuDstarFilter.cxx
6 // Description:
7 //
8 // Allows the user to search for (mu D*) combinations
9 // with both same and opposite charges.
10 // For D*+-, the decay D*+ -> D0 pi_s+ is selected
11 // with D0 in the nominal decay mode, D0 -> K- pi+ (Br = 3.947%),
12 // and, if "D0Kpi_only=false", 14 main other decays to 2 or 3 particles (except nu_e and nu_mu)
13 // in case they can immitate the nominal decay:
14 // D0 -> K-mu+nu, K-e+nu, pi-mu+nu, pi-e+nu,
15 // K-mu+nu pi0, K-e+nu pi0, pi-mu+nu pi0, pi-e+nu pi0,
16 // pi-pi+, K-K+, K-pi+pi0, K-K+pi0, pi-pi+pi0, K-pi+gamma
17 // Doubly Cabbibo supressed modes are also considered
18 // Requirements for non-nominal decay modes:
19 // D*+ -> ("K-" "pi+") pi_s+ (+c.c.) charges
20 // mKpiMin < m("K" "pi") < mKpiMax
21 // m("K" "pi" pi_s) - m("K" "pi") < delta_m_Max
22 //
23 // AuthorList:
24 // L K Gladilin (gladilin@mail.cern.ch) March 2023
25 // Versions:
26 // 1.0.0; 2023.04.01 (no jokes) - baseline version
27 // 1.0.1; 2023.04.14 - allow photons for "m_D0Kpi_only" for a case of Photos
28 // 1.0.2; 2023.04.17 - skip history Photos lines (status=3)
29 //
30 // Header for this module:-
31 
33 
34 // Framework Related Headers:-
35 #include "GaudiKernel/MsgStream.h"
36 
37 
38 // Other classes used by this class:-
39 #include <math.h>
40 #include <limits>
41 
42 #include "CLHEP/Vector/LorentzVector.h"
43 
44 #include "TLorentzVector.h"
45 
47 
48 //--------------------------------------------------------------------------
49 MuDstarFilter::MuDstarFilter(const std::string & name,
50  ISvcLocator * pSvcLocator): GenFilter(name, pSvcLocator) {
51  //--------------------------------------------------------------------------
52  declareProperty("PtMinMuon", m_PtMinMuon = 2500.);
53  declareProperty("PtMaxMuon", m_PtMaxMuon = 1e9);
54  declareProperty("EtaRangeMuon", m_EtaRangeMuon = 2.7);
55  //
56  declareProperty("PtMinDstar", m_PtMinDstar = 4500.);
57  declareProperty("PtMaxDstar", m_PtMaxDstar = 1e9);
58  declareProperty("EtaRangeDstar", m_EtaRangeDstar = 2.7);
59  declareProperty("RxyMinDstar", m_RxyMinDstar = -1e9);
60  //
61  declareProperty("PtMinPis", m_PtMinPis = 450.);
62  declareProperty("PtMaxPis", m_PtMaxPis = 1e9);
63  declareProperty("EtaRangePis", m_EtaRangePis = 2.7);
64  //
65  declareProperty("PtMinKpi", m_PtMinKpi = 900.);
66  declareProperty("PtMaxKpi", m_PtMaxKpi = 1e9);
67  declareProperty("EtaRangeKpi", m_EtaRangeKpi = 2.7);
68  //
69  declareProperty("D0Kpi_only", m_D0Kpi_only = false);
70  //
71  declareProperty("mKpiMin", m_mKpiMin = 1665.);
72  declareProperty("mKpiMax", m_mKpiMax = 2065.);
73  //
74  declareProperty("delta_m_Max", m_delta_m_Max = 220.);
75  //
76  declareProperty("DstarMu_m_Max", m_DstarMu_m_Max = 12000.);
77 }
78 
79 //--------------------------------------------------------------------------
81  //--------------------------------------------------------------------------
82 
83 }
84 
85 //---------------------------------------------------------------------------
87  //---------------------------------------------------------------------------
88  ATH_MSG_INFO( "MuDstarFilter v1.02 INITIALISING " );
89 
90  ATH_MSG_INFO( "PtMinMuon = " << m_PtMinMuon );
91  ATH_MSG_INFO( "PtMaxMuon = " << m_PtMaxMuon );
92  ATH_MSG_INFO( "EtaRangeMuon = " << m_EtaRangeMuon );
93 
94  ATH_MSG_INFO( "PtMinDstar = " << m_PtMinDstar );
95  ATH_MSG_INFO( "PtMaxDstar = " << m_PtMaxDstar );
96  ATH_MSG_INFO( "EtaRangeDstar = " << m_EtaRangeDstar );
97 
98  ATH_MSG_INFO( "RxyMinDstar = " << m_RxyMinDstar );
99 
100  ATH_MSG_INFO( "PtMinPis = " << m_PtMinPis );
101  ATH_MSG_INFO( "PtMaxPis = " << m_PtMaxPis );
102  ATH_MSG_INFO( "EtaRangePis = " << m_EtaRangePis );
103 
104  ATH_MSG_INFO( "D0Kpi_only = " << m_D0Kpi_only );
105  ATH_MSG_INFO( "PtMinKpi = " << m_PtMinKpi );
106  ATH_MSG_INFO( "PtMaxKpi = " << m_PtMaxKpi );
107  ATH_MSG_INFO( "EtaRangeKpi = " << m_EtaRangeKpi );
108 
109  ATH_MSG_INFO( "mKpiMin = " << m_mKpiMin );
110  ATH_MSG_INFO( "mKpiMax = " << m_mKpiMax );
111 
112  ATH_MSG_INFO( "delta_m_Max = " << m_delta_m_Max );
113 
114  ATH_MSG_INFO( "DstarMu_m_Max = " << m_DstarMu_m_Max );
115 
116  return StatusCode::SUCCESS;
117 }
118 
119 //---------------------------------------------------------------------------
121  //---------------------------------------------------------------------------
122  return StatusCode::SUCCESS;
123 }
124 
125 //---------------------------------------------------------------------------
127  //---------------------------------------------------------------------------
128 
129  // Loop over all events in McEventCollection
131 
132  ATH_MSG_DEBUG(" MuDstarFilter filtering ");
133 
134  for (itr = events_const() -> begin(); itr != events_const() -> end(); ++itr) {
135 
136  double primx = 0.;
137  double primy = 0.;
138 
139  #ifdef HEPMC3
140  HepMC::ConstGenVertexPtr vprim = * (( * itr) -> vertices().begin());
141  #else
142  HepMC::GenVertexPtr vprim = * (( * itr) -> vertices_begin());
143  #endif
144 
145  primx = vprim -> position().x();
146  primy = vprim -> position().y();
147 
148  ATH_MSG_DEBUG("MuDstarFilter: PV x, y = " << primx << " , " << primy);
149 
150  int NumMuons = 0;
151  int NumDstars = 0;
152 
153  std::vector < HepMC::ConstGenParticlePtr > Muons;
154 
155  // Loop over all particles in the event
156  const HepMC::GenEvent * genEvt = ( * itr);
157  for (const auto & pitr: * genEvt) {
158 
159  if (!MC::isPhysical(pitr)) continue; // photos history line
160 
161  // muons
162  if (MC::isMuon(pitr)) {
163  if ((pitr -> momentum().perp() >= m_PtMinMuon) &&
164  (pitr -> momentum().perp() < m_PtMaxMuon) &&
165  (std::abs(pitr -> momentum().eta()) < m_EtaRangeMuon)) {
166  ++NumMuons;
167  Muons.push_back(pitr);
168  }
169  }
170  }
171 
172  ATH_MSG_DEBUG( "MuDstarFilter: NumMuons = " << NumMuons);
173 
174  if (NumMuons == 0) break;
175  for (const auto & pitr: * genEvt) {
176 
177  if (!MC::isPhysical(pitr)) continue; // photos history line
178 
179  // Dstars
180  if (std::abs(pitr -> pdg_id()) == MC::DSTAR) {
181  if ((pitr -> momentum().perp() >= m_PtMinDstar) &&
182  (pitr -> momentum().perp() < m_PtMaxDstar) &&
183  (std::abs(pitr -> momentum().eta()) < m_EtaRangeDstar)) {
184 
185  //Check if has end_vertex
186  if (!(pitr -> end_vertex())) continue;
187 
188  double Rxy = std::sqrt(pow(pitr -> end_vertex() -> position().x() - primx, 2) + pow(pitr -> end_vertex() -> position().y() - primy, 2));
189 
190  ATH_MSG_DEBUG("MuDstarFilter: R_xy(Dstar) = " << Rxy);
191 
192  if (Rxy < m_RxyMinDstar) continue;
193 
194  // Child
195  #ifdef HEPMC3
196  auto firstChild = pitr -> end_vertex() -> particles_out().begin();
197  auto endChild = pitr -> end_vertex() -> particles_out().end();
198  auto thisChild = firstChild;
199  #else
200  HepMC::GenVertex::particle_iterator firstChild =
201  pitr -> end_vertex() -> particles_begin(HepMC::children);
202  HepMC::GenVertex::particle_iterator endChild =
203  pitr -> end_vertex() -> particles_end(HepMC::children);
204  HepMC::GenVertex::particle_iterator thisChild = firstChild;
205  #endif
206 
207  if (( * firstChild) -> pdg_id() == pitr -> pdg_id()) continue;
208 
209  TLorentzVector p4_K;
210  TLorentzVector p4_pi;
211  TLorentzVector p4_pis;
212  TLorentzVector p4_D0;
213  TLorentzVector p4_Dstar;
214  TLorentzVector p4_Mu;
215  TLorentzVector p4_DstarMu;
216 
217  int pis_pdg = 0;
218  int K_pdg = 0;
219 
220  int NumD0 = 0;
221  int NumPis = 0;
222 
223  for (; thisChild != endChild; ++thisChild) {
224 
225  if (!MC::isPhysical(*thisChild)) continue; // photos history line
226 
227  if (std::abs(( * thisChild) -> pdg_id()) == MC::PIPLUS) {
228  if ((( * thisChild) -> momentum().perp() >= m_PtMinPis) &&
229  (( * thisChild) -> momentum().perp() < m_PtMaxPis) &&
230  (std::abs(( * thisChild) -> momentum().eta()) < m_EtaRangePis)) {
231  ++NumPis;
232 
233  p4_pis.SetPtEtaPhiM(( * thisChild) -> momentum().perp(), ( * thisChild) -> momentum().eta(), ( * thisChild) -> momentum().phi(), m_PionMass);
234  pis_pdg = ( * thisChild) -> pdg_id();
235  }
236  }
237  } // thisChild
238 
239  ATH_MSG_DEBUG("MuDstarFilter: NumPis = " << NumPis );
240 
241  if (NumPis != 1) continue;
242 
243  HepMC::ConstGenParticlePtr D0Child1 = 0;
244  HepMC::ConstGenParticlePtr D0Child2 = 0;
245  HepMC::ConstGenParticlePtr D0ChildMu = 0;
246 
247  int NumChildD0 = 0;
248  int NumChildD0Charged = 0;
249  int NumChildD0neutrinos = 0;
250  int NumChildD0gammas = 0;
251  int ChargeD0Child1 = 0;
252  int ChargeD0Child2 = 0;
253  int NumChildD0K = 0;
254  int NumChildD0pi = 0;
255  int NumChildD0mu = 0;
256 
257  for (thisChild = firstChild; thisChild != endChild; ++thisChild) {
258 
259  if (!MC::isPhysical(*thisChild)) continue; // photos history line
260 
261  if (std::abs(( * thisChild) -> pdg_id()) == MC::D0) {
262  if ((( * thisChild) -> end_vertex())) {
263  #ifdef HEPMC3
264  auto firstChild1 = ( * thisChild) -> end_vertex() -> particles_out().begin();
265  auto endChild1 = ( * thisChild) -> end_vertex() -> particles_out().end();
266  auto thisChild1 = firstChild1;
267  #else
268  HepMC::GenVertex::particle_iterator firstChild1 =
269  ( * thisChild) -> end_vertex() -> particles_begin(HepMC::children);
270  HepMC::GenVertex::particle_iterator endChild1 =
271  ( * thisChild) -> end_vertex() -> particles_end(HepMC::children);
272  HepMC::GenVertex::particle_iterator thisChild1 = firstChild1;
273 
274  #endif
275 
276  for (; thisChild1 != endChild1; ++thisChild1) {
277 
278  if (!MC::isPhysical(*thisChild1)) continue; // photos history line
279 
280  if (std::abs(( * thisChild1) -> pdg_id()) == MC::ELECTRON || std::abs(( * thisChild1) -> pdg_id()) == MC::MUON ||
281  std::abs(( * thisChild1) -> pdg_id()) == MC::PIPLUS || std::abs(( * thisChild1) -> pdg_id()) == MC::KPLUS) {
282 
283  ++NumChildD0;
284 
285  if ((( * thisChild1) -> momentum().perp() >= m_PtMinKpi) &&
286  (( * thisChild1) -> momentum().perp() < m_PtMaxKpi) &&
287  (std::abs(( * thisChild1) -> momentum().eta()) < m_EtaRangeKpi)) {
288  ++NumChildD0Charged;
289 
290  if (NumChildD0Charged == 1) {
291  D0Child1 = ( * thisChild1);
292  if (( * thisChild1) -> pdg_id() == -MC::ELECTRON || ( * thisChild1) -> pdg_id() == -MC::MUON ||
293  ( * thisChild1) -> pdg_id() == MC::PIPLUS || ( * thisChild1) -> pdg_id() == MC::KPLUS) {
294  ChargeD0Child1 = +1;
295  } else {
296  ChargeD0Child1 = -1;
297  }
298  }
299  if (NumChildD0Charged == 2) {
300  D0Child2 = ( * thisChild1);
301  if (( * thisChild1) -> pdg_id() == -MC::ELECTRON || ( * thisChild1) -> pdg_id() == -MC::MUON ||
302  ( * thisChild1) -> pdg_id() == MC::PIPLUS || ( * thisChild1) -> pdg_id() == MC::KPLUS) {
303  ChargeD0Child2 = +1;
304  } else {
305  ChargeD0Child2 = -1;
306  }
307  }
308  if (std::abs(( * thisChild1) -> pdg_id()) == MC::MUON) {
309  ++NumChildD0mu;
310  D0ChildMu = ( * thisChild1);
311  }
312  if (std::abs(( * thisChild1) -> pdg_id()) == MC::PIPLUS) ++NumChildD0pi;
313  if (std::abs(( * thisChild1) -> pdg_id()) == MC::KPLUS) {
314  ++NumChildD0K;
315  K_pdg = ( * thisChild1) -> pdg_id();
316  }
317  }
318 
319  } else if (std::abs(( * thisChild1) -> pdg_id()) == MC::PI0) {
320  ++NumChildD0;
321  } else if (std::abs(( * thisChild1) -> pdg_id()) == MC::NU_E || std::abs(( * thisChild1) -> pdg_id()) == MC::NU_MU) {
322  ++NumChildD0neutrinos;
323  } else if (std::abs(( * thisChild1) -> pdg_id()) == MC::PHOTON) {
324  ++NumChildD0gammas;
325  } else if (std::abs(( * thisChild1) -> pdg_id()) == MC::K0 || std::abs(( * thisChild1) -> pdg_id()) == MC::K0L ||
326  std::abs(( * thisChild1) -> pdg_id()) == MC::K0S) {
327  ++NumChildD0;
328  ++NumChildD0;
329  } else if ((( * thisChild1) -> end_vertex())) {
330 
331  //
332  #ifdef HEPMC3
333  auto firstChild2 = ( * thisChild1) -> end_vertex() -> particles_out().begin();
334  auto endChild2 = ( * thisChild1) -> end_vertex() -> particles_out().end();
335  auto thisChild2 = firstChild2;
336  #else
337  HepMC::GenVertex::particle_iterator firstChild2 =
338  ( * thisChild1) -> end_vertex() -> particles_begin(HepMC::children);
339  HepMC::GenVertex::particle_iterator endChild2 =
340  ( * thisChild1) -> end_vertex() -> particles_end(HepMC::children);
341  HepMC::GenVertex::particle_iterator thisChild2 = firstChild2;
342 
343  #endif
344  for (; thisChild2 != endChild2; ++thisChild2) {
345 
346  if (!MC::isPhysical(*thisChild2)) continue; // photos history line
347 
348  if (std::abs(( * thisChild2) -> pdg_id()) == MC::ELECTRON || std::abs(( * thisChild2) -> pdg_id()) == MC::MUON ||
349  std::abs(( * thisChild2) -> pdg_id()) == MC::PIPLUS || std::abs(( * thisChild2) -> pdg_id()) == MC::KPLUS) {
350 
351  ++NumChildD0;
352 
353  if ((( * thisChild2) -> momentum().perp() >= m_PtMinKpi) &&
354  (( * thisChild2) -> momentum().perp() < m_PtMaxKpi) &&
355  (std::abs(( * thisChild2) -> momentum().eta()) < m_EtaRangeKpi)) {
356  ++NumChildD0Charged;
357 
358  if (NumChildD0Charged == 1) {
359  D0Child1 = ( * thisChild2);
360  if (( * thisChild2) -> pdg_id() == -MC::ELECTRON || ( * thisChild2) -> pdg_id() == -MC::MUON ||
361  ( * thisChild2) -> pdg_id() == MC::PIPLUS || ( * thisChild2) -> pdg_id() == MC::KPLUS) {
362  ChargeD0Child1 = +1;
363  } else {
364  ChargeD0Child1 = -1;
365  }
366  }
367  if (NumChildD0Charged == 2) {
368  D0Child2 = ( * thisChild2);
369  if (( * thisChild2) -> pdg_id() == -MC::ELECTRON || ( * thisChild2) -> pdg_id() == -MC::MUON ||
370  ( * thisChild2) -> pdg_id() == MC::PIPLUS || ( * thisChild2) -> pdg_id() == MC::KPLUS) {
371  ChargeD0Child2 = +1;
372  } else {
373  ChargeD0Child2 = -1;
374  }
375  }
376  if (std::abs(( * thisChild2) -> pdg_id()) == MC::MUON) {
377  ++NumChildD0mu;
378  D0ChildMu = ( * thisChild2);
379  }
380  }
381 
382  } else if (std::abs(( * thisChild2) -> pdg_id()) == MC::PI0) {
383  ++NumChildD0;
384  } else if (std::abs(( * thisChild2) -> pdg_id()) == MC::NU_E || std::abs(( * thisChild2) -> pdg_id()) == MC::NU_MU) {
385  ++NumChildD0neutrinos;
386  } else if (std::abs(( * thisChild2) -> pdg_id()) == MC::PHOTON) {
387  ++NumChildD0gammas;
388  } else if (std::abs(( * thisChild2) -> pdg_id()) == MC::K0 || std::abs(( * thisChild2) -> pdg_id()) == MC::K0L ||
389  std::abs(( * thisChild2) -> pdg_id()) == MC::K0S) {
390  ++NumChildD0;
391  ++NumChildD0;
392  } else if ((( * thisChild2) -> end_vertex())) {
393  ++NumChildD0;
394  ++NumChildD0;
395  } else {
396  ++NumChildD0;
397  ++NumChildD0;
398  ATH_MSG_DEBUG("MuDstarFilter: unexpected D0 granddaughter = " << (*thisChild2)->pdg_id() );
399  }
400  } // thisChild2
401 
402  } else {
403  ++NumChildD0;
404  ++NumChildD0;
405  ATH_MSG_DEBUG("MuDstarFilter: unexpected D0 daughter = " << (*thisChild1)->pdg_id() );
406  }
407  } // thisChild1
408 
409  ATH_MSG_DEBUG("MuDstarFilter: NumChildD0, NumChildD0Charged = " << NumChildD0 << " , " << NumChildD0Charged);
410 
411  if (NumChildD0 <= 3 && NumChildD0Charged == 2 && ChargeD0Child1 * ChargeD0Child2 < 0) {
412  if (m_D0Kpi_only) {
413  if (NumChildD0 == 2 && NumChildD0K == 1 && NumChildD0pi == 1) ++NumD0;
414  } else {
415  ++NumD0;
416  }
417  }
418  } // enVertex
419  } // D0
420 
421  ATH_MSG_DEBUG("MuDstarFilter: NumD0 = " << NumD0 );
422 
423  if (NumD0 == 1) {
424 
425  if (pis_pdg * ChargeD0Child1 < 0) {
426 
427  p4_K.SetPtEtaPhiM(D0Child1 -> momentum().perp(), D0Child1 -> momentum().eta(), D0Child1 -> momentum().phi(), m_KaonMass);
428  p4_pi.SetPtEtaPhiM(D0Child2 -> momentum().perp(), D0Child2 -> momentum().eta(), D0Child2 -> momentum().phi(), m_PionMass);
429 
430  } else {
431 
432  p4_K.SetPtEtaPhiM(D0Child2 -> momentum().perp(), D0Child2 -> momentum().eta(), D0Child2 -> momentum().phi(), m_KaonMass);
433  p4_pi.SetPtEtaPhiM(D0Child1 -> momentum().perp(), D0Child1 -> momentum().eta(), D0Child1 -> momentum().phi(), m_PionMass);
434  }
435 
436  p4_D0 = p4_K + p4_pi;
437  double mKpi = p4_D0.M();
438 
439  ATH_MSG_DEBUG("MuDstarFilter: mKpi = " << mKpi );
440 
441  if (mKpi >= m_mKpiMin && mKpi <= m_mKpiMax) {
442 
443  p4_Dstar = p4_D0 + p4_pis;
444 
445  double delta_m = p4_Dstar.M() - mKpi;
446 
447  ATH_MSG_DEBUG("MuDstarFilter: delta_m = " << delta_m );
448 
449  if (delta_m <= m_delta_m_Max) {
450  ++NumDstars;
451 
452  ATH_MSG_DEBUG("MuDstarFilter: NumDstars = " << NumDstars );
453 
454  for (size_t i = 0; i < Muons.size(); ++i) {
455 
456  if (NumChildD0mu == 1) {
457  ATH_MSG_DEBUG("MuDstarFilter: Mu(pT), D0Mu(pT) = " << Muons[i]->momentum().perp() << " , " << D0ChildMu->momentum().perp());
458  if (std::fabs(Muons[i] -> momentum().perp() - D0ChildMu -> momentum().perp())< std::numeric_limits<double>::epsilon()) continue;
459  ATH_MSG_DEBUG("MuDstarFilter: Mu(pT), D0Mu(pT) = " << Muons[i]->momentum().perp() << " , " << D0ChildMu->momentum().perp() ) ;
460  }
461 
462  p4_Mu.SetPtEtaPhiM(Muons[i] -> momentum().perp(), Muons[i] -> momentum().eta(), Muons[i] -> momentum().phi(), m_MuonMass);
463 
464  p4_DstarMu = p4_Dstar + p4_Mu;
465 
466  ATH_MSG_DEBUG("MuDstarFilter: p4_DstarMu.M() = " << p4_DstarMu.M());
467 
468  if (p4_DstarMu.M() <= m_DstarMu_m_Max) {
469 
470  ATH_MSG_INFO("MuDstarFilter: MuDstar candidate found" );
471  ATH_MSG_INFO("MuDstarFilter: p4_DstarMu.M() = " << p4_DstarMu.M() );
472  ATH_MSG_INFO("MuDstarFilter: NumChildD0, NumChildD0Charged = " << NumChildD0 << " , " << NumChildD0Charged );
473  ATH_MSG_INFO("MuDstarFilter: NumChildD0K, NumChildD0pi, NumChildD0mu = " << NumChildD0K << " , " << NumChildD0pi << " , " << NumChildD0mu );
474 
475  if (NumChildD0mu == 1) {
476 
477  ATH_MSG_DEBUG("MuDstarFilter: Mu(pT), D0Mu(pT) = " << Muons[i]->momentum().perp() << " , " << D0ChildMu->momentum().perp() );
478  }
479 
480  ATH_MSG_INFO("MuDstarFilter: NumChildD0neutrinos, NumChildD0gammas = " << NumChildD0neutrinos << " , " << NumChildD0gammas );
481  ATH_MSG_INFO("MuDstarFilter: pis_pdg, K_pdg, ChargeD0Child1, ChargeD0Child2 = " << pis_pdg << " , " << K_pdg << " , " << ChargeD0Child1 << " , " << ChargeD0Child2 );
482 
483  setFilterPassed(true);
484  return StatusCode::SUCCESS;
485  }
486  } // for i
487 
488  } //delta_m
489  } // mKpi
490  } // NumD0
491  } // thisChild
492  } // PtMinDstar
493  } // DSTAR
494  } // pitr
495  } // itr
496 
497  //
498  // if we get here we have failed
499  //
500  setFilterPassed(false);
501  return StatusCode::SUCCESS;
502 }
503 
HepMC::GenVertexPtr
HepMC::GenVertex * GenVertexPtr
Definition: GenVertex.h:59
MuDstarFilter::filterInitialize
virtual StatusCode filterInitialize()
Definition: MuDstarFilter.cxx:86
MuDstarFilter::m_PtMinKpi
double m_PtMinKpi
Definition: MuDstarFilter.h:62
DataModel_detail::const_iterator
Const iterator class for DataVector/DataList.
Definition: DVLIterator.h:82
MuDstarFilter::m_PtMaxKpi
double m_PtMaxKpi
Definition: MuDstarFilter.h:63
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:64
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
perp
Scalar perp() const
perp method - perpenticular length
Definition: AmgMatrixBasePlugin.h:35
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:79
GenBase::events_const
const McEventCollection * events_const() const
Access the current event's McEventCollection (const)
Definition: GenBase.h:96
AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
conifer::pow
constexpr int pow(int x)
Definition: conifer.h:20
MuDstarFilter::m_PtMaxMuon
double m_PtMaxMuon
Definition: MuDstarFilter.h:50
PlotCalibFromCool.begin
begin
Definition: PlotCalibFromCool.py:94
MuDstarFilter::m_PionMass
const double m_PionMass
Definition: MuDstarFilter.h:78
MuDstarFilter::filterEvent
virtual StatusCode filterEvent()
Definition: MuDstarFilter.cxx:126
MuDstarFilter.h
x
#define x
MuDstarFilter::m_delta_m_Max
double m_delta_m_Max
Definition: MuDstarFilter.h:71
mergePhysValFiles.end
end
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:93
MuDstarFilter::MuDstarFilter
MuDstarFilter(const std::string &name, ISvcLocator *pSvcLocator)
Definition: MuDstarFilter.cxx:49
MC::isPhysical
bool isPhysical(const T &p)
Definition: HepMCHelpers.h:32
MuDstarFilter::m_PtMinMuon
double m_PtMinMuon
Definition: MuDstarFilter.h:49
GenFilter
Base class for event generator filtering modules.
Definition: GenFilter.h:30
MuDstarFilter::m_EtaRangePis
double m_EtaRangePis
Definition: MuDstarFilter.h:60
MuDstarFilter::m_DstarMu_m_Max
double m_DstarMu_m_Max
Definition: MuDstarFilter.h:73
MuDstarFilter::filterFinalize
virtual StatusCode filterFinalize()
Definition: MuDstarFilter.cxx:120
ParticleGun_EoverP_Config.momentum
momentum
Definition: ParticleGun_EoverP_Config.py:63
lumiFormat.i
int i
Definition: lumiFormat.py:92
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
MuDstarFilter::~MuDstarFilter
virtual ~MuDstarFilter()
Definition: MuDstarFilter.cxx:80
MUON
xAOD::Muon MUON
D3PD INCLUDES.
Definition: TileCellFillerTool.h:37
MuDstarFilter::m_D0Kpi_only
bool m_D0Kpi_only
Definition: MuDstarFilter.h:66
MuDstarFilter::m_KaonMass
const double m_KaonMass
Definition: MuDstarFilter.h:79
MuDstarFilter::m_MuonMass
const double m_MuonMass
Definition: MuDstarFilter.h:77
MuDstarFilter::m_EtaRangeDstar
double m_EtaRangeDstar
Definition: MuDstarFilter.h:55
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
MuDstarFilter::m_PtMaxPis
double m_PtMaxPis
Definition: MuDstarFilter.h:59
y
#define y
MuDstarFilter::m_mKpiMax
double m_mKpiMax
Definition: MuDstarFilter.h:69
python.DecayParser.children
children
Definition: DecayParser.py:32
MuDstarFilter::m_EtaRangeKpi
double m_EtaRangeKpi
Definition: MuDstarFilter.h:64
MuDstarFilter::m_mKpiMin
double m_mKpiMin
Definition: MuDstarFilter.h:68
HepMC::ConstGenVertexPtr
const HepMC::GenVertex * ConstGenVertexPtr
Definition: GenVertex.h:60
MuDstarFilter::m_PtMinDstar
double m_PtMinDstar
Definition: MuDstarFilter.h:53
MuDstarFilter::m_RxyMinDstar
double m_RxyMinDstar
Definition: MuDstarFilter.h:56
MuDstarFilter::m_EtaRangeMuon
double m_EtaRangeMuon
Definition: MuDstarFilter.h:51
MuDstarFilter::m_PtMinPis
double m_PtMinPis
Definition: MuDstarFilter.h:58
HepMCHelpers.h
MuDstarFilter::m_PtMaxDstar
double m_PtMaxDstar
Definition: MuDstarFilter.h:54
isMuon
bool isMuon(const T &p)
Definition: AtlasPID.h:145