ATLAS Offline Software
Loading...
Searching...
No Matches
TTbarWToLeptonFilter Class Reference

Require the event to contain at least one charged lepton (from W decay, which should come from top) with pt at or above Ptcut. More...

#include <TTbarWToLeptonFilter.h>

Inheritance diagram for TTbarWToLeptonFilter:
Collaboration diagram for TTbarWToLeptonFilter:

Public Member Functions

 TTbarWToLeptonFilter (const std::string &name, ISvcLocator *pSvcLocator)
virtual StatusCode filterEvent ()
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const
Event loop algorithm methods: not to be overloaded
StatusCode initialize ()
StatusCode execute ()
StatusCode finalize ()
Gen-specific event loop methods: to be overloaded!
virtual StatusCode filterInitialize ()
virtual StatusCode filterFinalize ()
Counters
int nPassed () const
int nFailed () const
int nNeeded () const
Event collection accessors (const and non-const)
HepMC::GenEvent *event ATLAS_NOT_CONST_THREAD_SAFE ()
 Access the current signal event (first in the McEventCollection)
McEventCollection *events ATLAS_NOT_CONST_THREAD_SAFE ()
 Access the current event's McEventCollection.
const HepMC::GenEvent * event_const () const
 Access the current signal event (const)
const McEventCollectionevents_const () const
 Access the current event's McEventCollection (const)
const McEventCollectionevents_const (const EventContext &ctx) const
Particle data accessors
const ServiceHandle< IPartPropSvc > partPropSvc () const
 Access the particle property service.
const HepPDT::ParticleDataTable & particleTable () const
 Get a particle data table.
const HepPDT::ParticleDataTable & pdt () const
 Shorter alias to get a particle data table.
const HepPDT::ParticleData * particleData (int pid) const
 Access an element in the particle data table.

Public Attributes

 Ptcut
 Add this filter to the algs required to be successful for streaming if "TTbarWToLeptonFilter" not in StreamEVGEN.RequireAlgs: StreamEVGEN.RequireAlgs += ["TTbarWToLeptonFilter"].

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
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)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Protected Attributes

Counters and requirements
int m_nPass
int m_nFail
int m_nNeeded
Properties
std::string m_mcEventKey {}
 StoreGate key for the MC event collection (defaults to GEN_EVENT)
BooleanProperty m_mkMcEvent {this, "MakeMcEvent", false, "Create a new MC event collection if it doesn't exist"}
 Flag to determine if a new MC event collection should be made if it doesn't exist.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

double m_Ptmin
int m_numLeptons
bool m_fourTopsFilter
bool m_SSMLFilter
DataObjIDColl m_extendedExtraObjects
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Utility event-mangling functions

Todo
Replace with HepMC units when available
ServiceHandle< IPartPropSvc > m_ppSvc {this, "PartPropSvc", "PartPropSvc"}
 Handle on the particle property service.
SG::ReadHandleKey< McEventCollectionm_mcevents_const { this, "McEventKey", "GEN_EVENT", "StoreGate key of the MC event collection" }
 Const handle to the MC event collection.
void GeVToMeV (HepMC::GenEvent *evt)
 Scale event energies/momenta by x 1000.
void MeVToGeV (HepMC::GenEvent *evt)
 Scale event energies/momenta by x 1/1000.
void cmTomm (HepMC::GenEvent *evt)
 Scale event lengths by x 10.
void mmTocm (HepMC::GenEvent *evt)
 Scale event lengths by x 1/10.

Detailed Description

Require the event to contain at least one charged lepton (from W decay, which should come from top) with pt at or above Ptcut.

Events that do not contain t AND t bar quarks are rejected. Only tops decaying to W X are analyzed and counted in this algorithm.

Author
Gia Khoriauli, June 2008
Andy Buckley, extension to specific lepton multiplicities, April 2012

Definition at line 16 of file TTbarWToLeptonFilter.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< Algorithm > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ TTbarWToLeptonFilter()

