ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::FourLeptonVertexingAlgorithm Class Reference

#include <FourLeptonVertexingAlgorithm.h>

Inheritance diagram for DerivationFramework::FourLeptonVertexingAlgorithm:
Collaboration diagram for DerivationFramework::FourLeptonVertexingAlgorithm:

Public Types

using LeptonQuadruplet = std::array<const xAOD::IParticle*, 4>

Public Member Functions

 FourLeptonVertexingAlgorithm (const std::string &n, ISvcLocator *p)
StatusCode initialize () override
StatusCode execute (const EventContext &ctx) const override
StatusCode finalize () override
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual bool isClonable () const override
 Specify if the algorithm is clonable.
virtual unsigned int cardinality () const override
 Cardinality (Maximum number of clones that can exist) special value 0 means that algorithm is reentrant.
virtual StatusCode sysExecute (const EventContext &ctx) override
 Execute an algorithm.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
virtual bool filterPassed (const EventContext &ctx) const
virtual void setFilterPassed (bool state, const EventContext &ctx) const
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

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.

Private Types

using MuonTrk = xAOD::Muon::TrackParticleType
typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

bool passSelection (const xAOD::Electron *elec) const
bool passSelection (const xAOD::Muon *muon) const
bool passSelection (const LeptonQuadruplet &quad) const
std::vector< LeptonQuadrupletbuildAllQuadruplets (const EventContext &ctx) const
std::unique_ptr< xAOD::VertexfitQuadruplet (const EventContext &ctx, const LeptonQuadruplet &quad) const
const xAOD::TrackParticletrackParticle (const xAOD::IParticle *part) const
int charge (const xAOD::IParticle *part) const
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

SG::ReadHandleKey< xAOD::MuonContainerm_muonKey
 Input containers.
SG::ReadHandleKey< xAOD::ElectronContainerm_elecKey
SG::ReadHandleKey< xAOD::EventInfom_evtKey {this, "EventInfoKey", "EventInfo"}
ToolHandle< CP::IMuonSelectionToolm_muonSelTool
ToolHandle< IAsgElectronLikelihoodToolm_elecSelTool
ToolHandle< Trk::IVertexFitterm_fitter {this, "VertexFitter", "", "Vertex fitter tool"}
Gaudi::Property< float > m_minMuonPt
Gaudi::Property< float > m_minElecPt {this, "MinElecPt", 4.5 * Gaudi::Units::GeV, " Minimum pt cut applied on the electron"}
Gaudi::Property< int > m_muonTrkProp
MuonTrk m_muonTrk {MuonTrk::Primary}
Gaudi::Property< bool > m_elecUseGSF
Gaudi::Property< float > m_z0Cut
Gaudi::Property< float > m_lowSFOS_Cut {this, "LowSFOSMass", 3.5 * Gaudi::Units::GeV, "Minimal mass for the lower SFOS pair"}
Gaudi::Property< float > m_highPair_Cut {this, "HighMassPair", 40.*Gaudi::Units::GeV, "Minimum mass for the largest lepton pair"}
Gaudi::Property< float > m_LeadPtCut {this, "LeadingLeptonPt", 15. *Gaudi::Units::GeV,"Minimal momentum of the leading lepton in the quad "}
Gaudi::Property< float > m_SubLeadPtCut {this, "SubLeadingLeptonPt", 10. *Gaudi::Units::GeV,"Minimal momentum of the leading lepton in the quad "}
Gaudi::Property< float > m_VtxChi2Cut {this, "VertexChi2Cut", 20, "Maximal chii n.D.o.F cut on the reconstructed vertex" }
SG::WriteHandleKey< xAOD::VertexContainerm_vtxKey {this, "OutVertexContainer", "FourLeptonVertices", "Output vertex container"}
Gaudi::Property< bool > m_pruneCov {this, "PruneCovariance", true, "Clears the covariance vector"}
Gaudi::Property< bool > m_pruneWeight {this, "PruneTrackWeights", true, "Clears the track weight vector"}
std::array< std::atomic< unsigned int >, 8 > m_num_lep ATLAS_THREAD_SAFE {}
 Simple counter to evaluate the number of leptons per event.
std::array< std::atomic< unsigned int >, 8 > m_num_ele ATLAS_THREAD_SAFE {}
std::array< std::atomic< unsigned int >, 8 > m_num_muo ATLAS_THREAD_SAFE {}
std::atomic< unsigned int > m_tried_fits {0}
 How many vertex fits were attempted.
