ATLAS Offline Software
MLLSelector.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3  */
4 
6 #include "TopEvent/EventTools.h"
7 
8 namespace top {
9  MLLSelector::MLLSelector(const std::string& params) :
10  SignValueSelector("MLL", params) {
11  }
12 
13  bool MLLSelector::apply(const top::Event& event) const {
14  const xAOD::IParticle* lepton0 = 0;
15  const xAOD::IParticle* lepton1 = 0;
16 
17  if (event.m_electrons.size() == 2) {
18  lepton0 = event.m_electrons[0];
19  lepton1 = event.m_electrons[1];
20  } else if (event.m_muons.size() == 2) {
21  lepton0 = event.m_muons[0];
22  lepton1 = event.m_muons[1];
23  } else if (event.m_electrons.size() == 1 && event.m_muons.size() == 1) {
24  lepton0 = event.m_electrons[0];
25  lepton1 = event.m_muons[0];
26  } else {
27  throw std::runtime_error("MLLSelector::apply: Need two charge leptons");
28  }
29 
30  const double mll = top::invariantMass(*lepton0, *lepton1);
31  return checkFloat(mll, value());
32  }
33 
35  // If any of the required collections is a nullptr (i.e. has not been
36  // loaded) return false.
37  if (not event.m_electrons
38  or not event.m_muons) {
39  return false;
40  }
41 
42  const xAOD::IParticle* lepton0 = 0;
43  const xAOD::IParticle* lepton1 = 0;
44 
45  if (event.m_electrons->size() == 2) {
46  lepton0 = (*event.m_electrons)[0];
47  lepton1 = (*event.m_electrons)[1];
48  } else if (event.m_muons->size() == 2) {
49  lepton0 = (*event.m_muons)[0];
50  lepton1 = (*event.m_muons)[1];
51  } else if (event.m_electrons->size() == 1 && event.m_muons->size() == 1) {
52  lepton0 = (*event.m_electrons)[0];
53  lepton1 = (*event.m_muons)[0];
54  } else {
55  throw std::runtime_error("MLLSelector::applyParticleLevel: Need two charged leptons");
56  }
57 
58  const double mll = top::invariantMass(*lepton0, *lepton1);
59  return checkFloat(mll, value());
60  }
61 }
top
TopConfig A simple configuration that is NOT a singleton.
Definition: AnalysisTrackingHelper.cxx:58
top::invariantMass
double invariantMass(const xAOD::IParticle &p1, const xAOD::IParticle &p2)
Calculate the invariant mass of two particles.
Definition: EventTools.cxx:48
top::MLLSelector::apply
bool apply(const top::Event &event) const override
This does stuff based on the information in an event.
Definition: MLLSelector.cxx:13
xAOD::IParticle
Class providing the definition of the 4-vector interface.
Definition: Event/xAOD/xAODBase/xAODBase/IParticle.h:40
MLLSelector.h
top::MLLSelector::MLLSelector
MLLSelector(const std::string &params)
Definition: MLLSelector.cxx:9
top::SignValueSelector
Many of the tools need a sign (>=) and a value (2).
Definition: SignValueSelector.h:16
EventTools.h
A few functions for doing operations on particles / events. Currently holds code for dR,...
top::SignValueSelector::checkFloat
bool checkFloat(double value, double cut) const
Compare a cut supplied by the user with the value calculated in the event.
Definition: SignValueSelector.cxx:133
event
POOL::TEvent event(POOL::TEvent::kClassAccess)
top::ParticleLevelEvent
Definition: ParticleLevelEvent.h:24
top::MLLSelector::applyParticleLevel
bool applyParticleLevel(const top::ParticleLevelEvent &event) const override
This does stuff based on the information in a particle level event.
Definition: MLLSelector.cxx:34
top::SignValueSelector::value
double value() const
Get the cut value assigned in the constructor.
Definition: SignValueSelector.cxx:94
top::Event
Very simple class to hold event data after reading from a file.
Definition: Event.h:49
PowhegControl_ttFCNC_NLO.params
params
Definition: PowhegControl_ttFCNC_NLO.py:226