TTbarWToLeptonFilter::TTbarWToLeptonFilter ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 7 of file TTbarWToLeptonFilter.cxx.

8 : GenFilter(name,pSvcLocator)
9{
10 declareProperty("Ptcut", m_Ptmin=200000.);
11 declareProperty("NumLeptons", m_numLeptons=-1); // Negative for >0, positive integers for the specific number
12 declareProperty("fourTopsFilter", m_fourTopsFilter=false);//four top filter or not
13 declareProperty("SSMLFilter", m_SSMLFilter=false);//Same sign multilepton filter or not
14
15}
Base class for event generator filtering modules.
Definition GenFilter.h:30

Member Function Documentation

◆ ATLAS_NOT_CONST_THREAD_SAFE() [1/2]

HepMC::GenEvent *event GenBase::ATLAS_NOT_CONST_THREAD_SAFE ( )
inlineinherited

Access the current signal event (first in the McEventCollection)

Note
This function will make a new McEventCollection if there is not already a valid one and MakeMcEvent=True.

Definition at line 76 of file GenBase.h.

76 {
77 if (events()->empty())
78 ATH_MSG_ERROR("McEventCollection is empty during first event access");
79 return *(events()->begin());
80 }
#define ATH_MSG_ERROR(x)
static const Attributes_t empty

◆ ATLAS_NOT_CONST_THREAD_SAFE() [2/2]

McEventCollection *events GenBase::ATLAS_NOT_CONST_THREAD_SAFE ( )
inherited

Access the current event's McEventCollection.

Note
This function will make a new McEventCollection if there is not already a valid one and MakeMcEvent=True.

◆ cmTomm()

void GenBase::cmTomm ( HepMC::GenEvent * evt)
protectedinherited

Scale event lengths by x 10.

Definition at line 78 of file GenBase.cxx.

78 {
79 for (HepMC::GenEvent::vertex_iterator vtx = evt->vertices_begin(); vtx != evt->vertices_end(); ++vtx) {
80 const HepMC::FourVector fv((*vtx)->position().x() * 10,
81 (*vtx)->position().y() * 10,
82 (*vtx)->position().z() * 10,
83 (*vtx)->position().t() * 10);
84 (*vtx)->set_position(fv);
85 }
86}

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ event_const()

const HepMC::GenEvent * GenBase::event_const ( ) const
inlineinherited

Access the current signal event (const)

Definition at line 83 of file GenBase.h.

83 {
84 if (events_const()->empty())
85 ATH_MSG_ERROR("Const McEventCollection is empty during first event access");
86 return *(events_const()->begin());
87 }
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
const McEventCollection * events_const() const
Access the current event's McEventCollection (const)
Definition GenBase.h:96

◆ events_const() [1/2]

const McEventCollection * GenBase::events_const ( ) const
inlineinherited

Access the current event's McEventCollection (const)

Definition at line 96 of file GenBase.h.

96 {
97 return events_const( getContext() );
98 }

◆ events_const() [2/2]

const McEventCollection * GenBase::events_const ( const EventContext & ctx) const
inlineinherited

Definition at line 99 of file GenBase.h.

99 {
100 SG::ReadHandle<McEventCollection> ret = SG::makeHandle(m_mcevents_const, ctx);
101 if (!ret.isValid())
102 ATH_MSG_ERROR("No McEventCollection found in StoreGate with key " << m_mcevents_const.key());
103 return ret.cptr();
104 }
SG::ReadHandleKey< McEventCollection > m_mcevents_const
Const handle to the MC event collection.
Definition GenBase.h:163
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode GenFilter::execute ( )
virtualinherited
Todo
Probably the filter should only look at the first event... right?

Reimplemented from GenBase.

Definition at line 29 of file GenFilter.cxx.