std::atomic< unsigned int > m_good_fits {0}
 How many vertex fits were actually successful.
DataObjIDColl m_extendedExtraObjects
 Extra output dependency collection, extended by AthAlgorithmDHUpdate to add symlinks.
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

Detailed Description

Definition at line 26 of file FourLeptonVertexingAlgorithm.h.

Member Typedef Documentation

◆ LeptonQuadruplet

◆ MuonTrk

using DerivationFramework::FourLeptonVertexingAlgorithm::MuonTrk = xAOD::Muon::TrackParticleType
private

Definition at line 65 of file FourLeptonVertexingAlgorithm.h.

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ FourLeptonVertexingAlgorithm()

DerivationFramework::FourLeptonVertexingAlgorithm::FourLeptonVertexingAlgorithm ( const std::string & n,
ISvcLocator * p )

Definition at line 43 of file FourLeptonVertexingAlgorithm.cxx.

43: AthReentrantAlgorithm(n, p) {}

Member Function Documentation

◆ buildAllQuadruplets()

std::vector< LeptonQuadruplet > DerivationFramework::FourLeptonVertexingAlgorithm::buildAllQuadruplets ( const EventContext & ctx) const
private

Retrieve the container from the store gate

Dump the leptons in a common vector to build all possible quadruplets

Increment the lepton counter

None of the leptons has enough momentum to be selected in the offline analysis

The leading lepton has too little momentum. All subsequent combinations will have even less. Bail out of the loop

Definition at line 95 of file FourLeptonVertexingAlgorithm.cxx.

95 {
96 std::vector<LeptonQuadruplet> to_ret{};
97
99 SG::ReadHandle<xAOD::MuonContainer> muonContainer{m_muonKey, ctx};
100 if (!muonContainer.isValid()) {
101 ATH_MSG_FATAL("Failed to retrieve the muon container for input " << m_muonKey.fullKey());
102 throw std::runtime_error("Invalid key access");
103 }
104 SG::ReadHandle<xAOD::ElectronContainer> elecContainer{m_elecKey, ctx};
105 if (!elecContainer.isValid()) {
106 ATH_MSG_FATAL("Failed to retrieve the electron container for input " << m_elecKey.fullKey());
107 throw std::runtime_error("Invalid key access");
108 }
110 using PartVec = std::vector<const xAOD::IParticle*>;
111 PartVec selected_lep{};
112 const size_t num_lep = elecContainer->size() + muonContainer->size();
113 if (num_lep < 4) {
114 ATH_MSG_DEBUG("Less than four objects reconstructed");
115 ++m_num_lep[num_lep];
116 ++m_num_ele[elecContainer->size()];
117 ++m_num_muo[muonContainer->size()];
118 return to_ret;
119 }
120 selected_lep.reserve(num_lep);
121 size_t n_ele{0},n_muo{0};
122 for (const xAOD::Muon* muon : *muonContainer) {
123 if (passSelection(muon)) {
124 ATH_MSG_DEBUG("Add " << muon);
125 selected_lep.push_back(muon);
126 ++n_ele;
127 }
128 }
129 for (const xAOD::Electron* elec : *elecContainer) {
130 if (passSelection(elec)) {
131 ATH_MSG_DEBUG("Add " << elec);
132 selected_lep.push_back(elec);
133 ++n_muo;
134 }
135 }
137 ++m_num_lep[std::min(m_num_lep.size() -1 , n_ele + n_muo)];
138 ++m_num_ele[std::min(m_num_ele.size() -1 , n_ele)];
139 ++m_num_muo[std::min(m_num_muo.size() -1 , n_muo)];
140
141 if (selected_lep.size() < 4) {
142 ATH_MSG_DEBUG("Less than four leptons survived.");
143 return to_ret;
144 }
145 std::sort(selected_lep.begin(),selected_lep.end(),
146 [](const xAOD::IParticle* a, const xAOD::IParticle* b) {
147 return a->pt() > b->pt();});
149 if (selected_lep[0]->pt() < m_LeadPtCut || selected_lep[1]->pt() < m_SubLeadPtCut) return to_ret;
150
151 to_ret.reserve(std::pow(selected_lep.size(), 4));
152 for (PartVec::const_iterator itr1 = selected_lep.begin() + 3; itr1 != selected_lep.end(); ++itr1) {
153 for (PartVec::const_iterator itr2 = selected_lep.begin() + 2; itr2 != itr1; ++itr2) {
154 for (PartVec::const_iterator itr3 = selected_lep.begin() + 1; itr3 != itr2; ++itr3) {
155 if ( (*itr3)->pt() < m_SubLeadPtCut) return to_ret;
156 for (PartVec::const_iterator itr4 = selected_lep.begin(); itr4 != itr3; ++itr4) {
157 LeptonQuadruplet quad{*itr4, *itr3, *itr2, *itr1};
160 if (quad[0]->pt() < m_LeadPtCut) return to_ret;
161 ATH_MSG_DEBUG("Create new quaduplet " << quad);
162 to_ret.push_back(std::move(quad));
163 }
164 }
165 }
166 }
167 return to_ret;
168 }
#define ATH_MSG_FATAL(x)
#define ATH_MSG_DEBUG(x)
static Double_t a
xAOD::MuonContainer * muonContainer
size_type size() const noexcept
Returns the number of elements in the collection.
SG::ReadHandleKey< xAOD::ElectronContainer > m_elecKey
std::array< const xAOD::IParticle *, 4 > LeptonQuadruplet
SG::ReadHandleKey< xAOD::MuonContainer > m_muonKey
Input containers.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
Muon_v1 Muon
Reference the current persistent version:
Electron_v1 Electron
Definition of the current "egamma version".

◆ cardinality()

unsigned int AthCommonReentrantAlgorithm< Gaudi::Algorithm >::cardinality ( ) const
overridevirtualinherited

Cardinality (Maximum number of clones that can exist) special value 0 means that algorithm is reentrant.

Override this to return 0 for reentrant algorithms.

Definition at line 75 of file AthCommonReentrantAlgorithm.cxx.

64{
65 return 0;
66}

◆ charge()

int DerivationFramework::FourLeptonVertexingAlgorithm::charge ( const xAOD::IParticle * part) const
private

Definition at line 176 of file FourLeptonVertexingAlgorithm.cxx.

176 {
177 static const SG::AuxElement::ConstAccessor<float> acc_charge{"charge"};
178 if (acc_charge.isAvailable(*part)) return acc_charge(*part);
179 return trackParticle(part)->charge();
180 }
const xAOD::TrackParticle * trackParticle(const xAOD::IParticle *part) const
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
float charge() const
Returns the charge.

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Gaudi::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< Gaudi::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< Gaudi::Algorithm > >::detStore ( ) const
inlineinherited

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

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Gaudi::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 DerivationFramework::FourLeptonVertexingAlgorithm::execute ( const EventContext & ctx) const
override

Setup the output container

Definition at line 74 of file FourLeptonVertexingAlgorithm.cxx.

74 {
76 SG::WriteHandle<xAOD::VertexContainer> vtxContainer{m_vtxKey, ctx};
77 ATH_CHECK(vtxContainer.record(std::make_unique<xAOD::VertexContainer>(), std::make_unique<xAOD::VertexAuxContainer>()));
78 xAOD::VertexContainer* out_container = vtxContainer.ptr();
79
80 std::vector<LeptonQuadruplet> all_quads = buildAllQuadruplets(ctx);
81 for (const LeptonQuadruplet& quad : all_quads) {
82 std::unique_ptr<xAOD::Vertex> candidate = fitQuadruplet(ctx, quad);
83 if (candidate) out_container->push_back(std::move(candidate));
84 }
85 return StatusCode::SUCCESS;
86 }
#define ATH_CHECK
Evaluate an expression and check for errors.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
SG::WriteHandleKey< xAOD::VertexContainer > m_vtxKey
std::vector< LeptonQuadruplet > buildAllQuadruplets(const EventContext &ctx) const
std::unique_ptr< xAOD::Vertex > fitQuadruplet(const EventContext &ctx, const LeptonQuadruplet &quad) const
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
pointer_type ptr()
Dereference the pointer.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Gaudi::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 & AthCommonReentrantAlgorithm< Gaudi::Algorithm >::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

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

Definition at line 94 of file AthCommonReentrantAlgorithm.cxx.

90{
91 // If we didn't find any symlinks to add, just return the collection
92 // from the base class. Otherwise, return the extended collection.
93 if (!m_extendedExtraObjects.empty()) {
95 }
97}
An algorithm that can be simultaneously executed in multiple threads.

◆ filterPassed()