29 {
30 if (events_const()->empty()) {
31 ATH_MSG_ERROR("No events found in McEventCollection");
32 return StatusCode::FAILURE;
33 } else if (events_const()->size() > 1) {
35 ATH_MSG_WARNING("More than one event in current McEventCollection -- which is valid?");
36 }
38#ifdef HEPMC3
39 if (filterPassed() || m_keepAll ) {
40#else
41 if (filterPassed() ) {
42#endif
43 ATH_MSG_DEBUG("Event passed filter");
44 m_nPass += 1;
45 } else {
46 ATH_MSG_DEBUG("Event failed filter");
47 m_nFail += 1;
48 }
49 // Bail out once we have enough events
50 if (m_nPass >= m_nNeeded && m_nNeeded > 0)
51 sc = StatusCode::FAILURE;
52 return sc;
53}
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static Double_t sc
int m_nPass
Definition GenFilter.h:65
int m_nNeeded
Definition GenFilter.h:67
int m_nFail
Definition GenFilter.h:66
virtual StatusCode filterEvent()=0
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ filterEvent()

StatusCode TTbarWToLeptonFilter::filterEvent ( )
virtual

Implements GenFilter.

Definition at line 18 of file TTbarWToLeptonFilter.cxx.

18 {
19 int N_quark_t = 0;
20 int N_quark_tbar = 0;
21 int N_quark_t_all = 0;
22 int N_quark_tbar_all = 0;
23 int N_pt_above_cut = 0;
24 int N_pt_above_cut_plus = 0;
25 int N_pt_above_cut_minus = 0;
26
27 int foundlepton[6] = {0, 0, 0, 0, 0, 0}; // e+, e-, mu+, mu-, tau+, tau-
28 int count_found_leptons = 1; // In ttbar each charged lepton flavour is counted at most once
29 if(m_fourTopsFilter) count_found_leptons = 2; // In four tops, one can have the same charged lepton flavour twice
30
31 for (McEventCollection::const_iterator itr = events()->begin(); itr!=events()->end(); ++itr) {
32 const HepMC::GenEvent* genEvt = (*itr);
33#ifdef HEPMC3
34 for (const auto& pitr: *genEvt) {
35 if (std::abs(pitr->pdg_id()) != 6) continue;
36 if ( pitr->pdg_id() == 6 ) N_quark_t_all++;
37 if ( pitr->pdg_id() == -6 ) N_quark_tbar_all++;
38 auto decayVtx = pitr->end_vertex();
39 // Verify if we got a valid pointer and retrieve the number of daughters
40 if (!decayVtx) continue;
41 // For this analysis we are not interested in t->t MC structures, only in decays
42 if (decayVtx->particles_out().size() < 2) continue;
43 for (const auto& child_mcpart: decayVtx->particles_out()) {
44 // Implicitly assume that tops always decay to W X
45 if (std::abs(child_mcpart->pdg_id()) != 24) continue;
46 if ( pitr->pdg_id() == 6 ) N_quark_t++;
47 if ( pitr->pdg_id() == -6 ) N_quark_tbar++;
48
49 bool useNextVertex = false;
50 auto w_decayVtx = child_mcpart->end_vertex();
51 while (w_decayVtx) {
52 useNextVertex = false;
53 for (const auto& grandchild_mcpart: *w_decayVtx) {
54 int grandchild_pid = grandchild_mcpart->pdg_id();
55 ATH_MSG_DEBUG("W (t/tbar) has " << w_decayVtx->particles_out().size() << " children and the pdg_id of the next is " << grandchild_pid);
56 // Check if the W's child is W again. If yes, then move to its next decay vertex in a decay tree
57 if (std::abs(grandchild_pid) == 24) {
58 w_decayVtx = grandchild_mcpart->end_vertex();
59 // If something wrong comes from truth...
60 if (!w_decayVtx) {
61 ATH_MSG_ERROR("A stable W is found... ");
62 break;
63 }
64 useNextVertex = true;
65 break;
66 }
67
68 // use brute force to use only leptons that have not been found already
69 if (grandchild_pid == -11 && foundlepton[0] < count_found_leptons) {
70 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
71 foundlepton[0]++;
72 N_pt_above_cut++;
73 N_pt_above_cut_minus++;
74 }
75 }
76 if (grandchild_pid == 11 && foundlepton[1] < count_found_leptons) {
77 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
78 foundlepton[1]++;
79 N_pt_above_cut++;
80 N_pt_above_cut_plus++;
81
82 }
83 }
84 if (grandchild_pid == -13 && foundlepton[2] < count_found_leptons) {
85 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
86 foundlepton[2]++;
87 N_pt_above_cut++;
88 N_pt_above_cut_minus++;
89 }
90 }
91 if (grandchild_pid == 13 && foundlepton[3] < count_found_leptons) {
92 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
93 foundlepton[3]++;
94 N_pt_above_cut++;
95 N_pt_above_cut_plus++;
96 }
97 }
98 if (grandchild_pid == -15 && foundlepton[4] < count_found_leptons) {
99 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
100 foundlepton[4]++;
101 N_pt_above_cut++;
102 N_pt_above_cut_minus++;
103 }
104 }
105 if (grandchild_pid == 15 && foundlepton[5] < count_found_leptons) {
106 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
107 foundlepton[5]++;
108 N_pt_above_cut++;
109 N_pt_above_cut_plus++;
110 }
111 }
112 }
113 // If investigation of W's next decay vertex is not required then finish looking for leptons
114 if (!useNextVertex) break;
115 }
116 }
117 }
118#else
119 for (HepMC::GenEvent::particle_const_iterator pitr = genEvt->particles_begin(); pitr != genEvt->particles_end(); ++pitr) {
120 if (std::abs((*pitr)->pdg_id()) == 6) {
121 if ( (*pitr)->pdg_id() == 6 ) N_quark_t_all++;
122 if ( (*pitr)->pdg_id() == -6 ) N_quark_tbar_all++;
123
124 int n_daughters = 0;
125
126 HepMC::GenParticle * mcpart = (*pitr);
127 const HepMC::GenVertex * decayVtx = mcpart->end_vertex();
128
129 // Verify if we got a valid pointer and retrieve the number of daughters
130 if (decayVtx != 0) n_daughters = decayVtx->particles_out_size();
131
132 // For this analysis we are not interested in t->t MC structures, only in decays
133 if (n_daughters >= 2) {
134 HepMC::GenVertex::particles_in_const_iterator child_mcpartItr = decayVtx->particles_out_const_begin();
135 HepMC::GenVertex::particles_in_const_iterator child_mcpartItrE = decayVtx->particles_out_const_end();
136 for (; child_mcpartItr != child_mcpartItrE; ++child_mcpartItr) {
137 HepMC::GenParticle * child_mcpart = (*child_mcpartItr);
138
139 // Implicitly assume that tops always decay to W X
140 if (std::abs(child_mcpart->pdg_id()) == 24) {
141 if ( (*pitr)->pdg_id() == 6 ) N_quark_t++;
142 if ( (*pitr)->pdg_id() == -6 ) N_quark_tbar++;
143
144 bool useNextVertex = false;
145 const HepMC::GenVertex * w_decayVtx = child_mcpart->end_vertex();
146
147 while (w_decayVtx) {
148
149 useNextVertex = false;
150 int mcpart_n_particles_out = w_decayVtx->particles_out_size();
151
152 HepMC::GenVertex::particles_out_const_iterator grandchild_mcpartItr = w_decayVtx->particles_out_const_begin();
153 HepMC::GenVertex::particles_out_const_iterator grandchild_mcpartItrE = w_decayVtx->particles_out_const_end();
154
155 for (; grandchild_mcpartItr != grandchild_mcpartItrE; ++grandchild_mcpartItr) {
156 HepMC::GenParticle * grandchild_mcpart = (*grandchild_mcpartItr);
157 int grandchild_pid = grandchild_mcpart->pdg_id();
158
159 ATH_MSG_DEBUG("W (t/tbar) has " << mcpart_n_particles_out << " children and the pdg_id of the next is " << grandchild_pid);
160
161 // Check if the W's child is W again. If yes, then move to its next decay vertex in a decay tree
162 if (std::abs(grandchild_pid) == 24) {
163 w_decayVtx = grandchild_mcpart->end_vertex();
164
165 // If something wrong comes from truth...
166 if (!w_decayVtx) {
167 ATH_MSG_ERROR("A stable W is found... ");
168 break;
169 }
170
171 useNextVertex = true;
172 break;
173 }
174
175 // use brute force to use only leptons that have not been found already
176 if (grandchild_pid == -11 && foundlepton[0] < count_found_leptons) {
177 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
178 foundlepton[0]++;
179 N_pt_above_cut++;
180 N_pt_above_cut_minus++;
181 }
182 }
183 if (grandchild_pid == 11 && foundlepton[1] < count_found_leptons) {
184 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
185 foundlepton[1]++;
186 N_pt_above_cut++;
187 N_pt_above_cut_plus++;
188
189 }
190 }
191 if (grandchild_pid == -13 && foundlepton[2] < count_found_leptons) {
192 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
193 foundlepton[2]++;
194 N_pt_above_cut++;
195 N_pt_above_cut_minus++;
196 }
197 }
198 if (grandchild_pid == 13 && foundlepton[3] < count_found_leptons) {
199 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
200 foundlepton[3]++;
201 N_pt_above_cut++;
202 N_pt_above_cut_plus++;
203 }
204 }
205 if (grandchild_pid == -15 && foundlepton[4] < count_found_leptons) {
206 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
207 foundlepton[4]++;
208 N_pt_above_cut++;
209 N_pt_above_cut_minus++;
210 }
211 }
212 if (grandchild_pid == 15 && foundlepton[5] < count_found_leptons) {
213 if (grandchild_mcpart->momentum().perp() >= m_Ptmin) {
214 foundlepton[5]++;
215 N_pt_above_cut++;
216 N_pt_above_cut_plus++;
217 }
218 }
219
220 }
221
222 // If investigation of W's next decay vertex is not required then finish looking for leptons
223 if (!useNextVertex) break;
224 }
225 }
226 }
227 }
228 }
229 }
230#endif
231 }
232
233 ATH_MSG_INFO("Found " << N_quark_t_all << " t quarks in event record");
234 ATH_MSG_INFO("Found " << N_quark_tbar_all << " tbar quarks in event record");
235 ATH_MSG_INFO("Found " << N_quark_t << " t -> W X decays");
236 ATH_MSG_INFO("Found " << N_quark_tbar << " tbar -> W X decays");
237 ATH_MSG_INFO("Num leptons from W decays from tops with lepton pt above " << m_Ptmin/1000.0 << " Gaudi::Units::GeV/c = " << N_pt_above_cut);
238
239 int count_tops = 1;
240 if(m_fourTopsFilter) count_tops = 2;
241 if (N_quark_t_all < count_tops || N_quark_tbar_all < count_tops) {
242
243 ATH_MSG_ERROR("No t or tbar quarks were found in a (presumably) ttbar event! Event is rejected.");
244 setFilterPassed(false);
245 return StatusCode::SUCCESS;
246 }
247
248 if (N_quark_t < count_tops || N_quark_tbar < count_tops) {
249
250 ATH_MSG_ERROR("No t or tbar quarks were found decaying to W in a (presumably) ttbar event! Event is rejected. Event dump follows.");
251 int event = 0;
252 for (McEventCollection::const_iterator itr = events()->begin(); itr!=events()->end(); ++itr) {
253 event++;
254 const HepMC::GenEvent* genEvt = (*itr);
255 int part=0 ;
256 for (const auto& mcpart: *genEvt) {
257 part++;
258 int pid = mcpart->pdg_id();
259 ATH_MSG_ERROR("In event (from MC collection) " << event << " particle number " << part << " has pdg_id = " << pid);
260 // retrieve decay vertex
261 auto decayVtx = mcpart->end_vertex();
262 // verify if we got a valid pointer
263 if ( !decayVtx ) continue;
264 int part_child=0;
265 for (const auto& child_mcpart: *decayVtx) {
266 part_child++;
267 int child_pid = child_mcpart->pdg_id();
268 ATH_MSG_ERROR(" child " << part_child << " with pdg_id = " << child_pid);
269 }
270 }
271 }
272 setFilterPassed(false);
273 return StatusCode::SUCCESS;
274 }
275
276 // If NumLeptons is negative (default), accept if Nlep > 0, otherwise only accept an exact match
277 if (m_numLeptons < 0) {
278 if ( (N_quark_t >= 2 || N_quark_tbar >= 2) && !m_fourTopsFilter) {
279 ATH_MSG_WARNING("More than one t -> W X or tbar -> W X decays found. Event is accepted anyway.");
280 }
281 if ( (N_quark_t > 2 || N_quark_tbar > 2) && m_fourTopsFilter) {
282 ATH_MSG_WARNING("More than one t -> W X or tbar -> W X decays found. Event is accepted anyway.");
283 }
284 setFilterPassed(N_pt_above_cut > 0);
285 } else {
286 if(m_fourTopsFilter){
287 if(m_SSMLFilter) setFilterPassed( (N_pt_above_cut >= m_numLeptons) && (N_pt_above_cut_plus >= 2 || N_pt_above_cut_minus >= 2));
288 else setFilterPassed(N_pt_above_cut >= m_numLeptons);}
289 else setFilterPassed(N_pt_above_cut == m_numLeptons);
290 }
291
292 return StatusCode::SUCCESS;
293}
#define ATH_MSG_INFO(x)
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838