virtual bool AthCommonReentrantAlgorithm< Gaudi::Algorithm >::filterPassed ( const EventContext & ctx) const
inlinevirtualinherited

Definition at line 96 of file AthCommonReentrantAlgorithm.h.

96 {
97 return execState( ctx ).filterPassed();
98 }
virtual bool filterPassed(const EventContext &ctx) const

◆ finalize()

StatusCode DerivationFramework::FourLeptonVertexingAlgorithm::finalize ( )
override

Definition at line 62 of file FourLeptonVertexingAlgorithm.cxx.

62 {
63 unsigned int n_processed = std::accumulate(m_num_lep.begin(), m_num_lep.end(), 0, [](const std::atomic<unsigned int>& l, unsigned int n){
64 return n +l;
65 });
66 ATH_MSG_INFO("Proccessed in total "<<n_processed<<" events. Lepton multiplicities observed:");
67 for (unsigned int l = 0 ; l < m_num_lep.size(); ++l) {
68 ATH_MSG_INFO(" ---- Events with "<<l<<" leptons: "<<m_num_lep[l]<<" (elec/muon) "<<m_num_ele[l]<<"/"<<m_num_muo[l]<<" ("<<(100.*m_num_lep[l] / std::max(1u,n_processed))<<"%)" );
69 }
70 const unsigned int n_trials = m_tried_fits;
71 ATH_MSG_INFO("---> Attempted fits "<<m_tried_fits<<" out of which "<<m_good_fits<<" ("<<(100.* m_good_fits / std::max(1u, n_trials) )<<"%) succeeded." );
72 return StatusCode::SUCCESS;
73 }
#define ATH_MSG_INFO(x)
std::atomic< unsigned int > m_good_fits
How many vertex fits were actually successful.
std::atomic< unsigned int > m_tried_fits
How many vertex fits were attempted.
l
Printing final latex table to .tex output file.

◆ fitQuadruplet()

std::unique_ptr< xAOD::Vertex > DerivationFramework::FourLeptonVertexingAlgorithm::fitQuadruplet ( const EventContext & ctx,
const LeptonQuadruplet & quad ) const
private

Use the beam spot as an initial constraint

Definition at line 181 of file FourLeptonVertexingAlgorithm.cxx.

181 {
182 if (!passSelection(quad)) return nullptr;
183 ++m_tried_fits;
184 std::vector<const xAOD::TrackParticle*> trks{};
185 trks.reserve(4);
186 for (const xAOD::IParticle* lep : quad) { trks.push_back(trackParticle(lep)); }
188 SG::ReadHandle<xAOD::EventInfo> evtInfo{m_evtKey, ctx};
189 if (!evtInfo.isValid()) {
190 ATH_MSG_FATAL("Failed to retrieve the event info " << m_evtKey.fullKey());
191 throw std::runtime_error("Invalid key access");
192 }
193
194 const Amg::Vector3D beam_spot{evtInfo->beamPosX(), evtInfo->beamPosY(), evtInfo->beamPosZ()};
195 std::unique_ptr<xAOD::Vertex> common_vtx = m_fitter->fit(ctx, trks, beam_spot);
196 if (!common_vtx) {
197 ATH_MSG_DEBUG("Fit from " << quad << " failed.");
198 return nullptr;
199 }
200 const float redChi2 = (common_vtx->chiSquared() / common_vtx->numberDoF());
201 ATH_MSG_DEBUG("Fit from " << quad << " gave a vertex with position at " << common_vtx->position() << " with chi2 "
202 << redChi2 << " nDoF: " << common_vtx->numberDoF());
203 if (redChi2 > m_VtxChi2Cut) {
204 ATH_MSG_DEBUG("Reject due to bad chi2");
205 return nullptr;
206 }
207 static const std::vector<float> empty_vec{};
208 if (m_pruneWeight) common_vtx->setTrackWeights(empty_vec);
209 if (m_pruneCov) common_vtx->setCovariance(empty_vec);
210 std::vector<TrkLink> track_links{};
211 for (const xAOD::TrackParticle* trk : trks) {
212 const xAOD::TrackParticleContainer* cont = dynamic_cast<const xAOD::TrackParticleContainer*>(trk->container());
213 TrkLink link{*cont, trk->index(), ctx};
214 link.toPersistent();
215 track_links.push_back(std::move(link));
216 }
217 common_vtx->setTrackParticleLinks(track_links);
218 ++m_good_fits;
219 return common_vtx;
220 }
Eigen::Matrix< double, 3, 1 > Vector3D
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".