◆ filterFinalize()

◆ filterInitialize()

virtual StatusCode GenFilter::filterInitialize ( )
inlinevirtualinherited

Reimplemented in BoostedHadTopAndTopPair, ChargedTracksWeightFilter, DecaysFinalStateFilter, DecayTimeFilter, DiBjetFilter, DiLeptonMassFilter, DiPhotonFilter, DirectPhotonFilter, FourLeptonInvMassFilter, FourLeptonMassFilter, HeavyFlavorHadronFilter, HTFilter, JetFilter, LeptonFilter, LeptonPairFilter, M4MuIntervalFilter, MassRangeFilter, MuDstarFilter, MultiBjetFilter, MultiCjetFilter, MultiParticleFilter, ParentChildFilter, ParentChildwStatusFilter, ParentTwoChildrenFilter, ParticleDecayFilter, ParticleFilter, QCDTruthJetFilter, QCDTruthMultiJetFilter, SameParticleHardScatteringFilter, TauFilter, TrimuMassRangeFilter, TruthJetFilter, TTbarPlusHeavyFlavorFilter, TTbarWithJpsimumuFilter, VBFForwardJetsFilter, VBFHbbEtaSortingFilter, VBFMjjIntervalFilter, xAODBSignalFilter, xAODChargedTracksFilter, xAODChargedTracksWeightFilter, xAODDecaysFinalStateFilter, xAODDecayTimeFilter, xAODDiLeptonMassFilter, xAODDiPhotonFilter, xAODDirectPhotonFilter, xAODElectronFilter, xAODForwardProtonFilter, xAODFourLeptonMassFilter, xAODHeavyFlavorHadronFilter, xAODHTFilter, xAODJetFilter, xAODLeptonFilter, xAODLeptonPairFilter, xAODM4MuIntervalFilter, xAODMETFilter, xAODMuDstarFilter, xAODMultiBjetFilter, xAODMultiCjetFilter, xAODMultiElecMuTauFilter, xAODMultiElectronFilter, xAODMultiLeptonFilter, xAODMultiMuonFilter, xAODMuonFilter, xAODParentChildFilter, xAODParentTwoChildrenFilter, xAODParticleDecayFilter, xAODParticleFilter, xAODPhotonFilter, xAODSameParticleHardScatteringFilter, xAODSplitPhotonFilter, xAODTauFilter, xAODTTbarWithJpsimumuFilter, xAODTTbarWToLeptonFilter, xAODVBFForwardJetsFilter, xAODVBFMjjIntervalFilter, xAODXtoVVDecayFilterExtended, and XtoVVDecayFilterExtended.