◆ initialize()

StatusCode DerivationFramework::FourLeptonVertexingAlgorithm::initialize ( )
override

Definition at line 45 of file FourLeptonVertexingAlgorithm.cxx.

45 {
46 ATH_CHECK(m_muonKey.initialize());
47 ATH_CHECK(m_elecKey.initialize());
48 ATH_CHECK(m_evtKey.initialize());
49 ATH_CHECK(m_vtxKey.initialize());
50 ATH_CHECK(m_fitter.retrieve());
51
52 if (!m_elecSelTool.empty()) ATH_CHECK(m_elecSelTool.retrieve());
53 if (!m_muonSelTool.empty()) ATH_CHECK(m_muonSelTool.retrieve());
54
55 if (m_muonTrkProp < MuonTrk::Primary || m_muonTrk > MuonTrk::MSOnlyExtrapolatedMuonSpectrometerTrackParticle) {
56 ATH_MSG_FATAL("A bogous muon track particle has been picked " << m_muonTrkProp);
57 return StatusCode::FAILURE;
58 }
59 m_muonTrk = static_cast<MuonTrk>(m_muonTrkProp.value());
60 return StatusCode::SUCCESS;
61 }
ToolHandle< IAsgElectronLikelihoodTool > m_elecSelTool

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Gaudi::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.

◆ isClonable()

◆ msg()

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

Definition at line 24 of file AthCommonMsg.h.

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

◆ msgLvl()

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

Definition at line 30 of file AthCommonMsg.h.

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

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Gaudi::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.

◆ passSelection() [1/3]

bool DerivationFramework::FourLeptonVertexingAlgorithm::passSelection ( const LeptonQuadruplet & quad) const
private

Check the di lepton masses in the quad

If the pair is SFOS check the mass

The high mass pair must at least satisfy a certain threshold

Definition at line 221 of file FourLeptonVertexingAlgorithm.cxx.

221 {
222 float max_m{-1.};
223 for (unsigned int i = 1; i < quad.size(); ++i) {
224 for (unsigned int j = 0; j < i; ++j) {
225 const xAOD::TrackParticle* trk_a = trackParticle(quad[i]);
226 const xAOD::TrackParticle* trk_b = trackParticle(quad[j]);
227
228 if (std::abs(trk_b->z0() - trk_a->z0()) > m_z0Cut) return false;
229
231 const float di_m1 = (quad[i]->p4() + quad[j]->p4()).M();
232 const bool isSFOS1 = charge(quad[i])*charge(quad[j]) < 0 && quad[i]->type() == quad[j]->type();
234 if (isSFOS1 && di_m1 < m_lowSFOS_Cut) return false;
235 max_m = std::max(di_m1, max_m);
236 for (unsigned int k = 1 ; k < quad.size(); ++k) {
237 if (k == i || k == j) continue;
238 for (unsigned int l = 0; l < k ; ++l) {
239 if (l == i || l == j) continue;
240 const float di_m2 = (quad[k]->p4() + quad[l]->p4()).M();
241 const bool isSFOS2 = charge(quad[k])*charge(quad[l]) < 0 && quad[k]->type() == quad[l]->type();
242 if (isSFOS2 && di_m2 < m_lowSFOS_Cut) return false;
243 max_m = std::max(max_m, di_m2);
244 }
245 }
246 }
247 }
249 if (m_highPair_Cut > max_m) return false;
250
251 return true;
252 }
float z0() const
Returns the parameter.

◆ passSelection() [2/3]

bool DerivationFramework::FourLeptonVertexingAlgorithm::passSelection ( const xAOD::Electron * elec) const
private

Definition at line 88 of file FourLeptonVertexingAlgorithm.cxx.

88 {
89 return elec->pt() >= m_minElecPt && (m_elecSelTool.empty() || m_elecSelTool->accept(elec)) &&
91 }
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Definition Egamma_v1.cxx:66
std::vector< const xAOD::TrackParticle * > getTrackParticlesVec(const xAOD::Egamma *eg, bool useBremAssoc=true, bool allParticles=true)
Return a list of all or only the best TrackParticle associated to the object.

◆ passSelection() [3/3]

bool DerivationFramework::FourLeptonVertexingAlgorithm::passSelection ( const xAOD::Muon * muon) const
private

Definition at line 92 of file FourLeptonVertexingAlgorithm.cxx.

92 {
93 return muon->pt() >= m_minMuonPt && (m_muonSelTool.empty() || m_muonSelTool->accept(*muon)) && muon->trackParticle(m_muonTrk);
94 }

◆ 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< Gaudi::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< Gaudi::Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ setFilterPassed()

virtual void AthCommonReentrantAlgorithm< Gaudi::Algorithm >::setFilterPassed ( bool state,
const EventContext & ctx ) const
inlinevirtualinherited

Definition at line 100 of file AthCommonReentrantAlgorithm.h.

100 {
102 }
virtual void setFilterPassed(bool state, const EventContext &ctx) const

◆ sysExecute()

StatusCode AthCommonReentrantAlgorithm< Gaudi::Algorithm >::sysExecute ( const EventContext & ctx)
overridevirtualinherited

Execute an algorithm.

We override this in order to work around an issue with the Algorithm base class storing the event context in a member variable that can cause crashes in MT jobs.

Definition at line 85 of file AthCommonReentrantAlgorithm.cxx.

77{
78 return BaseAlg::sysExecute (ctx);
79}

◆ sysInitialize()

StatusCode AthCommonReentrantAlgorithm< Gaudi::Algorithm >::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< Gaudi::Algorithm > >.

Reimplemented in HypoBase, and InputMakerBase.

Definition at line 61 of file AthCommonReentrantAlgorithm.cxx.

107 {
109
110 if (sc.isFailure()) {
111 return sc;
112 }
113
114 ServiceHandle<ICondSvc> cs("CondSvc",name());
115 for (auto h : outputHandles()) {
116 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
117 // do this inside the loop so we don't create the CondSvc until needed
118 if ( cs.retrieve().isFailure() ) {
119 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
120 return StatusCode::SUCCESS;
121 }
122 if (cs->regHandle(this,*h).isFailure()) {
124 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
125 << " with CondSvc");
126 }
127 }
128 }
129 return sc;
130}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Gaudi::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.

◆ trackParticle()

const xAOD::TrackParticle * DerivationFramework::FourLeptonVertexingAlgorithm::trackParticle ( const xAOD::IParticle * part) const
private

Definition at line 169 of file FourLeptonVertexingAlgorithm.cxx.

169 {
170 if (part->type() == xAOD::Type::Muon)
171 return static_cast<const xAOD::Muon*>(part)->trackParticle(m_muonTrk);
172 else if (part->type() == xAOD::Type::Electron)
173 return xAOD::EgammaHelpers::getTrackParticlesVec(static_cast<const xAOD::Electron*>(part), !m_elecUseGSF)[0];
174 return dynamic_cast<const xAOD::TrackParticle*>(part);
175 }
@ Muon
The object is a muon.
Definition ObjectType.h:48
@ Electron
The object is an electron.
Definition ObjectType.h:46

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< Gaudi::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 }

Member Data Documentation

◆ ATLAS_THREAD_SAFE [1/3]

std::array<std::atomic<unsigned int>, 8> m_num_lep DerivationFramework::FourLeptonVertexingAlgorithm::ATLAS_THREAD_SAFE {}
mutableprivate

Simple counter to evaluate the number of leptons per event.

Definition at line 91 of file FourLeptonVertexingAlgorithm.h.

91{};

◆ ATLAS_THREAD_SAFE [2/3]

std::array<std::atomic<unsigned int>, 8> m_num_ele DerivationFramework::FourLeptonVertexingAlgorithm::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 92 of file FourLeptonVertexingAlgorithm.h.

92{};

◆ ATLAS_THREAD_SAFE [3/3]

std::array<std::atomic<unsigned int>, 8> m_num_muo DerivationFramework::FourLeptonVertexingAlgorithm::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 93 of file FourLeptonVertexingAlgorithm.h.

93{};

◆ m_detStore

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

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_elecKey

SG::ReadHandleKey<xAOD::ElectronContainer> DerivationFramework::FourLeptonVertexingAlgorithm::m_elecKey
private
Initial value:
{this, "ElectronContainer", "Electrons",
"Location of the electron container in the store gate"}

Definition at line 50 of file FourLeptonVertexingAlgorithm.h.

50 {this, "ElectronContainer", "Electrons",
51 "Location of the electron container in the store gate"};

◆ m_elecSelTool