Definition at line 45 of file GenFilter.h.

45{ return StatusCode::SUCCESS; }

◆ finalize()

StatusCode GenFilter::finalize ( )
inherited

Definition at line 56 of file GenFilter.cxx.

56 {
57 ATH_MSG_INFO("Events passed = " << m_nPass << " Events failed = " << m_nFail);
59 return StatusCode::SUCCESS;
60}
#define CHECK(...)
Evaluate an expression and check for errors.
virtual StatusCode filterFinalize()
Definition GenFilter.h:47

◆ GeVToMeV()

void GenBase::GeVToMeV ( HepMC::GenEvent * evt)
protectedinherited

Scale event energies/momenta by x 1000.

Todo
Add HepMC units awareness and do it differently when HepMC provides this functionality directly (and reference-based FourVector accessors)

Definition at line 58 of file GenBase.cxx.

58 {
59 for (HepMC::GenEvent::particle_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p) {
60 const HepMC::FourVector fv((*p)->momentum().px() * 1000,
61 (*p)->momentum().py() * 1000,
62 (*p)->momentum().pz() * 1000,
63 (*p)->momentum().e() * 1000);
64 (*p)->set_momentum(fv);
65 (*p)->set_generated_mass(1000 * (*p)->generated_mass());
66 }
67}