ToolHandle<IAsgElectronLikelihoodTool> DerivationFramework::FourLeptonVertexingAlgorithm::m_elecSelTool
private
Initial value:
{this, "ElectronSelectionTool", "",
"Applied preselection on the electrons. Will be ignored if empty"}

Definition at line 56 of file FourLeptonVertexingAlgorithm.h.

56 {this, "ElectronSelectionTool", "",
57 "Applied preselection on the electrons. Will be ignored if empty"};

◆ m_elecUseGSF

Gaudi::Property<bool> DerivationFramework::FourLeptonVertexingAlgorithm::m_elecUseGSF
private
Initial value:
{this, "PickElecGsfTrk", true,
"If this property is set to true it will pick up the electron GSF track"}

Definition at line 70 of file FourLeptonVertexingAlgorithm.h.

70 {this, "PickElecGsfTrk", true,
71 "If this property is set to true it will pick up the electron GSF track"};

◆ m_evtKey

SG::ReadHandleKey<xAOD::EventInfo> DerivationFramework::FourLeptonVertexingAlgorithm::m_evtKey {this, "EventInfoKey", "EventInfo"}
private

Definition at line 52 of file FourLeptonVertexingAlgorithm.h.

52{this, "EventInfoKey", "EventInfo"};

◆ m_evtStore

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

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthCommonReentrantAlgorithm< Gaudi::Algorithm >::m_extendedExtraObjects
privateinherited

Extra output dependency collection, extended by AthAlgorithmDHUpdate to add symlinks.

Empty if no symlinks were found.

Definition at line 114 of file AthCommonReentrantAlgorithm.h.

◆ m_fitter

ToolHandle<Trk::IVertexFitter> DerivationFramework::FourLeptonVertexingAlgorithm::m_fitter {this, "VertexFitter", "", "Vertex fitter tool"}
private

Definition at line 59 of file FourLeptonVertexingAlgorithm.h.

59{this, "VertexFitter", "", "Vertex fitter tool"};

◆ m_good_fits

std::atomic<unsigned int> DerivationFramework::FourLeptonVertexingAlgorithm::m_good_fits {0}
mutableprivate

How many vertex fits were actually successful.

Definition at line 98 of file FourLeptonVertexingAlgorithm.h.

98{0};

◆ m_highPair_Cut

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_highPair_Cut {this, "HighMassPair", 40.*Gaudi::Units::GeV, "Minimum mass for the largest lepton pair"}
private

Definition at line 78 of file FourLeptonVertexingAlgorithm.h.

78{this, "HighMassPair", 40.*Gaudi::Units::GeV, "Minimum mass for the largest lepton pair"};

◆ m_LeadPtCut

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_LeadPtCut {this, "LeadingLeptonPt", 15. *Gaudi::Units::GeV,"Minimal momentum of the leading lepton in the quad "}
private

Definition at line 80 of file FourLeptonVertexingAlgorithm.h.

80{this, "LeadingLeptonPt", 15. *Gaudi::Units::GeV,"Minimal momentum of the leading lepton in the quad "};

◆ m_lowSFOS_Cut

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_lowSFOS_Cut {this, "LowSFOSMass", 3.5 * Gaudi::Units::GeV, "Minimal mass for the lower SFOS pair"}
private

Definition at line 77 of file FourLeptonVertexingAlgorithm.h.

77{this, "LowSFOSMass", 3.5 * Gaudi::Units::GeV, "Minimal mass for the lower SFOS pair"};

◆ m_minElecPt

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_minElecPt {this, "MinElecPt", 4.5 * Gaudi::Units::GeV, " Minimum pt cut applied on the electron"}
private

Definition at line 63 of file FourLeptonVertexingAlgorithm.h.

63{this, "MinElecPt", 4.5 * Gaudi::Units::GeV, " Minimum pt cut applied on the electron"};

◆ m_minMuonPt

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_minMuonPt
private
Initial value:
{this, "MinMuonPt", 3. * Gaudi::Units::GeV,
"Minimum pt cut applied on the muon -- Unfortunately not handled by the MST"}

Definition at line 61 of file FourLeptonVertexingAlgorithm.h.

61 {this, "MinMuonPt", 3. * Gaudi::Units::GeV,
62 "Minimum pt cut applied on the muon -- Unfortunately not handled by the MST"};

◆ m_muonKey

SG::ReadHandleKey<xAOD::MuonContainer> DerivationFramework::FourLeptonVertexingAlgorithm::m_muonKey
private
Initial value:
{this, "MuonContainer", "Muons",
"Location of the input muon container in the StoreGate"}

Input containers.

Definition at line 48 of file FourLeptonVertexingAlgorithm.h.

48 {this, "MuonContainer", "Muons",
49 "Location of the input muon container in the StoreGate"};

◆ m_muonSelTool

ToolHandle<CP::IMuonSelectionTool> DerivationFramework::FourLeptonVertexingAlgorithm::m_muonSelTool
private
Initial value:
{this, "MuonSelectionTool", "",
"Applied preselection on the muons. Will be ignored if empty"}

Definition at line 54 of file FourLeptonVertexingAlgorithm.h.

54 {this, "MuonSelectionTool", "",
55 "Applied preselection on the muons. Will be ignored if empty"};

◆ m_muonTrk

MuonTrk DerivationFramework::FourLeptonVertexingAlgorithm::m_muonTrk {MuonTrk::Primary}
private

Definition at line 68 of file FourLeptonVertexingAlgorithm.h.

68{MuonTrk::Primary};

◆ m_muonTrkProp

Gaudi::Property<int> DerivationFramework::FourLeptonVertexingAlgorithm::m_muonTrkProp
private
Initial value:
{this, "PickMuonTrk", MuonTrk::InnerDetectorTrackParticle,
"Pick the proper track particle from the muon"}

Definition at line 66 of file FourLeptonVertexingAlgorithm.h.

66 {this, "PickMuonTrk", MuonTrk::InnerDetectorTrackParticle,
67 "Pick the proper track particle from the muon"};

◆ m_pruneCov

Gaudi::Property<bool> DerivationFramework::FourLeptonVertexingAlgorithm::m_pruneCov {this, "PruneCovariance", true, "Clears the covariance vector"}
private

Definition at line 87 of file FourLeptonVertexingAlgorithm.h.

87{this, "PruneCovariance", true, "Clears the covariance vector"};

◆ m_pruneWeight

Gaudi::Property<bool> DerivationFramework::FourLeptonVertexingAlgorithm::m_pruneWeight {this, "PruneTrackWeights", true, "Clears the track weight vector"}
private

Definition at line 88 of file FourLeptonVertexingAlgorithm.h.

88{this, "PruneTrackWeights", true, "Clears the track weight vector"};

◆ m_SubLeadPtCut

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_SubLeadPtCut {this, "SubLeadingLeptonPt", 10. *Gaudi::Units::GeV,"Minimal momentum of the leading lepton in the quad "}
private

Definition at line 81 of file FourLeptonVertexingAlgorithm.h.

81{this, "SubLeadingLeptonPt", 10. *Gaudi::Units::GeV,"Minimal momentum of the leading lepton in the quad "};

◆ m_tried_fits

std::atomic<unsigned int> DerivationFramework::FourLeptonVertexingAlgorithm::m_tried_fits {0}
mutableprivate

How many vertex fits were attempted.

Definition at line 96 of file FourLeptonVertexingAlgorithm.h.

96{0};

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.

◆ m_VtxChi2Cut

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_VtxChi2Cut {this, "VertexChi2Cut", 20, "Maximal chii n.D.o.F cut on the reconstructed vertex" }
private

Definition at line 83 of file FourLeptonVertexingAlgorithm.h.

83{this, "VertexChi2Cut", 20, "Maximal chii n.D.o.F cut on the reconstructed vertex" };

◆ m_vtxKey

SG::WriteHandleKey<xAOD::VertexContainer> DerivationFramework::FourLeptonVertexingAlgorithm::m_vtxKey {this, "OutVertexContainer", "FourLeptonVertices", "Output vertex container"}
private

Definition at line 86 of file FourLeptonVertexingAlgorithm.h.

86{this, "OutVertexContainer", "FourLeptonVertices", "Output vertex container"};

◆ m_z0Cut

Gaudi::Property<float> DerivationFramework::FourLeptonVertexingAlgorithm::m_z0Cut
private
Initial value:
{this, "DeltaZ0Cut", 7.5 * Gaudi::Units::mm,
"All leptons in the quadruplet have to be seperated from each other by maximum"}

Definition at line 73 of file FourLeptonVertexingAlgorithm.h.

73 {this, "DeltaZ0Cut", 7.5 * Gaudi::Units::mm,
74 "All leptons in the quadruplet have to be seperated from each other by maximum"};

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