◆ initialize()

StatusCode GenFilter::initialize ( )
virtualinherited

Reimplemented from GenBase.

Definition at line 20 of file GenFilter.cxx.

20 {
22 m_nPass = 0;
23 m_nFail = 0;
25 return StatusCode::SUCCESS;
26}
virtual StatusCode initialize() override
Definition GenBase.cxx:17
virtual StatusCode filterInitialize()
Definition GenFilter.h:45

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ MeVToGeV()

void GenBase::MeVToGeV ( HepMC::GenEvent * evt)
protectedinherited

Scale event energies/momenta by x 1/1000.

Definition at line 68 of file GenBase.cxx.

68 {
69 for (HepMC::GenEvent::particle_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p) {
70 const HepMC::FourVector fv((*p)->momentum().px() / 1000,
71 (*p)->momentum().py() / 1000,
72 (*p)->momentum().pz() / 1000,
73 (*p)->momentum().e() / 1000);
74 (*p)->set_momentum(fv);
75 (*p)->set_generated_mass((*p)->generated_mass() / 1000);
76 }
77}

◆ mmTocm()

void GenBase::mmTocm ( HepMC::GenEvent * evt)
protectedinherited

Scale event lengths by x 1/10.

Definition at line 87 of file GenBase.cxx.

87 {
88 for (HepMC::GenEvent::vertex_iterator vtx = evt->vertices_begin(); vtx != evt->vertices_end(); ++vtx) {
89 const HepMC::FourVector fv((*vtx)->position().x() / 10,
90 (*vtx)->position().y() / 10,
91 (*vtx)->position().z() / 10,
92 (*vtx)->position().t() / 10);
93 (*vtx)->set_position(fv);
94 }
95}

◆ msg()

MsgStream & AthCommonMsg< Algorithm >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< Algorithm >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ nFailed()

int GenFilter::nFailed ( ) const
inlineinherited

Definition at line 53 of file GenFilter.h.

53{ return m_nFail; }

◆ nNeeded()

int GenFilter::nNeeded ( ) const
inlineinherited

Definition at line 54 of file GenFilter.h.

54{ return m_nNeeded; }

◆ nPassed()

int GenFilter::nPassed ( ) const
inlineinherited

Definition at line 52 of file GenFilter.h.

52{ return m_nPass; }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ particleData()

const HepPDT::ParticleData * GenBase::particleData ( int pid) const
inlineinherited

Access an element in the particle data table.

Definition at line 126 of file GenBase.h.

126 {
127 return pdt().particle(HepPDT::ParticleID(std::abs(pid)));
128 }
const HepPDT::ParticleDataTable & pdt() const
Shorter alias to get a particle data table.
Definition GenBase.h:123

◆ particleTable()

const HepPDT::ParticleDataTable & GenBase::particleTable ( ) const
inlineinherited

Get a particle data table.

Definition at line 118 of file GenBase.h.

118 {
119 return *(m_ppSvc->PDT());
120 }
ServiceHandle< IPartPropSvc > m_ppSvc
Handle on the particle property service.
Definition GenBase.h:160

◆ partPropSvc()

const ServiceHandle< IPartPropSvc > GenBase::partPropSvc ( ) const
inlineinherited

Access the particle property service.

Definition at line 113 of file GenBase.h.

113 {
114 return m_ppSvc;
115 }

◆ pdt()

const HepPDT::ParticleDataTable & GenBase::pdt ( ) const
inlineinherited

Shorter alias to get a particle data table.

Definition at line 123 of file GenBase.h.

123{ return particleTable(); }
const HepPDT::ParticleDataTable & particleTable() const
Get a particle data table.
Definition GenBase.h:118

◆ renounce()

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 > AthCommonDataStore< AthCommonMsg< Algorithm > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
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)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_fourTopsFilter

bool TTbarWToLeptonFilter::m_fourTopsFilter
private

Definition at line 26 of file TTbarWToLeptonFilter.h.

◆ m_mcEventKey

std::string GenBase::m_mcEventKey {}
protectedinherited

StoreGate key for the MC event collection (defaults to GEN_EVENT)

Definition at line 137 of file GenBase.h.

137{};

◆ m_mcevents_const

SG::ReadHandleKey<McEventCollection> GenBase::m_mcevents_const { this, "McEventKey", "GEN_EVENT", "StoreGate key of the MC event collection" }
privateinherited

Const handle to the MC event collection.

Definition at line 163 of file GenBase.h.

163{ this, "McEventKey", "GEN_EVENT", "StoreGate key of the MC event collection" };

◆ m_mkMcEvent

BooleanProperty GenBase::m_mkMcEvent {this, "MakeMcEvent", false, "Create a new MC event collection if it doesn't exist"}
protectedinherited

Flag to determine if a new MC event collection should be made if it doesn't exist.

Definition at line 139 of file GenBase.h.

139{this, "MakeMcEvent", false, "Create a new MC event collection if it doesn't exist"};

◆ m_nFail

int GenFilter::m_nFail
protectedinherited

Definition at line 66 of file GenFilter.h.

◆ m_nNeeded

int GenFilter::m_nNeeded
protectedinherited

Definition at line 67 of file GenFilter.h.

◆ m_nPass

int GenFilter::m_nPass
protectedinherited

Definition at line 65 of file GenFilter.h.

◆ m_numLeptons

int TTbarWToLeptonFilter::m_numLeptons
private

Definition at line 25 of file TTbarWToLeptonFilter.h.

◆ m_ppSvc

ServiceHandle<IPartPropSvc> GenBase::m_ppSvc {this, "PartPropSvc", "PartPropSvc"}
privateinherited

Handle on the particle property service.

Definition at line 160 of file GenBase.h.

160{this, "PartPropSvc", "PartPropSvc"};

◆ m_Ptmin

double TTbarWToLeptonFilter::m_Ptmin
private

Definition at line 24 of file TTbarWToLeptonFilter.h.

◆ m_SSMLFilter

bool TTbarWToLeptonFilter::m_SSMLFilter
private

Definition at line 27 of file TTbarWToLeptonFilter.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< Algorithm > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< Algorithm > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.

◆ Ptcut

TTbarWToLeptonFilter.Ptcut

Add this filter to the algs required to be successful for streaming if "TTbarWToLeptonFilter" not in StreamEVGEN.RequireAlgs: StreamEVGEN.RequireAlgs += ["TTbarWToLeptonFilter"].

Default cut params

Definition at line 15 of file TTbarWToLeptonFilter.py.


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