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

#include <TrigMultiTrkComboHypo.h>

Inheritance diagram for TrigMultiTrkComboHypo:
Collaboration diagram for TrigMultiTrkComboHypo:

Public Member Functions

 TrigMultiTrkComboHypo (const std::string &name, ISvcLocator *pSvcLocator)
 TrigMultiTrkComboHypo ()=delete
virtual StatusCode initialize () override
virtual StatusCode execute (const EventContext &context) const override
template<typename T>
StatusCode mergeLeptonsFromDecisions (TrigMultiTrkState< T > &state) const
template<typename T>
StatusCode findMultiLeptonCandidates (TrigMultiTrkState< T > &state) const
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

const SG::ReadHandleKeyArray< TrigCompositeUtils::DecisionContainer > & decisionsInput () const
const SG::WriteHandleKeyArray< TrigCompositeUtils::DecisionContainer > & decisionsOutput () const
const Combo::MultiplicityReqMap & triggerMultiplicityMap () const
const Combo::LegMap & legToInputCollectionMap () const
ToolHandleArray< ComboHypoToolBase > & hypoTools ()
const ToolHandleArray< ComboHypoToolBase > & hypoTools () const
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

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

StatusCode mergeTracksFromViews (TrigMultiTrkStateBase &) const
 Go through state.previousDecisions() and fetch xAOD::TrackParticle objects associated with the nearest SG::View.
template<typename CONTAINER>
StatusCode mergeTracksFromDecisions (TrigMultiTrkStateBase &) const
 Go through state.previousDecisions(), fetch xAOD::Muons/xAODElectron objects attached to decisions and save links to the their xAOD::TrackParticle objects in state.tracks().
template<typename CONTAINER>
StatusCode mergeLeptonsFromDecisions (TrigMultiTrkState< CONTAINER > &) const
 Go through state.previousDecisions(), fetch xAOD::Muons/xAODElectron objects attached to decisions and save links to them in state.leptons().
template<typename CONTAINER>
StatusCode findMultiLeptonCandidates (TrigMultiTrkState< CONTAINER > &) const
 Make all possible combinations from state.leptons(), fit tracks to the common vertex, create xAOD::TrigBphys objects and put them into state.trigBphysCollection()
StatusCode filterTrackCombinations (TrigMultiTrkStateBase &) const
 Make all possible combinations from state.tracks(), fit tracks to the common vertex and set state.isEventAccepted to true if at least one good combination is found; no xAOD::TrigBphys objects will be created.
StatusCode processMergedElectrons (TrigMultiTrkState< xAOD::ElectronContainer > &) const
 Make all possible combinations from electrons originating from the same BPH-0DR3-EM7J15 cluster, only one track should pass e5_lhvloose requirements; initialRoI word will be saved into xAOD::TrigBphys objects.
StatusCode findMuTrkCandidates (TrigMultiTrkState< xAOD::MuonContainer > &) const
 Build J/psi candidates from muon from SG::View and tracks from the same SG::View, to be used for Tag-and-Probe PEB chains (similar to HLT_mu6_bJpsimutrk_MuonTrkPEB_L1MU5VF)
StatusCode copyDecisionObjects (TrigMultiTrkStateBase &) const
 All appropriate decisions from state.previousDecisions() will be copied to state.decisions() if flag state.isEventAccepted() is set by filterTrackCombinations() method; to be used only in Streamer mode.
StatusCode copyAdditionalDecisionObjects (TrigMultiTrkStateBase &) const
 For chains from CombinedSlice (similar to 'HLT_e9_lhvloose_e5_lhvloose_bBeeM6000_mu4_L1BPH-0M9-EM7-EM5_MU6') we have decisionsInput().size() > 1, so that we should copy decisions created by EmptySteps from muon part.
StatusCode createDecisionObjects (TrigMultiTrkStateBase &) const
 Create a decision for each xAOD::TrigBphys object from state.trigBphysCollection() and save it to state.decisions(); use hypoTools() to assign correct decisionIDs, not compatible with Streamer mode.
std::unique_ptr< xAOD::Vertexfit (const std::vector< ElementLink< xAOD::TrackParticleContainer > > &trackParticleLinks, const std::vector< double > &particleMasses, Trk::IVKalState &fitterState) const
 Perform a vertex fit on selected tracks.
xAOD::TrigBphysmakeTrigBPhys (const xAOD::Vertex &vertex, const std::vector< double > &particleMasses, const xAOD::Vertex &beamSpot, const Trk::IVKalState &fitterState) const
 Construct the trigger object that may be stored for debugging or matching.
bool isIdenticalTracks (const xAOD::TrackParticle *lhs, const xAOD::TrackParticle *rhs) const
 Attempts to identify identical tracks by selection on DeltaR.
bool isIdenticalTracks (const xAOD::Muon *lhs, const xAOD::Muon *rhs) const
bool isIdenticalTracks (const xAOD::Electron *lhs, const xAOD::Electron *rhs) const
bool isInMassRange (double mass, size_t idx) const
bool passedDeltaRcut (const std::vector< xAOD::TrackParticle::GenVecFourMom_t > &momenta) const
StatusCode copyDecisions (const Combo::LegDecisionsMap &passingLegs, const EventContext &context) const
 iterates over the inputs and for every object (no filtering) crates output object linked to input moving the decisions that are mentioned in the passing set
StatusCode extractFeatureAndRoI (const HLT::Identifier &chainLegId, const ElementLink< TrigCompositeUtils::DecisionContainer > &EL, SG::sgkey_t &featureKey, TrigCompositeUtils::Decision::index_type &featureIndex, SG::sgkey_t &roiKey, TrigCompositeUtils::Decision::index_type &roiIndex, bool &roiFullscan, bool &objectRequestsNoMultiplicityCheck, SG::SGKeyMap< std::set< uint32_t > > &priorFeaturesMap) const
 For a given Decision node from a HypoAlg, extracts type-less identification data on the node's Feature and seeding ROI.
StatusCode fillDecisionsMap (Combo::LegDecisionsMap &dmap, const EventContext &context) const
 iterates over all inputs, associating inputs to legs
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::TrackParticleContainerm_trackParticleContainerKey {this, "TrackCollectionKey", "Tracks", "input TrackParticle container name"}
SG::WriteHandleKey< xAOD::TrigBphysContainerm_trigBphysContainerKey {this, "TrigBphysCollectionKey", "TrigBphysContainer", "output TrigBphysContainer name"}
SG::ReadCondHandleKey< InDet::BeamSpotDatam_beamSpotKey {this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot"}
Gaudi::Property< std::vector< unsigned int > > m_nTrk
Gaudi::Property< std::vector< int > > m_nTrkCharge
Gaudi::Property< std::vector< std::vector< double > > > m_trkMass
Gaudi::Property< std::vector< std::vector< double > > > m_trkPt
Gaudi::Property< std::vector< std::pair< double, double > > > m_massRange
Gaudi::Property< bool > m_applyOverlapRemoval
Gaudi::Property< bool > m_combineInputDecisionCollections
Gaudi::Property< bool > m_useLeptonMomentum
Gaudi::Property< bool > m_checkMultiplicity
Gaudi::Property< float > m_deltaR
Gaudi::Property< float > m_deltaRMax
Gaudi::Property< float > m_deltaRMin
Gaudi::Property< float > m_chi2
Gaudi::Property< bool > m_isStreamer
Gaudi::Property< bool > m_isMuTrkMode
Gaudi::Property< bool > m_doElectrons
Gaudi::Property< std::string > m_trigLevel
Gaudi::Property< std::vector< std::string > > m_mergedElectronChains
Gaudi::Property< double > m_caloClusterEtThreshold
ToolHandle< InDet::VertexPointEstimatorm_vertexPointEstimator {this, "VertexPointEstimator", "", "tool to find starting point for the vertex fitter"}
ToolHandle< Trk::TrkVKalVrtFitterm_vertexFitter {this, "VertexFitter", "", "VKalVrtFitter tool to fit tracks into the common vertex"}
ToolHandle< Trk::V0Toolsm_v0Tools {this, "V0Tools", "", "tool to calculate Lxy/LxyError of dimuon candidate wrt beam spot"}
ToolHandle< GenericMonitoringToolm_monTool {this, "MonTool", "", "monitoring tool"}
TrigCompositeUtils::DecisionIDContainer m_allowedIDs
TrigCompositeUtils::DecisionIDContainer m_resolvedElectronIDs
TrigCompositeUtils::DecisionIDContainer m_mergedElectronIDs
double m_trkPtMin = 0.0
SG::ReadHandleKeyArray< TrigCompositeUtils::DecisionContainerm_inputs { this, "HypoInputDecisions", {}, "Input Decisions" }
SG::WriteHandleKeyArray< TrigCompositeUtils::DecisionContainerm_outputs { this, "HypoOutputDecisions", {}, "Output Decisions" }
Gaudi::Property< bool > m_requireUniqueROI
Gaudi::Property< Combo::MultiplicityReqMap > m_multiplicitiesReqMap
Gaudi::Property< Combo::LegMap > m_legToInputCollectionMap
Gaudi::Property< bool > m_checkMultiplicityMap
ToolHandleArray< ComboHypoToolBasem_hypoTools {this, "ComboHypoTools", {}, "Tools to perform selection"}
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 138 of file TrigMultiTrkComboHypo.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ TrigMultiTrkComboHypo() [1/2]

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

Definition at line 67 of file TrigMultiTrkComboHypo.cxx.

68 : ::ComboHypo(name, pSvcLocator) {}
ComboHypo(const std::string &name, ISvcLocator *pSvcLocator)
Definition ComboHypo.cxx:13

◆ TrigMultiTrkComboHypo() [2/2]

TrigMultiTrkComboHypo::TrigMultiTrkComboHypo ( )
delete

Member Function Documentation

◆ 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}

◆ copyAdditionalDecisionObjects()

StatusCode TrigMultiTrkComboHypo::copyAdditionalDecisionObjects ( TrigMultiTrkStateBase & state) const
private

For chains from CombinedSlice (similar to 'HLT_e9_lhvloose_e5_lhvloose_bBeeM6000_mu4_L1BPH-0M9-EM7-EM5_MU6') we have decisionsInput().size() > 1, so that we should copy decisions created by EmptySteps from muon part.

Definition at line 898 of file TrigMultiTrkComboHypo.cxx.

898 {
899
900 if (decisionsInput().size() > 1) {
901 ATH_MSG_DEBUG( "Found more than one decision input key, decisionsInput().size = " << decisionsInput().size() );
902 for (size_t i = 1; i < decisionsInput().size(); ++i) {
903 ATH_MSG_DEBUG( "Copying decisions from " << decisionsInput().at(i).key() << " to " << decisionsOutput().at(i).key() );
904 auto previousDecisionsHandle = SG::makeHandle(decisionsInput().at(i), state.context());
905 CHECK( previousDecisionsHandle.isValid() );
906 ATH_MSG_DEBUG( "Running with "<< previousDecisionsHandle->size() << " previous decisions" );
907 SG::WriteHandle<DecisionContainer> outputDecisionsHandle = TrigCompositeUtils::createAndStore(decisionsOutput().at(i), state.context());
908 for (const Decision* previousDecision : *previousDecisionsHandle) {
909 if (!TrigCompositeUtils::isAnyIDPassing(previousDecision, m_allowedIDs)) continue;
910
911 DecisionIDContainer previousDecisionIDs;
912 TrigCompositeUtils::decisionIDs(previousDecision, previousDecisionIDs);
914 std::set_intersection(previousDecisionIDs.begin(), previousDecisionIDs.end(), m_allowedIDs.begin(), m_allowedIDs.end(),
915 std::inserter(decisionIDs, decisionIDs.begin()));
916
918 TrigCompositeUtils::linkToPrevious(decision, previousDecision, state.context());
920 }
921 }
922 }
923 return StatusCode::SUCCESS;
924}
#define ATH_MSG_DEBUG(x)
#define CHECK(...)
Evaluate an expression and check for errors.
void decisionIDs(const Decision *d, DecisionIDContainer &id)
Extracts DecisionIDs stored in the Decision object.
const SG::WriteHandleKeyArray< TrigCompositeUtils::DecisionContainer > & decisionsOutput() const
Definition ComboHypo.h:42
const SG::ReadHandleKeyArray< TrigCompositeUtils::DecisionContainer > & decisionsInput() const
Definition ComboHypo.h:41
const EventContext & context() const
pointer_type ptr()
Dereference the pointer.
TrigCompositeUtils::DecisionIDContainer m_allowedIDs
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
void insertDecisionIDs(const Decision *src, Decision *dest)
Appends the decision IDs of src to the dest decision object.
Decision * newDecisionIn(DecisionContainer *dc, const std::string &name)
Helper method to create a Decision object, place it in the container and return a pointer to it.
const std::string & comboHypoAlgNodeName()
std::set< DecisionID > DecisionIDContainer
SG::WriteHandle< DecisionContainer > createAndStore(const SG::WriteHandleKey< DecisionContainer > &key, const EventContext &ctx)
Creates and right away records the DecisionContainer with the key.
void linkToPrevious(Decision *d, const std::string &previousCollectionKey, size_t previousIndex)
Links to the previous object, location of previous 'seed' decision supplied by hand.
void decisionIDs(const Decision *d, DecisionIDContainer &destination)
Extracts DecisionIDs stored in the Decision object.
bool isAnyIDPassing(const Decision *d, const DecisionIDContainer &required)
Checks if any of the DecisionIDs passed in arg required is availble in Decision object.

◆ copyDecisionObjects()

StatusCode TrigMultiTrkComboHypo::copyDecisionObjects ( TrigMultiTrkStateBase & state) const
private

All appropriate decisions from state.previousDecisions() will be copied to state.decisions() if flag state.isEventAccepted() is set by filterTrackCombinations() method; to be used only in Streamer mode.

Definition at line 834 of file TrigMultiTrkComboHypo.cxx.

834 {
835
836 if (state.isEventAccepted()) {
837 ATH_MSG_DEBUG( "Copying decisions from " << decisionsInput().at(0).key() << " to " << decisionsOutput().at(0).key() );
838 for (const Decision* previousDecision : state.previousDecisions()) {
839 if (!TrigCompositeUtils::isAnyIDPassing(previousDecision, m_allowedIDs)) continue;
840
841 DecisionIDContainer previousDecisionIDs;
842 TrigCompositeUtils::decisionIDs(previousDecision, previousDecisionIDs);
844 std::set_intersection(previousDecisionIDs.begin(), previousDecisionIDs.end(), m_allowedIDs.begin(), m_allowedIDs.end(),
845 std::inserter(decisionIDs, decisionIDs.begin()));
846
848 TrigCompositeUtils::linkToPrevious(decision, previousDecision, state.context());
850 }
851
852 // copy additional decisions for combined chains, as 'HLT_e9_lhvloose_e5_lhvloose_bBeeM6000_mu4_L1BPH-0M9-EM7-EM5_MU6'
854 }
855 return StatusCode::SUCCESS;
856}
#define ATH_CHECK
Evaluate an expression and check for errors.
TrigCompositeUtils::DecisionContainer & decisions()
const TrigCompositeUtils::DecisionContainer & previousDecisions() const
Gaudi::Property< bool > m_combineInputDecisionCollections
StatusCode copyAdditionalDecisionObjects(TrigMultiTrkStateBase &) const
For chains from CombinedSlice (similar to 'HLT_e9_lhvloose_e5_lhvloose_bBeeM6000_mu4_L1BPH-0M9-EM7-EM...

◆ copyDecisions()

StatusCode ComboHypo::copyDecisions ( const Combo::LegDecisionsMap & passingLegs,
const EventContext & context ) const
privateinherited

iterates over the inputs and for every object (no filtering) crates output object linked to input moving the decisions that are mentioned in the passing set

Definition at line 86 of file ComboHypo.cxx.

86 {
87 DecisionIDContainer passing;
88 for (auto const& element : passingLegs) {
89 passing.insert(element.first);
90 }
91
92 ATH_MSG_DEBUG( "Copying "<<passing.size()<<" positive decision IDs to outputs");
93
94 for ( size_t input_counter = 0; input_counter < m_inputs.size(); ++input_counter ) {
95 // new output decisions
96 SG::WriteHandle<DecisionContainer> outputHandle = createAndStore(m_outputs.at(input_counter), context );
97 auto outDecisions = outputHandle.ptr();
98 auto inputHandle = SG::makeHandle( m_inputs.at(input_counter), context );
99 if ( inputHandle.isValid() ) {
100
101 for (const Decision* inputDecision : *inputHandle) {
102 auto thisEL = TrigCompositeUtils::decisionToElementLink(inputDecision, context);
103
104 // from all positive decision in the input only the ones that survived counting are passed over
105 const DecisionIDContainer& common = passedDecisionIDs(inputDecision, passing);
106
107 // check if this EL is in the combination map for the passing decIDs:
108 ATH_MSG_DEBUG("Searching this element in the map: ("<<thisEL.dataID() << " , " << thisEL.index()<<")");
109 DecisionIDContainer finalIds;
110 for (const DecisionID c : common){
111 const HLT::Identifier cID = HLT::Identifier(c);
112 // add the decID only if this candidate passed the combination selection
113 const std::vector<ElementLink<DecisionContainer>>& Comb=passingLegs.at(c);
114 if(std::find(Comb.begin(), Comb.end(), thisEL) == Comb.end()) {
115 continue;
116 }
117 ATH_MSG_DEBUG(" Adding "<< cID <<" because EL is found in the passingLegs map");
118 finalIds.insert( cID.numeric() ); // all Ids used by the Filter, including legs
119 if (TrigCompositeUtils::isLegId ( cID )){
120 const HLT::Identifier mainChain = TrigCompositeUtils::getIDFromLeg( cID );
121 finalIds.insert( mainChain.numeric() );
122 ATH_MSG_DEBUG(" Adding "<< mainChain <<" consequently");
123 }
124 }
125
126 Decision* newDec = newDecisionIn( outDecisions, inputDecision, comboHypoAlgNodeName(), context );
127 ATH_MSG_DEBUG("New decision (Container Index:" << input_counter << ", Element Index:"<< newDec->index() <<") has "
129 << " valid initialRoI, "<< TrigCompositeUtils::getLinkToPrevious(newDec).size() <<" previous decisions and "<<finalIds.size()<<" decision IDs") ;
130 insertDecisionIDs( finalIds, newDec );
131
132 }
133 }
134
135 if (msgLvl(MSG::DEBUG)) {
136 ATH_MSG_DEBUG("Output Handle " << m_outputs.at(input_counter).key() << " with " << outputHandle->size() <<" Decision objects");
137 for (const Decision* d : *outputHandle){
138 DecisionIDContainer objDecisions;
139 decisionIDs( d, objDecisions );
140 ATH_MSG_DEBUG(" Decision object #" << d->index() << " with " << objDecisions.size()<<" positive decision IDs");
141 for (const TrigCompositeUtils::DecisionID id : objDecisions ) {
142 ATH_MSG_DEBUG(" --- Passes: " << HLT::Identifier( id ));
143 }
144 }
145 }
146 }
147
148 return StatusCode::SUCCESS;
149}
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition AtlasPID.h:878
bool msgLvl(const MSG::Level lvl) const
SG::WriteHandleKeyArray< TrigCompositeUtils::DecisionContainer > m_outputs
Definition ComboHypo.h:51
SG::ReadHandleKeyArray< TrigCompositeUtils::DecisionContainer > m_inputs
Definition ComboHypo.h:50
TrigCompositeUtils::DecisionID numeric() const
numeric ID
size_t index() const
Return the index of this element within its container.
unsigned int DecisionID
const std::vector< ElementLink< DecisionContainer > > getLinkToPrevious(const Decision *d)
returns links to previous decision object 'seed'
DecisionIDContainer passedDecisionIDs(const Decision *d, const T &required)
return DecisionIDs in Decision object that match the required ones
HLT::Identifier getIDFromLeg(const HLT::Identifier &legIdentifier)
Generate the HLT::Identifier which corresponds to the chain name from the leg name.
LinkInfo< T > findLink(const Decision *start, const std::string &linkName, const bool suppressMultipleLinksWarning=false)
Perform a recursive search for ElementLinks of type T and name 'linkName', starting from Decision obj...
const std::string & initialRoIString()
ElementLink< DecisionContainer > decisionToElementLink(const Decision *d, const EventContext &ctx)
Takes a raw pointer to a Decision and returns an ElementLink to the Decision.
bool isLegId(const HLT::Identifier &legIdentifier)
Recognise whether the chain ID is a leg ID.

◆ createDecisionObjects()

StatusCode TrigMultiTrkComboHypo::createDecisionObjects ( TrigMultiTrkStateBase & state) const
private

Create a decision for each xAOD::TrigBphys object from state.trigBphysCollection() and save it to state.decisions(); use hypoTools() to assign correct decisionIDs, not compatible with Streamer mode.

Definition at line 859 of file TrigMultiTrkComboHypo.cxx.

859 {
860
861 size_t idx = 0;
862 for (const xAOD::TrigBphys* triggerObject : state.trigBphysCollection()) {
863 ATH_MSG_DEBUG( "Found xAOD::TrigBphys: mass / chi2 = " << triggerObject->mass() << " / " << triggerObject->fitchi2() );
864
865 auto triggerObjectEL = ElementLink<xAOD::TrigBphysContainer>(state.trigBphysCollection(), triggerObject->index());
866 ATH_CHECK( triggerObjectEL.isValid() );
867
868 // create a new output Decision object, backed by the 'decisions' container.
870
871 std::vector<const DecisionIDContainer*> previousDecisionIDs;
872 for (const size_t& i : state.getTrigBphysLegIndices(idx)) {
873
874 // attach all previous decisions: if the same previous decision is called twice, that's fine - internally takes care of that
875 // we already have an array of links to the previous decisions, so there is no need to use TrigCompositeUtils::linkToPrevious()
877 previousDecisionIDs.push_back(&state.getDecisionIDs(i));
878 }
879
880 // set mandatory feature ElementLink to xAOD::TrigBphys object
882 decision->setDetail<int32_t>("noCombo", 1);
883
884 for (const auto& tool : hypoTools()) {
885 ATH_MSG_DEBUG( "Go to " << tool );
886 if (!m_checkMultiplicity || state.checkMultiplicity(tool->legMultiplicity(), tool->legDecisionIds())) ATH_CHECK( tool->decideOnSingleObject(decision, previousDecisionIDs) );
887 }
888 ++idx;
889 }
890
891 // copy additional decisions from Empty muon step for Combined chains, as 'HLT_e9_lhvloose_e5_lhvloose_bBeeM6000_mu4_L1BPH-0M9-EM7-EM5_MU6'
893
894 return StatusCode::SUCCESS;
895}
ToolHandleArray< ComboHypoToolBase > & hypoTools()
Definition ComboHypo.h:45
xAOD::TrigBphysContainer & trigBphysCollection()
Gaudi::Property< bool > m_checkMultiplicity
virtual bool checkMultiplicity(const std::vector< int > &legMultiplicity, const std::vector< HLT::Identifier > &legDecisionIDs) const =0
std::vector< size_t > & getTrigBphysLegIndices(size_t idx)
virtual TrigCompositeUtils::DecisionIDContainer & getDecisionIDs(size_t)=0
virtual std::vector< ElementLink< TrigCompositeUtils::DecisionContainer > > & getDecisionLinks(size_t)=0
bool setObjectLink(const std::string &name, const ElementLink< CONTAINER > &link)
Set the link to an object.
bool setDetail(const std::string &name, const TYPE &value)
Set an TYPE detail on the object.
bool addObjectCollectionLinks(const std::string &collectionName, const std::vector< ElementLink< CONTAINER > > &links)
Add links to multiple objects within a collection. Performs de-duplication.
const std::string & featureString()
const std::string & seedString()
TrigBphys_v1 TrigBphys
Definition TrigBphys.h:18
TrigBphysContainer_v1 TrigBphysContainer

◆ decisionsInput()

const SG::ReadHandleKeyArray< TrigCompositeUtils::DecisionContainer > & ComboHypo::decisionsInput ( ) const
inlineprotectedinherited

Definition at line 41 of file ComboHypo.h.

41{ return m_inputs; }

◆ decisionsOutput()

const SG::WriteHandleKeyArray< TrigCompositeUtils::DecisionContainer > & ComboHypo::decisionsOutput ( ) const
inlineprotectedinherited

Definition at line 42 of file ComboHypo.h.

42{ return m_outputs; }

◆ 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 TrigMultiTrkComboHypo::execute ( const EventContext & context) const
overridevirtual

Reimplemented from ComboHypo.

Definition at line 210 of file TrigMultiTrkComboHypo.cxx.

210 {
211
212 ATH_MSG_DEBUG( "TrigMultiTrkHypo::execute() starts" );
213
214 ATH_MSG_DEBUG( "decision input key: " << decisionsInput().at(0).key() );
215 auto previousDecisionsHandle = SG::makeHandle(decisionsInput().at(0), context);
216 CHECK( previousDecisionsHandle.isValid() );
217 ATH_MSG_DEBUG( "Running with "<< previousDecisionsHandle->size() << " previous decisions" );
218
219 SG::WriteHandle<DecisionContainer> outputDecisionsHandle = TrigCompositeUtils::createAndStore(decisionsOutput().at(0), context);
220
221 const InDet::BeamSpotData* beamSpotData = nullptr;
222 if (!m_isStreamer) {
223 SG::ReadCondHandle<InDet::BeamSpotData> beamSpotHandle {m_beamSpotKey, context};
224 ATH_CHECK( beamSpotHandle.isValid() );
225 beamSpotData = *beamSpotHandle;
226 }
227
228 std::unique_ptr<TrigMultiTrkState<xAOD::MuonContainer>> muonState;
229 std::unique_ptr<TrigMultiTrkState<xAOD::ElectronContainer>> electronState;
230 TrigMultiTrkStateBase* commonState = nullptr;
231 if (m_doElectrons) {
232 electronState = std::make_unique<TrigMultiTrkState<xAOD::ElectronContainer>>(context, *previousDecisionsHandle, *outputDecisionsHandle, nullptr, beamSpotData);
233 commonState = electronState.get();
234 }
235 else {
236 muonState = std::make_unique<TrigMultiTrkState<xAOD::MuonContainer>>(context, *previousDecisionsHandle, *outputDecisionsHandle, nullptr, beamSpotData);
237 commonState = muonState.get();
238 }
239
240 if (m_isStreamer) {
241 if (m_trigLevel == "L2") {
242 ATH_CHECK( mergeTracksFromViews(*commonState) );
243 }
244 else if (m_trigLevel == "L2IO" || m_trigLevel == "L2MT") {
246 }
247 else if (m_trigLevel == "EF") {
249 }
250 ATH_CHECK( filterTrackCombinations(*commonState) );
251 ATH_CHECK( copyDecisionObjects(*commonState) );
252 }
253 else {
254 auto trigBphysHandle = SG::makeHandle(m_trigBphysContainerKey, context);
255 ATH_CHECK( trigBphysHandle.record(std::make_unique<xAOD::TrigBphysContainer>(),
256 std::make_unique<xAOD::TrigBphysAuxContainer>()) );
257 commonState->setTrigBphysCollection(trigBphysHandle.ptr());
258
259 if (m_doElectrons) {
260 ATH_CHECK( mergeLeptonsFromDecisions(*electronState) );
261 ATH_CHECK( findMultiLeptonCandidates(*electronState) );
262 ATH_CHECK( processMergedElectrons(*electronState) );
263 }
264 else if (m_isMuTrkMode) {
265 ATH_CHECK( findMuTrkCandidates(*muonState) );
266 }
267 else {
270 }
271 ATH_CHECK( createDecisionObjects(*commonState) );
272 }
273
274 ATH_MSG_DEBUG( "TrigMultiTrkHypo::execute() terminates with StatusCode::SUCCESS" );
275 return StatusCode::SUCCESS;
276}
StatusCode mergeLeptonsFromDecisions(TrigMultiTrkState< CONTAINER > &) const
Go through state.previousDecisions(), fetch xAOD::Muons/xAODElectron objects attached to decisions an...
SG::WriteHandleKey< xAOD::TrigBphysContainer > m_trigBphysContainerKey
StatusCode findMuTrkCandidates(TrigMultiTrkState< xAOD::MuonContainer > &) const
Build J/psi candidates from muon from SG::View and tracks from the same SG::View, to be used for Tag-...
StatusCode filterTrackCombinations(TrigMultiTrkStateBase &) const
Make all possible combinations from state.tracks(), fit tracks to the common vertex and set state....
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
Gaudi::Property< bool > m_isStreamer
Gaudi::Property< bool > m_doElectrons
StatusCode findMultiLeptonCandidates(TrigMultiTrkState< CONTAINER > &) const
Make all possible combinations from state.leptons(), fit tracks to the common vertex,...
StatusCode copyDecisionObjects(TrigMultiTrkStateBase &) const
All appropriate decisions from state.previousDecisions() will be copied to state.decisions() if flag ...
StatusCode processMergedElectrons(TrigMultiTrkState< xAOD::ElectronContainer > &) const
Make all possible combinations from electrons originating from the same BPH-0DR3-EM7J15 cluster,...
Gaudi::Property< std::string > m_trigLevel
Gaudi::Property< bool > m_isMuTrkMode
StatusCode mergeTracksFromViews(TrigMultiTrkStateBase &) const
Go through state.previousDecisions() and fetch xAOD::TrackParticle objects associated with the neares...
StatusCode createDecisionObjects(TrigMultiTrkStateBase &) const
Create a decision for each xAOD::TrigBphys object from state.trigBphysCollection() and save it to sta...
StatusCode mergeTracksFromDecisions(TrigMultiTrkStateBase &) const
Go through state.previousDecisions(), fetch xAOD::Muons/xAODElectron objects attached to decisions an...

◆ extractFeatureAndRoI()

StatusCode ComboHypo::extractFeatureAndRoI ( const HLT::Identifier & chainLegId,
const ElementLink< TrigCompositeUtils::DecisionContainer > & EL,
SG::sgkey_t & featureKey,
TrigCompositeUtils::Decision::index_type & featureIndex,
SG::sgkey_t & roiKey,
TrigCompositeUtils::Decision::index_type & roiIndex,
bool & roiFullscan,
bool & objectRequestsNoMultiplicityCheck,
SG::SGKeyMap< std::set< uint32_t > > & priorFeaturesMap ) const
privateinherited

For a given Decision node from a HypoAlg, extracts type-less identification data on the node's Feature and seeding ROI.

Parameters
[in]chainLegIdThe HLT::Identifer of the chain (leg) we're extracting features for.
[in]ELThe Decision node from the HypoAlg, expected to have a "feature" link attached to it. Expected to be able to locate a "initialRoI" in its history if RequireUniqueROI=True.
[out]featureKeyType-less SG Key hash of the collection hosting the Decision node's feature .
[out]featureIndexIndex inside the featureKey collection.
[out]roiKeyType-less SG Key hash of the collection hosting the Decision node's initial ROI collection.
[out]roiIndexIndex inside the roiKey collection.
[out]roiIsFullscanFlag indicating if the located initial ROI has the FullScan flag enabled. Triggers special behaviour allowing the ROI to satisfy arbitrary multiplicities in an arbitrary number of legs.
[out]objectRequestsNoMultiplicityCheckFlag indicating of the DecisionObject requested not be included in the multiplicity computation. Triggers special behaviour allowing the DecisionObject to satisfy arbitrary multiplicities in an arbitrary number of legs.
[in,out]priorFeaturesMapData structure collating for a given feature (key) what the prior features were integrated over all previous steps (value set).

Definition at line 378 of file ComboHypo.cxx.

387{
388 // Return collections for the findLinks call.
389 // While we will be focusing on the most recent feature, for tag-and-probe we need to keep a record of the features from the prior steps too.
390
391 // Construct a sub-graph following just this leg back through the nav
392 DecisionIDContainer chainLegIdSet = {chainLegId.numeric()};
393 TrigCompositeUtils::NavGraph subGraph;
394 recursiveGetDecisions((*dEL), subGraph, chainLegIdSet, /*enforceDecisionOnStartNode =*/ true);
395
396 if (subGraph.finalNodes().size() != 1) {
397 ATH_MSG_ERROR("We are only expecting to search from a single navigation node in extractFeatureAndRoI");
398 return StatusCode::FAILURE;
399 }
400 const NavGraphNode* start = *(subGraph.finalNodes().begin());
401
402 std::vector<SG::sgkey_t> keys;
403 std::vector<uint32_t> clids; // We don't care about the class ID. This part gets ignored.
404 std::vector<Decision::index_type> indicies;
405 std::vector<const Decision*> sources;
406
407 std::set<const Decision*> fullyExploredFrom; // This is a cache which typelessFindLinks will use to avoid re-visiting already explored regions of the graph
408 // Note: This call to typelessFindLinks is exploring from a NavGraphNode* rather than a Decision*,
409 // this indicates that the search is restricted to a sub-graph (specifically, only following one chain-leg)
410 const bool foundFeature = typelessFindLinks(start, featureString(), keys, clids, indicies, sources, TrigDefs::allFeaturesOfType, &fullyExploredFrom);
411
412 const Decision* featureSource = nullptr;
413 // The "most recent" feature (from the step just run) is the one we find first. Hence it's at index 0
414 if (foundFeature) {
415 featureKey = keys.at(0);
416 featureIndex = indicies.at(0);
417 featureSource = sources.at(0);
418 }
419
420 objectRequestsNoMultiplicityCheck = (featureSource and featureSource->hasDetail<int32_t>("noCombo") and featureSource->getDetail<int32_t>("noCombo") == 1);
421
422 if (foundFeature and priorFeaturesMap.count(featureKey + featureIndex) == 0) {
423 const std::string* key_str = evtStore()->keyToString(featureKey);
424 ATH_MSG_DEBUG("Note: Will use feature hash " << featureKey + featureIndex << ", for " << (key_str ? *key_str : "UNKNOWN") << " index=" << featureIndex);
425 // Use the deep-search data to look further back than .at(0)
426 // Here's where we keep the record of the features in previous steps. Step ordering is unimportant, we can use a set.
427 if (keys.size() > 1) {
428 for (size_t i = 1; i < keys.size(); ++i) { // Skip the 1st entry, this will be equal to featureKey and featureIndex from typelessFindLink above.
429 // featureKey + featureIndex should be considered as equivalent to a per-feature hash (featureKey is a real hash, featureIndex is an offset index)
430 if (featureKey + featureIndex == keys.at(i) + indicies.at(i)) {
431 continue; // Do not add the case where a feature is re-attached to more than one step.
432 }
433 priorFeaturesMap[featureKey + featureIndex].insert(keys.at(i) + indicies.at(i));
434 }
435 } else { // Exactly one feature. Make a note of this by inserting an empty set, such that we don't do this search again.
436 priorFeaturesMap.insert( std::pair<uint32_t, std::set<uint32_t>>(featureKey + featureIndex, std::set<uint32_t>()) );
437 }
438 }
439
440 // Try and get seeding ROI data too.
441 uint32_t roiClid{0}; // Unused
442 const Decision* roiSource{nullptr}; // Unused
443 const bool foundROI = typelessFindLink(subGraph, initialRoIString(), roiKey, roiClid, roiIndex, roiSource);
444 if (foundROI) {
445 ElementLink<TrigRoiDescriptorCollection> roiEL(roiKey, roiIndex);
446 ATH_CHECK( roiEL.isValid() );
447 roiIsFullscan = (*(roiEL))->isFullscan();
448 if (!foundFeature) {
449 const std::string* roi_str = evtStore()->keyToString(roiKey);
450 ATH_MSG_DEBUG("Note: Located fallback-ROI, if used this will have feature hash =" << roiKey + roiIndex << ", for " << (roi_str ? *roi_str : "UNKNOWN") << " index=" << roiIndex);
451 }
452 }
453
454 if (!foundFeature && !foundROI) {
455 ATH_MSG_WARNING("Did not find the feature or initialRoI for " << dEL.dataID() << " index " << dEL.index());
456 }
457
458 return StatusCode::SUCCESS;
459}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
const std::vector< NavGraphNode * > & finalNodes() const
Get all final nodes.
Definition NavGraph.cxx:99
bool hasDetail(const std::string &name) const
Check if a given type of detail is available.
bool getDetail(const std::string &name, TYPE &value) const
Get an TYPE detail from the object.
void recursiveGetDecisions(const Decision *start, NavGraph &navGraph, const DecisionIDContainer &ids, const bool enforceDecisionOnStartNode)
Search back in time from "node" and locate all paths back through Decision objects for a given chain.
bool typelessFindLinks(const Decision *start, const std::string &linkName, std::vector< sgkey_t > &keyVec, std::vector< CLID > &clidVec, std::vector< Decision::index_type > &indexVec, std::vector< const Decision * > &sourceVec, const unsigned int behaviour, std::set< const Decision * > *fullyExploredFrom)
search back the TC links for the object of type T linked to the one of TC (recursively) Returns the l...
bool typelessFindLink(const Decision *start, const std::string &linkName, sgkey_t &key, CLID &clid, Decision::index_type &index, const Decision *&source, const bool suppressMultipleLinksWarning)
Perform a recursive search for ElementLinks of any time and name 'linkName', starting from Decision o...
static const unsigned int allFeaturesOfType
Run 3 "enum". Return all features along legs (still with type and container checks)
setEventNumber uint32_t

◆ 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.

◆ fillDecisionsMap()

StatusCode ComboHypo::fillDecisionsMap ( Combo::LegDecisionsMap & dmap,
const EventContext & context ) const
privateinherited

iterates over all inputs, associating inputs to legs

Definition at line 462 of file ComboHypo.cxx.

462 {
463 for ( size_t inputContainerIndex = 0; inputContainerIndex < m_inputs.size(); ++inputContainerIndex ) {
464 auto inputHandle = SG::makeHandle( m_inputs.at(inputContainerIndex), context );
465 if ( !inputHandle.isValid() ) {
466 ATH_MSG_ERROR( "No input ReadHandle from " << inputHandle.key() );
467 return StatusCode::FAILURE;
468 }
469 ATH_MSG_DEBUG( "-- Found ReadHandle from " << inputHandle.key() <<" with "<< inputHandle->size() << " elements:" );
470 for ( const Decision* decision : *inputHandle ) {
471 ATH_MSG_DEBUG( "-- -- Input Decision #"<< decision->index() <<" with "<< decisionIDs( decision ).size() << " active IDs. Populating the multiplicity map:" );
472 for ( const DecisionID id: decisionIDs( decision ) ) {
473
474 // Handles regular and multi-leg case
475 const auto& [chainName, chainLeg] = getNameAndIndexFromLeg(HLT::Identifier(id).name());
476
477 // We need to check if we are configured to accept DecisionObjects passing 'chainName' ...
478 Combo::LegMap::const_iterator it = m_legToInputCollectionMap.find(chainName);
479 if (it == m_legToInputCollectionMap.end()) {
480 ATH_MSG_VERBOSE("-- -- -- Ignoring the DecsionID " << id << " on leg " << chainLeg << " as it does not correspond to any of the " << m_legToInputCollectionMap.size() << " chains this Alg is processing.");
481 continue;
482 }
483
484 // ... and if so we need to further check that we are accepting passing IDs for chainLeg on the current inputContainerIndex
485 const std::vector<int>& legToInputCollectionIndex = it->second;
486 const size_t requiredInputContainerIndex = static_cast<size_t>(legToInputCollectionIndex.at(chainLeg));
487 if (requiredInputContainerIndex != inputContainerIndex) {
488 ATH_MSG_VERBOSE("-- -- -- Ignoring the DecisionID " << id << " on leg " << chainLeg << " as we are only permitted to accept passing objects on leg #" << chainLeg << " of " << chainName
489 << " which come from input collection index " << requiredInputContainerIndex << " (which is " << m_inputs.at(requiredInputContainerIndex).key() << ")"
490 << ". Not the current index " << inputContainerIndex << " (which is " << m_inputs.at(inputContainerIndex).key() << ")");
491 continue;
492 }
493
494 ATH_MSG_DEBUG( " ++++ " << HLT::Identifier( id ) );
495 dmap[id].push_back( TrigCompositeUtils::decisionToElementLink(decision, context) );
496 }
497 }
498 }
499
500 if (msgLvl(MSG::DEBUG)){
501 ATH_MSG_DEBUG( "Decision map filled :" );
502 size_t legCount = 0;
503 for (const auto& entry: dmap){
504 ATH_MSG_DEBUG("leg ["<<legCount<<"]: ");
505 const std::vector<ElementLink<DecisionContainer>>& decisions = entry.second;
506 ATH_MSG_DEBUG(" ++++ " << HLT::Identifier( entry.first ) <<" Number Decisions: "<< decisions.size());
507 for (const ElementLink<DecisionContainer>& d : decisions){
508 ATH_MSG_DEBUG(" Decision: (ContainerKey:"<<d.dataID()<<", DecisionElementIndex:"<<d.index()<<")");
509 }
510 legCount++;
511 }
512 }
513
514
515 return StatusCode::SUCCESS;
516}
#define ATH_MSG_VERBOSE(x)
Gaudi::Property< Combo::LegMap > m_legToInputCollectionMap
Definition ComboHypo.h:59
std::pair< std::string, int32_t > getNameAndIndexFromLeg(const std::string &name)
Extract the name and numeric index of a leg identifier.

◆ 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

◆ filterTrackCombinations()

StatusCode TrigMultiTrkComboHypo::filterTrackCombinations ( TrigMultiTrkStateBase & state) const
private

Make all possible combinations from state.tracks(), fit tracks to the common vertex and set state.isEventAccepted to true if at least one good combination is found; no xAOD::TrigBphys objects will be created.

Definition at line 492 of file TrigMultiTrkComboHypo.cxx.

492 {
493
494 const auto& tracks = state.tracks();
495 state.setEventAccepted(false);
496
497 // monitored variables
498 auto mon_nAcceptedTrk = Monitored::Scalar<int>("nAcceptedTrk", tracks.size());
499 auto mon_nVertexFitterCalls = Monitored::Scalar<int>("nVertexFitterCalls", 0);
500 auto mon_isEventAccepted = Monitored::Scalar<int>("acceptance", 0);
501 auto mon_timer = Monitored::Timer( "TIME_all" );
502
503 auto group = Monitored::Group(m_monTool,
504 mon_nAcceptedTrk, mon_nVertexFitterCalls, mon_isEventAccepted,
505 mon_timer);
506
507 for (size_t iTrk = 0; iTrk < m_nTrk.size(); ++iTrk) {
508 if (state.isEventAccepted()) break;
509 size_t nTrk = m_nTrk[iTrk];
510
511 if (tracks.size() < nTrk) {
512 ATH_MSG_DEBUG( "Could not build a subset of " << nTrk << " tracks from collection which contains only " << tracks.size() << " objects" );
513 continue;
514 }
515 ATH_MSG_DEBUG( "Consider combinations of " << nTrk << " tracks from collection which contains " << tracks.size() << " objects until find a good one" );
516
517 std::vector<ElementLink<xAOD::TrackParticleContainer>> tracklist(nTrk);
518 std::vector<xAOD::TrackParticle::GenVecFourMom_t> p(nTrk);
519
520 // tracks from the current combination will have non-zero value against their position in the 'idx' vector
521 // consider first nTrk tracks as an initial combination, then get the next one with std::prev_permutation()
522 std::vector<char> idx(tracks.size(), 0);
523 std::fill(idx.begin(), idx.begin() + nTrk, 1);
524 do {
525 // fill tracklist and momenta of tracks, also check that the track pT passes the threshold value
526 bool isValidCombination = true;
527 int totalCharge = 0;
528 size_t j = 0;
529 for (size_t i = 0; i < idx.size(); ++i) {
530 if (!idx[i]) continue;
531 const auto& trackEL = tracks[i];
532 tracklist[j] = trackEL;
533 const auto track = *trackEL;
534 p[j] = track->genvecP4();
535 p[j].SetM(m_trkMass[iTrk][j]);
536 totalCharge += static_cast<int>(track->charge());
537 if (p[j].Pt() < m_trkPt[iTrk][j]) {
538 isValidCombination = false;
539 break;
540 }
541 ++j;
542 }
543 if (!isValidCombination || (m_nTrkCharge[iTrk] >= 0 && totalCharge != m_nTrkCharge[iTrk]) || !passedDeltaRcut(p)) continue;
544
545 if (msgLvl(MSG::DEBUG)) {
546 ATH_MSG_DEBUG( "Dump found tracks before vertex fit: pT / eta / phi / charge" );
547 for (size_t i = 0; i < tracklist.size(); ++i) {
548 const auto track = *tracklist[i];
549 ATH_MSG_DEBUG( "track " << i + 1 << ": " << p[i].Pt() << " / " << p[i].Eta() << " / " << p[i].Phi() << " / " << track->charge() );
550 }
551 }
552
553 auto mass = (std::accumulate(p.begin(), p.end(), xAOD::TrackParticle::GenVecFourMom_t())).M();
554 ATH_MSG_DEBUG( "invariant mass: " << mass );
555
556 if (!isInMassRange(mass, iTrk)) continue;
557
558 auto fitterState = m_vertexFitter->makeState(state.context());
559 auto vertex = fit(tracklist, m_trkMass[iTrk], *fitterState);
560 ++mon_nVertexFitterCalls;
561 if (!vertex) continue;
562
563 ATH_MSG_DEBUG( "Filter found a subset of tracks which passed the rough selection: stop looking for other combinations" );
564 state.setEventAccepted(true);
565 break;
566
567 } while (prev_perm(idx));
568 }
569
570 if (!state.isEventAccepted()) {
571 ATH_MSG_DEBUG( "Filter could not find a good subset of tracks" );
572 }
573
574 mon_isEventAccepted = state.isEventAccepted();
575
576 return StatusCode::SUCCESS;
577}
@ Phi
Definition RPCdef.h:8
@ Eta
Definition RPCdef.h:8
ToolHandle< GenericMonitoringTool > m_monTool
std::unique_ptr< xAOD::Vertex > fit(const std::vector< ElementLink< xAOD::TrackParticleContainer > > &trackParticleLinks, const std::vector< double > &particleMasses, Trk::IVKalState &fitterState) const
Perform a vertex fit on selected tracks.
bool isInMassRange(double mass, size_t idx) const
Gaudi::Property< std::vector< std::vector< double > > > m_trkMass
Gaudi::Property< std::vector< int > > m_nTrkCharge
Gaudi::Property< std::vector< std::vector< double > > > m_trkPt
Gaudi::Property< std::vector< unsigned int > > m_nTrk
bool passedDeltaRcut(const std::vector< xAOD::TrackParticle::GenVecFourMom_t > &momenta) const
ToolHandle< Trk::TrkVKalVrtFitter > m_vertexFitter
std::vector< ElementLink< xAOD::TrackParticleContainer > > & tracks()
void setEventAccepted(bool flag=true)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzM4D< double > > GenVecFourMom_t
Base 4 Momentum type for TrackParticle.
float j(const xAOD::IParticle &, const xAOD::TrackMeasurementValidation &hit, const Eigen::Matrix3d &jab_inv)

◆ findMultiLeptonCandidates() [1/2]

template<typename CONTAINER>
StatusCode TrigMultiTrkComboHypo::findMultiLeptonCandidates ( TrigMultiTrkState< CONTAINER > & ) const
private

Make all possible combinations from state.leptons(), fit tracks to the common vertex, create xAOD::TrigBphys objects and put them into state.trigBphysCollection()

combination is a selection of items from a collection, such that the order of items does not matter leptons/tracks have already been sorted over pT, the algorithm will keep this order in created subsets of tracks, i.e. for nTrk = 2 combinations are [0, 1], [0, 2], ..., [1, 2], [1, 3], ... for nTrk = 3 combinations are [0, 1, 2], [0, 1, 3], ..., [0, 2, 3], [0, 2, 4], ..., [1, 2, 3], [1, 2, 4], ...

◆ findMultiLeptonCandidates() [2/2]

template<typename T>
StatusCode TrigMultiTrkComboHypo::findMultiLeptonCandidates ( TrigMultiTrkState< T > & state) const

Definition at line 581 of file TrigMultiTrkComboHypo.cxx.

581 {
582
583 state.trigBphysLegIndices().clear();
584 const auto& leptons = state.leptons();
585
586 // monitored variables
587 auto mon_nAcceptedTrk = Monitored::Scalar<int>("nAcceptedTrk", leptons.size());
588 auto mon_nCombination = Monitored::Scalar<int>("nCombination", 0);
589 auto mon_nCombinationBeforeFit = Monitored::Scalar<int>("nCombinationBeforeFit", 0);
590 auto mon_nBPhysObject = Monitored::Scalar<int>("nBPhysObject", 0);
591
592 std::vector<float> trkMassBeforeFit;
593 std::vector<float> bphysMass;
594 std::vector<float> d0track1, d0track2;
595 std::vector<float> pttrack1, pttrack2;
596 std::vector<float> etatrack1, etatrack2;
597 std::vector<int> bphysCharge;
598 auto mon_trkMassBeforeFit = Monitored::Collection("trkMassBeforeFit", trkMassBeforeFit);
599 auto mon_bphysChi2 = Monitored::Collection("bphysChi2", state.trigBphysCollection(), &xAOD::TrigBphys::fitchi2);
600 auto mon_bphysLxy = Monitored::Collection("bphysLxy", state.trigBphysCollection(), &xAOD::TrigBphys::lxy);
601 auto mon_bphysFitMass = Monitored::Collection("bphysFitMass", state.trigBphysCollection(), [](const xAOD::TrigBphys* x){ return x->fitmass()*0.001; });
602 auto mon_bphysMass = Monitored::Collection("bphysMass", bphysMass);
603 auto mon_d0track1 = Monitored::Collection("bphysd0_trk1", d0track1);
604 auto mon_d0track2 = Monitored::Collection("bphysd0_trk2", d0track2);
605 auto mon_pttrack1 = Monitored::Collection("bphysPt_trk1", pttrack1);
606 auto mon_pttrack2 = Monitored::Collection("bphysPt_trk2", pttrack2);
607 auto mon_etatrack1 = Monitored::Collection("bphysEtatrack1", etatrack1);
608 auto mon_etatrack2 = Monitored::Collection("bphysEtatrack2", etatrack2);
609 auto mon_bphysCharge = Monitored::Collection("bphysCharge", bphysCharge);
610
611 auto mon_timer = Monitored::Timer( "TIME_all" );
612
613 auto group = Monitored::Group(m_monTool,
614 mon_nAcceptedTrk, mon_nCombination, mon_nCombinationBeforeFit, mon_nBPhysObject,
615 mon_trkMassBeforeFit, mon_bphysChi2, mon_bphysLxy, mon_bphysFitMass, mon_bphysMass, mon_bphysCharge, mon_d0track1, mon_d0track2,
616 mon_timer, mon_pttrack1, mon_pttrack2, mon_etatrack1, mon_etatrack2);
617
618 for (size_t iTrk = 0; iTrk < m_nTrk.size(); ++iTrk) {
619 size_t nTrk = m_nTrk[iTrk];
620
621 if (leptons.size() < nTrk) {
622 ATH_MSG_DEBUG( "Could not build a subset of " << nTrk << " legs from collection which contains only " << leptons.size() << " objects" );
623 continue;
624 }
625 ATH_MSG_DEBUG( "Consider all combinations of " << nTrk << " legs from collection which contains " << leptons.size() << " objects" );
626
627 std::vector<size_t> leptonIndices(nTrk);
628 std::vector<ElementLink<xAOD::TrackParticleContainer>> tracklist(nTrk);
629 std::vector<xAOD::TrackParticle::GenVecFourMom_t> p(nTrk);
630
631 std::vector<char> combination(leptons.size(), 0);
632 std::fill(combination.begin(), combination.begin() + nTrk, 1);
633 do {
634 // fill tracklist and momenta of muons in the current subset
635 bool isValidCombination = true;
636 int charge = 0;
637 size_t j = 0;
638 for (size_t i = 0; i < combination.size(); ++i) {
639 if (!combination[i]) continue;
640 leptonIndices[j] = i;
641 auto leg = *leptons[i].link;
642 charge += static_cast<int>(lround(leg->charge()));
643 ElementLink<xAOD::TrackParticleContainer> trackEL;
644 if constexpr(std::is_same<T, xAOD::MuonContainer>::value) {
645 trackEL = leg->inDetTrackParticleLink();
646 }
647 else {
648 trackEL = leg->trackParticleLink();
649 }
650 tracklist[j] = trackEL;
652 p[j] = leg->genvecP4();
653 p[j].SetM(0.); // to keep consistency with TrigComboHypoTool::compute()
654 }
655 else {
656 p[j] = (*trackEL)->genvecP4();
657 p[j].SetM(m_trkMass[iTrk][j]);
658 }
659 if (p[j].Pt() < m_trkPt[iTrk][j]) {
660 isValidCombination = false;
661 break;
662 }
663 ++j;
664 }
665 if (!isValidCombination || (m_nTrkCharge[iTrk] >= 0 && charge != m_nTrkCharge[iTrk]) || !passedDeltaRcut(p)) continue;
666
667 if (msgLvl(MSG::DEBUG)) {
668 ATH_MSG_DEBUG( "Dump found leptons before vertex fit: pT / eta / phi / charge" );
669 for (size_t i = 0; i < tracklist.size(); ++i) {
670 const auto track = *tracklist[i];
671 ATH_MSG_DEBUG( "legs " << i + 1 << ": " << p[i].Pt() << " / " << p[i].Eta() << " / " << p[i].Phi() << " / " << track->charge() );
672 }
673 }
674
675 auto mass = (std::accumulate(p.begin(), p.end(), xAOD::TrackParticle::GenVecFourMom_t())).M();
676 ATH_MSG_DEBUG( "invariant mass: " << mass );
677
678 mon_nCombination++;
679 trkMassBeforeFit.push_back(mass * 0.001);
680 if (!isInMassRange(mass, iTrk)) continue;
681
682 mon_nCombinationBeforeFit++;
683 auto fitterState = m_vertexFitter->makeState(state.context());
684 auto vertex = fit(tracklist, m_trkMass[iTrk], *fitterState);
685 if (!vertex) continue;
686 xAOD::TrigBphys* trigBphys = makeTrigBPhys(*vertex, m_trkMass[iTrk], state.beamSpot(), *fitterState);
687 if (m_useLeptonMomentum) trigBphys->setMass(mass);
688 state.addTrigBphysObject(trigBphys, leptonIndices);
689
690 mon_nBPhysObject++;
691 bphysMass.push_back(mass * 0.001);
692 bphysCharge.push_back(charge);
693 d0track1.push_back((*tracklist[0])->d0());
694 d0track2.push_back((*tracklist[1])->d0());
695 pttrack1.push_back((*tracklist[0])->pt() * 0.001);
696 pttrack2.push_back((*tracklist[1])->pt() * 0.001);
697 etatrack1.push_back((*tracklist[0])->eta());
698 etatrack2.push_back((*tracklist[1])->eta());
699
700 } while (prev_perm(combination));
701 }
702 return StatusCode::SUCCESS;
703}
Scalar eta() const
pseudorapidity method
double charge(const T &p)
Definition AtlasPID.h:997
#define x
const xAOD::Vertex & beamSpot() const
xAOD::TrigBphys * makeTrigBPhys(const xAOD::Vertex &vertex, const std::vector< double > &particleMasses, const xAOD::Vertex &beamSpot, const Trk::IVKalState &fitterState) const
Construct the trigger object that may be stored for debugging or matching.
Gaudi::Property< bool > m_useLeptonMomentum
std::vector< std::vector< size_t > > & trigBphysLegIndices()
std::vector< LEPTON > & leptons()
virtual void addTrigBphysObject(xAOD::TrigBphys *trigBphysObject, const std::vector< size_t > &legIndices) override final
float fitchi2() const
accessor method: chi2 from vertex fit
float lxy() const
accessor method: lxy
void setMass(float)
Set the mass of the object.
ValuesCollection< T > Collection(std::string name, const T &collection)
Declare a monitored (double-convertible) collection.

◆ findMuTrkCandidates()

StatusCode TrigMultiTrkComboHypo::findMuTrkCandidates ( TrigMultiTrkState< xAOD::MuonContainer > & state) const
private

Build J/psi candidates from muon from SG::View and tracks from the same SG::View, to be used for Tag-and-Probe PEB chains (similar to HLT_mu6_bJpsimutrk_MuonTrkPEB_L1MU5VF)

Definition at line 766 of file TrigMultiTrkComboHypo.cxx.

766 {
767
768 ATH_MSG_DEBUG( "Try to find muon + track combinations from the same SG::View" );
769
770 auto& muons = state.leptons();
771 muons.clear();
772
773 const std::vector<double> particleMasses(2, PDG::mMuon);
774
775 for (const Decision* decision : state.previousDecisions()) {
776 if (!TrigCompositeUtils::isAnyIDPassing(decision, m_allowedIDs)) continue;
777
778 auto decisionEL = TrigCompositeUtils::decisionToElementLink(decision, state.context());
780 auto muonEL = decision->objectLink<xAOD::MuonContainer>(TrigCompositeUtils::featureString());
781 const auto muon = *muonEL;
782 if (!muon->trackParticle(xAOD::Muon::TrackParticleType::CombinedTrackParticle)) continue;
783 if (!muon->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle)) continue;
784 const auto muonInDetTrack = muon->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle);
785 auto muonMomentum = muonInDetTrack->genvecP4();
786 muonMomentum.SetM(PDG::mMuon);
787
788 // add muon from decision to state.leptons
791 muons.push_back({muonEL, std::vector<ElementLink<DecisionContainer>>(1, decisionEL), decisionIDs});
792
793 ATH_MSG_DEBUG( "Found muon (CombinedTrackParticle): " << muon->pt() << " / " << muon->eta() << " / " << muon->phi() << " / " << muon->charge() );
794
796 ATH_CHECK( viewLinkInfo.isValid() );
797 auto view = *viewLinkInfo.link;
798
799 auto tracksHandle = ViewHelper::makeHandle(view, m_trackParticleContainerKey, state.context());
800 ATH_CHECK( tracksHandle.isValid() );
801 ATH_MSG_DEBUG( "Tracks container " << m_trackParticleContainerKey << " size: " << tracksHandle->size() );
802
803 // try to fit muon and track into common vertex: first track is always muon, second tracks comes from the same SG::View
804 std::vector<ElementLink<xAOD::TrackParticleContainer>> tracklist(2);
805 tracklist[0] = muon->inDetTrackParticleLink();
806 for (size_t idx = 0; idx < tracksHandle->size(); ++idx) {
807 const xAOD::TrackParticle* track = tracksHandle->at(idx);
808
809 if (track->pt() < m_trkPt[0][1] || isIdenticalTracks(track, muonInDetTrack)) continue;
810 auto trackMomentum = track->genvecP4();
812 if (!isInMassRange((muonMomentum + trackMomentum).M(), 0)) continue;
813 if (m_nTrkCharge[0] >= 0 && muonInDetTrack->charge() * track->charge() > 0.) continue;
814
815 tracklist[1] = ViewHelper::makeLink<xAOD::TrackParticleContainer>(view, tracksHandle, idx);
816
817 ATH_MSG_DEBUG( "Dump found muon+track pair before vertex fit: pT / eta / phi / charge" << endmsg <<
818 " muon: " << muonMomentum.Pt() << " / " << muonMomentum.Eta() << " / " << muonMomentum.Phi() << " / " << muon->charge() << endmsg <<
819 " track: " << trackMomentum.Pt() << " / " << trackMomentum.Eta() << " / " << trackMomentum.Phi() << " / " << track->charge() );
820
821 auto fitterState = m_vertexFitter->makeState(state.context());
822 auto vertex = fit(tracklist, particleMasses, *fitterState);
823 if (!vertex) continue;
824 xAOD::TrigBphys* trigBphys = makeTrigBPhys(*vertex, particleMasses, state.beamSpot(), *fitterState);
825 // trigBphys->setRoiId(initialRoI->roiWord());
826 state.addTrigBphysObject(trigBphys, std::vector<size_t>(1, muons.size() - 1));
827 }
828 }
829
830 return StatusCode::SUCCESS;
831}
#define endmsg
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_trackParticleContainerKey
bool isIdenticalTracks(const xAOD::TrackParticle *lhs, const xAOD::TrackParticle *rhs) const
Attempts to identify identical tracks by selection on DeltaR.
const std::string & viewString()
ElementLink< T > makeLink(const SG::View *view, const SG::ReadHandle< T > &handle, size_t index)
Create EL to a collection in view.
Definition ViewHelper.h:309
auto makeHandle(const SG::View *view, const KEY &key, const EventContext &ctx)
Create a view handle from a handle key.
Definition ViewHelper.h:273
const std::string trackMomentum
TrackParticle_v1 TrackParticle
Reference the current persistent version:
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".

◆ fit()

std::unique_ptr< xAOD::Vertex > TrigMultiTrkComboHypo::fit ( const std::vector< ElementLink< xAOD::TrackParticleContainer > > & trackParticleLinks,
const std::vector< double > & particleMasses,
Trk::IVKalState & fitterState ) const
private

Perform a vertex fit on selected tracks.

Parameters
trackParticleLinksthe trackParticles to fit
trkMassTrack mass hypothesis for mass calculations
fitterStateThe temporary state object
Returns
The fitted vertex - null if fit fails or is very low quality

Definition at line 927 of file TrigMultiTrkComboHypo.cxx.

930 {
931
932 ATH_MSG_DEBUG( "Perform vertex fit" );
933 std::vector<const xAOD::TrackParticle*> tracklist(trackParticleLinks.size(), nullptr);
934 std::transform(trackParticleLinks.begin(), trackParticleLinks.end(), tracklist.begin(),
935 [](const ElementLink<xAOD::TrackParticleContainer>& link){ return *link; });
936
937 const Trk::Perigee& perigee1 = tracklist[0]->perigeeParameters();
938 const Trk::Perigee& perigee2 = tracklist[1]->perigeeParameters();
939 int flag = 0;
940 int errorcode = 0;
941 Amg::Vector3D startingPoint = m_vertexPointEstimator->getCirclesIntersectionPoint(&perigee1, &perigee2, flag, errorcode);
942 if (errorcode != 0) startingPoint = Amg::Vector3D::Zero(3);
943 ATH_MSG_DEBUG( "Starting point: (" << startingPoint(0) << ", " << startingPoint(1) << ", " << startingPoint(2) << ")" );
944
945 m_vertexFitter->setMassInputParticles(particleMasses, fitterState);
946 std::unique_ptr<xAOD::Vertex> vertex(m_vertexFitter->fit(tracklist, startingPoint, fitterState));
947 if (!vertex) {
948 ATH_MSG_DEBUG( "Vertex fit fails" );
949 return vertex;
950 }
951 if (vertex->chiSquared() > m_chi2) {
952 ATH_MSG_DEBUG( "Fit is successful, but vertex chi2 is too high, we are not going to save it (chi2 = " << vertex->chiSquared() << " > " << m_chi2.value() << ")" );
953 vertex.reset();
954 return vertex;
955 }
956 ATH_MSG_DEBUG( "Fit is successful" );
957 vertex->clearTracks();
958 vertex->setTrackParticleLinks(trackParticleLinks);
959 return vertex;
960}
Gaudi::Property< float > m_chi2
ToolHandle< InDet::VertexPointEstimator > m_vertexPointEstimator
Eigen::Matrix< double, 3, 1 > Vector3D
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
bool flag
Definition master.py:29
std::vector< ElementLink< xAOD::TrackParticleContainer > > trackParticleLinks(const xAOD::TauJet *tau, xAOD::TauJetParameters::TauTrackFlag flag=xAOD::TauJetParameters::TauTrackFlag::classifiedCharged)

◆ hypoTools() [1/2]

ToolHandleArray< ComboHypoToolBase > & ComboHypo::hypoTools ( )
inlineprotectedinherited

Definition at line 45 of file ComboHypo.h.

45{ return m_hypoTools; }
ToolHandleArray< ComboHypoToolBase > m_hypoTools
Definition ComboHypo.h:104

◆ hypoTools() [2/2]

const ToolHandleArray< ComboHypoToolBase > & ComboHypo::hypoTools ( ) const
inlineprotectedinherited

Definition at line 46 of file ComboHypo.h.

46{ return m_hypoTools; }

◆ initialize()

StatusCode TrigMultiTrkComboHypo::initialize ( )
overridevirtual

Reimplemented from ComboHypo.

Definition at line 71 of file TrigMultiTrkComboHypo.cxx.

71 {
72 ATH_MSG_DEBUG( "TrigMultiTrkHypo::initialize()" );
73
75
76 // check consistency of the properties
77 ATH_CHECK( !m_nTrk.empty() );
78
79 if (m_nTrkCharge.empty()) {
80 ATH_MSG_INFO( "totalCharge value is not specified, no charge selection for track combinations will be used" );
81 m_nTrkCharge = std::vector<int>(m_nTrk.size(), -1);
82 }
83
84 if (m_trkMass.empty()) {
85 ATH_MSG_INFO( "trackMasses value is not specified, muon/electron mass will be used" );
86 for (const auto& n : m_nTrk) {
87 m_trkMass.value().emplace_back(std::vector<double>(n, (m_doElectrons ? PDG::mElectron : PDG::mMuon)));
88 }
89 }
90 if (m_trkPt.empty()) {
91 ATH_MSG_INFO( "trackPtThresholds value is not specified" );
92 for (const auto& n : m_nTrk) {
93 m_trkPt.value().emplace_back(std::vector<double>(n, -100.));
94 }
95 }
96 m_trkPtMin = m_trkPt[0].back();
97 for (size_t i = 0; i < m_trkPt.size(); ++i) {
98 m_trkPtMin = std::min(m_trkPtMin, m_trkPt[i].back());
99 }
100
101 ATH_CHECK( m_trkMass.size() == m_nTrk.size() );
102 ATH_CHECK( m_trkPt.size() == m_nTrk.size() );
103
104 for (size_t i = 0; i < m_nTrk.size(); ++i) {
105 ATH_CHECK( m_trkMass[i].size() == m_nTrk[i] );
106 ATH_CHECK( m_trkPt[i].size() == m_nTrk[i] );
107 ATH_CHECK( std::is_sorted(m_trkPt[i].begin(), m_trkPt[i].end(), std::greater<>()) );
108 }
109
110 for (const auto& range : m_massRange.value()) {
111 ATH_CHECK( range.first < range.second );
112 }
113
114 // dump numerical values
115 if (msgLvl(MSG::DEBUG)) {
116 for (size_t i = 0; i < m_nTrk.size(); ++i) {
117 ATH_MSG_DEBUG( "vertex topology: nTrk = " << m_nTrk[i] );
118 for (size_t j = 0; j < m_nTrk[i]; ++j) {
119 ATH_MSG_DEBUG( " " << j + 1 << " trk: mass = " << m_trkMass[i][j] << ", Pt > " << m_trkPt[i][j] );
120 }
121 }
122 msg() << MSG::DEBUG << " mass range: {";
123 for (const auto& range : m_massRange.value()) {
124 msg() << MSG::DEBUG << " { " << range.first << ", " << range.second << " }";
125 }
126 msg() << MSG::DEBUG << " }" << std::endl;
127 }
128
129 if (m_isStreamer) {
130 ATH_MSG_DEBUG( "Configured to run in a streamer mode: no trigger objects will be created" );
131 }
132 if (!m_isStreamer && m_trigLevel != "EF") {
133 ATH_MSG_ERROR( "Could not create trigger objects from tracks or L2 CB muons, use the streamer mode for L2 step" );
134 return StatusCode::FAILURE;
135 }
136
137 ATH_CHECK( !((m_trigLevel == "L2IO" || m_trigLevel == "L2MT") && m_doElectrons) );
138
139 if (m_trigLevel == "L2" || (m_trigLevel == "EF" && m_isMuTrkMode)) {
142 }
143 else if (m_trigLevel == "L2IO" || m_trigLevel == "L2MT" || m_trigLevel == "EF") {
144 ATH_CHECK( m_trackParticleContainerKey.initialize(false) );
145 }
146 else {
147 ATH_MSG_ERROR( "trigLevel should be L2, L2IO, L2MT or EF, but " << m_trigLevel << " provided" );
148 return StatusCode::FAILURE;
149 }
150
151 ATH_CHECK( m_trigBphysContainerKey.initialize(!m_isStreamer.value()) );
152 ATH_CHECK( m_beamSpotKey.initialize(!m_isStreamer.value()) ); // need beamSpot only to create xAOD::TrigBphys object
153
154 ATH_CHECK( m_vertexFitter.retrieve() );
155 ATH_CHECK( m_vertexPointEstimator.retrieve() );
156 ATH_CHECK( m_v0Tools.retrieve() );
157
158 // allowed IDs to filter out incoming decisions at L2 level
159 for (const auto& item : triggerMultiplicityMap()) {
160 const HLT::Identifier id = HLT::Identifier(item.first);
161 m_allowedIDs.insert(id.numeric());
162 if (item.second.size() > 1) {
163 for (size_t i = 0; i < item.second.size(); i++) {
164 m_allowedIDs.insert(TrigCompositeUtils::createLegName(id, i).numeric());
165 }
166 }
167 }
168 if (msgLvl(MSG::DEBUG)) {
169 ATH_MSG_DEBUG( "Allowed decisions:" );
170 for (const DecisionID& id : m_allowedIDs) {
171 ATH_MSG_DEBUG( " +++ " << HLT::Identifier(id) );
172 }
173 }
174 if (m_doElectrons) {
175 for (const DecisionID& id : m_allowedIDs) {
176 std::string name = HLT::Identifier(id).name();
177 bool isMergedElectronChain = false;
178 for (size_t i = 0; i < m_mergedElectronChains.size(); ++i) {
179 if (name.find(m_mergedElectronChains.value().at(i)) != std::string::npos) {
180 isMergedElectronChain = true;
181 break;
182 }
183 }
184 (isMergedElectronChain ? m_mergedElectronIDs.insert(id) : m_resolvedElectronIDs.insert(id));
185 }
186 if (msgLvl(MSG::DEBUG)) {
187 ATH_MSG_DEBUG( "Well-separated electron decisions:" );
188 for (const DecisionID& id : m_resolvedElectronIDs) {
189 ATH_MSG_DEBUG( " +++ " << HLT::Identifier(id) );
190 }
191 ATH_MSG_DEBUG( "Overlapping electron decisions:" );
192 for (const DecisionID& id : m_mergedElectronIDs) {
193 ATH_MSG_DEBUG( " +++ " << HLT::Identifier(id) );
194 }
195 }
196 }
197
198 if (!m_monTool.empty()) {
199 ATH_CHECK( m_monTool.retrieve() );
200 ATH_MSG_DEBUG( "GenericMonitoringTool name:" << m_monTool );
201 }
202 else {
203 ATH_MSG_DEBUG( "No GenericMonitoringTool configured: no monitoring histograms will be available" );
204 }
205
206 return StatusCode::SUCCESS;
207}
#define ATH_MSG_INFO(x)
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)
const Combo::MultiplicityReqMap & triggerMultiplicityMap() const
Definition ComboHypo.h:43
virtual StatusCode initialize() override
Definition ComboHypo.cxx:22
TrigCompositeUtils::DecisionIDContainer m_mergedElectronIDs
ToolHandle< Trk::V0Tools > m_v0Tools
Gaudi::Property< std::vector< std::string > > m_mergedElectronChains
Gaudi::Property< std::vector< std::pair< double, double > > > m_massRange
TrigCompositeUtils::DecisionIDContainer m_resolvedElectronIDs
HLT::Identifier createLegName(const HLT::Identifier &chainIdentifier, size_t counter)
Generate the HLT::Identifier which corresponds to a specific leg of a given chain.

◆ 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()

◆ isIdenticalTracks() [1/3]

bool TrigMultiTrkComboHypo::isIdenticalTracks ( const xAOD::Electron * lhs,
const xAOD::Electron * rhs ) const
private

Definition at line 1027 of file TrigMultiTrkComboHypo.cxx.

1027 {
1028
1029 return isIdenticalTracks(*lhs->trackParticleLink(), *rhs->trackParticleLink());
1030}
const ElementLink< TrackParticleContainer > & trackParticleLink(size_t index=0) const
ElementLink to the xAOD::TrackParticle/s that match the electron candidate.

◆ isIdenticalTracks() [2/3]

bool TrigMultiTrkComboHypo::isIdenticalTracks ( const xAOD::Muon * lhs,
const xAOD::Muon * rhs ) const
private

◆ isIdenticalTracks() [3/3]

bool TrigMultiTrkComboHypo::isIdenticalTracks ( const xAOD::TrackParticle * lhs,
const xAOD::TrackParticle * rhs ) const
private

Attempts to identify identical tracks by selection on DeltaR.

Returns
true if 'identical', false otherwise

Definition at line 1015 of file TrigMultiTrkComboHypo.cxx.

1015 {
1016
1017 if (lhs->charge() * rhs->charge() < 0.) return false;
1018 return (ROOT::Math::VectorUtil::DeltaR(lhs->genvecP4(), rhs->genvecP4()) < m_deltaR);
1019}
Gaudi::Property< float > m_deltaR
GenVecFourMom_t genvecP4() const
The full 4-momentum of the particle : GenVector form.
float charge() const
Returns the charge.

◆ isInMassRange()

bool TrigMultiTrkComboHypo::isInMassRange ( double mass,
size_t idx ) const
private

Definition at line 1033 of file TrigMultiTrkComboHypo.cxx.

1033 {
1034
1035 const auto& range = m_massRange[idx];
1036 return (mass > range.first && mass < range.second);
1037}

◆ legToInputCollectionMap()

const Combo::LegMap & ComboHypo::legToInputCollectionMap ( ) const
inlineprotectedinherited

Definition at line 44 of file ComboHypo.h.

44{ return m_legToInputCollectionMap.value(); }

◆ makeTrigBPhys()

xAOD::TrigBphys * TrigMultiTrkComboHypo::makeTrigBPhys ( const xAOD::Vertex & vertex,
const std::vector< double > & particleMasses,
const xAOD::Vertex & beamSpot,
const Trk::IVKalState & fitterState ) const
private

Construct the trigger object that may be stored for debugging or matching.

Parameters
vertexFitted Vertex object
particleMassesTrack mass hypothesis for mass calculations
beamSpotThe beamspot vertex
fitterStateVertexer state object
Returns
Pointer to the TrigBphys object

Definition at line 963 of file TrigMultiTrkComboHypo.cxx.

967 {
968
969 double invariantMass = 0.;
970 double invariantMassError = 0.;
971 if (!m_vertexFitter->VKalGetMassError(invariantMass, invariantMassError, fitterState).isSuccess()) {
972 ATH_MSG_DEBUG( "Warning from TrkVKalVrtFitter: can not calculate uncertainties" );
973 invariantMass = -9999.;
974 }
975
977 for (size_t i = 0; i < vertex.nTrackParticles(); ++i) {
978 auto p = vertex.trackParticle(i)->genvecP4();
979 p.SetM(particleMasses[i]);
980 momentum += p;
981 }
982
984 result->makePrivateStore();
985 result->initialise(0, momentum.Eta(), momentum.Phi(), momentum.Pt(), xAOD::TrigBphys::MULTIMU, momentum.M(), xAOD::TrigBphys::EF);
986
987 if (m_doElectrons) result->setParticleType(xAOD::TrigBphys::JPSIEE);
988 result->setFitmass(invariantMass);
989 result->setFitchi2(vertex.chiSquared());
990 result->setFitndof(vertex.numberDoF());
991 result->setFitx(vertex.x());
992 result->setFity(vertex.y());
993 result->setFitz(vertex.z());
994 result->setTrackParticleLinks(vertex.trackParticleLinks());
995 result->setLxy(m_v0Tools->lxy(&vertex, &beamSpot));
996 result->setLxyError(m_v0Tools->lxyError(&vertex, &beamSpot));
997
999 "TrigBphys objects:\n\t " <<
1000 "roiId: " << result->roiId() << "\n\t " <<
1001 "particleType: " << result->particleType() << "\n\t " <<
1002 "level: " << result->level() << "\n\t " <<
1003 "eta: " << result->eta() << "\n\t " <<
1004 "phi: " << result->phi() << "\n\t " <<
1005 "mass: " << result->mass() << "\n\t " <<
1006 "fitmass: " << result->fitmass() << "\n\t " <<
1007 "chi2/NDF: " << result->fitchi2() << " / " << result->fitndof() << "\n\t " <<
1008 "vertex: (" << result->fitx() << ", " << result->fity() << ", " << result->fitz() << ")\n\t " <<
1009 "Lxy/LxyError: " << result->lxy() << " / " << result->lxyError() );
1010
1011 return result;
1012}

◆ mergeLeptonsFromDecisions() [1/2]

template<typename CONTAINER>
StatusCode TrigMultiTrkComboHypo::mergeLeptonsFromDecisions ( TrigMultiTrkState< CONTAINER > & ) const
private

Go through state.previousDecisions(), fetch xAOD::Muons/xAODElectron objects attached to decisions and save links to them in state.leptons().

◆ mergeLeptonsFromDecisions() [2/2]

template<typename T>
StatusCode TrigMultiTrkComboHypo::mergeLeptonsFromDecisions ( TrigMultiTrkState< T > & state) const

Definition at line 280 of file TrigMultiTrkComboHypo.cxx.

280 {
281
282 auto& leptons = state.leptons();
283 leptons.clear();
284
285 std::vector<const Decision*> previousDecisions(state.previousDecisions().begin(), state.previousDecisions().end());
286 std::map<const Decision*, int> decisionToInputCollectionIndexMap;
287 for (const auto& decision : previousDecisions) decisionToInputCollectionIndexMap.emplace(decision, 0);
289 for (size_t k = 1; k < decisionsInput().size(); ++k) {
290 auto previousDecisionsHandle = SG::makeHandle(decisionsInput().at(k), state.context());
291 CHECK( previousDecisionsHandle.isValid() );
292 ATH_MSG_DEBUG( "Adding " << previousDecisionsHandle->size() << " decisions from " << decisionsInput().at(k).key() );
293 for (const Decision* decision : *previousDecisionsHandle) {
294 previousDecisions.push_back(decision);
295 decisionToInputCollectionIndexMap.emplace(decision, static_cast<int>(k));
296 }
297 }
298 }
299
300 // all muons/electrons from views are already connected with previous decisions by TrigMuonEFHypoAlg
301 for (const Decision* decision : previousDecisions) {
303
304 ElementLink<T> leptonEL;
305 if (decision->hasObjectLink(TrigCompositeUtils::featureString(), ClassID_traits<T>::ID())) {
306 leptonEL = decision->objectLink<T>(TrigCompositeUtils::featureString());
307 }
308 else {
309 auto leptonLinkInfo = TrigCompositeUtils::findLink<T>(decision, TrigCompositeUtils::featureString(), true);
310 ATH_CHECK( leptonLinkInfo.isValid() );
311 leptonEL = leptonLinkInfo.link;
312 }
313
314 const auto lepton = *leptonEL;
315 if constexpr(std::is_same<T, xAOD::MuonContainer>::value) {
316 if (!lepton->trackParticle(xAOD::Muon::TrackParticleType::CombinedTrackParticle)) continue;
317 if (!lepton->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle)) continue;
318 }
319 else if constexpr(std::is_same<T, xAOD::ElectronContainer>::value) {
320 if (!lepton->trackParticle()) continue;
321 }
322 else {
323 ATH_MSG_ERROR( "mergeLeptonsFromDecisions(): no scenario for the provided CONTAINER is specified" );
324 return StatusCode::FAILURE;
325 }
326 auto decisionEL = TrigCompositeUtils::decisionToElementLink(decision, state.context());
327
328 auto itr = leptons.end();
330 itr = std::find_if(leptons.begin(), leptons.end(),
331 [this, lepton = lepton](const auto& x){ return this->isIdenticalTracks(lepton, *x.link); });
332 }
333 if (itr == leptons.end()) {
334 leptons.push_back({leptonEL, std::vector<ElementLink<DecisionContainer>>(1, decisionEL), DecisionIDContainer()});
335 }
336 else {
337 (*itr).decisionLinks.push_back(decisionEL);
338 }
339 }
340
341 // muon->pt() is equal to muon->trackParticle(xAOD::Muon::TrackParticleType::CombinedTrackParticle)->pt()
342 // and the later is used by TrigMuonEFHypoTool for classification of muEFCB candidates
343 std::sort(leptons.begin(), leptons.end(), [](const auto& lhs, const auto& rhs){ return ((*lhs.link)->pt() > (*rhs.link)->pt()); });
344
345 // for each muon we extract DecisionIDs stored in the associated Decision objects and copy them at muon.decisionIDs
346 const auto& legToInputCollectionIndexMap = legToInputCollectionMap();
347 for (auto& item : leptons) {
348 for (const ElementLink<xAOD::TrigCompositeContainer>& decisionEL : item.decisionLinks) {
350 auto decisionIndex = decisionToInputCollectionIndexMap[*decisionEL];
351 for (const auto& id : TrigCompositeUtils::decisionIDs(*decisionEL)) {
353 auto legIndex = static_cast<size_t>(TrigCompositeUtils::getIndexFromLeg(id));
354 std::string chain = TrigCompositeUtils::getIDFromLeg(id).name();
355 if (legToInputCollectionIndexMap.at(chain).at(legIndex) == decisionIndex) item.decisionIDs.insert(id);
356 }
357 }
358 else {
359 TrigCompositeUtils::decisionIDs(*decisionEL, item.decisionIDs);
360 }
361 }
362 }
363
364 if (msgLvl(MSG::DEBUG)) {
365 ATH_MSG_DEBUG( "Dump found leptons before vertex fit: " << leptons.size() << " candidates" );
366 for (const auto& item : leptons) {
367 const auto lepton = *item.link;
369 if constexpr(std::is_same<T, xAOD::MuonContainer>::value) {
370 track = lepton->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle);
371 }
372 else {
373 track = lepton->trackParticle();
374 }
375 ATH_MSG_DEBUG( " -- lepton InDetTrackParticle pt/eta/phi/q: " << track->pt() << " / " << track->eta() << " / " << track->phi() << " / " << track->charge() );
376 ATH_MSG_DEBUG( " lepton pt (muon: CombinedTrackParticle): " << lepton->pt() << " / " << lepton->eta() << " / " << lepton->phi() << " / " << lepton->charge() );
377 ATH_MSG_DEBUG( " allowed decisions:" );
378 for (const DecisionID& id : item.decisionIDs) {
379 ATH_MSG_DEBUG( " +++ " << HLT::Identifier(id) );
380 }
381 }
382 }
383
384 return StatusCode::SUCCESS;
385}
const Combo::LegMap & legToInputCollectionMap() const
Definition ComboHypo.h:44
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
std::string name() const
reports human redable name
Gaudi::Property< bool > m_applyOverlapRemoval
unsigned long long T
int32_t getIndexFromLeg(const HLT::Identifier &legIdentifier)
Extract the numeric index of a leg identifier.
bool passed(DecisionID id, const DecisionIDContainer &idSet)
checks if required decision ID is in the set of IDs in the container
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.

◆ mergeTracksFromDecisions()

template<typename CONTAINER>
StatusCode TrigMultiTrkComboHypo::mergeTracksFromDecisions ( TrigMultiTrkStateBase & state) const
private

Go through state.previousDecisions(), fetch xAOD::Muons/xAODElectron objects attached to decisions and save links to the their xAOD::TrackParticle objects in state.tracks().

Definition at line 439 of file TrigMultiTrkComboHypo.cxx.

439 {
440
441 auto& tracks = state.tracks();
442 tracks.clear();
443
444 // all muons/electrons from views are already connected with previous decisions by TrigMuonEFHypoAlg or by TrigEgammaPrecisionElectronHypoAlg
445 for (const Decision* decision : state.previousDecisions()) {
446 if (!TrigCompositeUtils::isAnyIDPassing(decision, m_allowedIDs)) continue;
447
449 auto leptonEL = decision->objectLink<CONTAINER>(TrigCompositeUtils::featureString());
450 const auto lepton = *leptonEL;
451
452 ElementLink<xAOD::TrackParticleContainer> trackEL;
453 if constexpr(std::is_same<CONTAINER, xAOD::MuonContainer>::value) {
454 if (!lepton->trackParticle(xAOD::Muon::TrackParticleType::CombinedTrackParticle)) continue;
455 if (!lepton->trackParticle(xAOD::Muon::TrackParticleType::InnerDetectorTrackParticle)) continue;
456 trackEL = lepton->inDetTrackParticleLink();
457 }
458 else if constexpr(std::is_same<CONTAINER, xAOD::L2CombinedMuonContainer>::value) {
459 if (!lepton->idTrack()) continue;
460 trackEL = lepton->idTrackLink();
461 }
462 else if constexpr(std::is_same<CONTAINER, xAOD::ElectronContainer>::value) {
463 if (!lepton->trackParticle()) continue;
464 trackEL = lepton->trackParticleLink();
465 }
466 else {
467 ATH_MSG_ERROR( "mergeTracksFromDecisions(): no scenario for the provided CONTAINER is specified" );
468 return StatusCode::FAILURE;
469 }
470
471 if (!trackEL.isValid()) continue;
473 std::find_if(tracks.begin(), tracks.end(),
474 [this, track = *trackEL](const auto& x){ return this->isIdenticalTracks(track, *x); }) == tracks.end()) {
475 tracks.emplace_back(trackEL);
476 }
477 }
478 std::sort(tracks.begin(), tracks.end(), [](const auto& lhs, const auto& rhs){ return ((*lhs)->pt() > (*rhs)->pt()); });
479
480 if (msgLvl(MSG::DEBUG)) {
481 ATH_MSG_DEBUG( "Dump found tracks before vertex fit: " << tracks.size() << " candidates" );
482 for (const auto& trackEL : tracks) {
483 const xAOD::TrackParticle* track = *trackEL;
484 ATH_MSG_DEBUG( " -- track pt/eta/phi/q: " << track->pt() << " / " << track->eta() << " / " << track->phi() << " / " << track->charge() );
485 }
486 }
487
488 return StatusCode::SUCCESS;
489}

◆ mergeTracksFromViews()

StatusCode TrigMultiTrkComboHypo::mergeTracksFromViews ( TrigMultiTrkStateBase & state) const
private

Go through state.previousDecisions() and fetch xAOD::TrackParticle objects associated with the nearest SG::View.

Enable overlap removal to get collection of unique objects at state.tracks().

Definition at line 388 of file TrigMultiTrkComboHypo.cxx.

388 {
389
390 auto& tracks = state.tracks();
391 tracks.clear();
392
393 std::set<const SG::View*> views;
394 for (const Decision* decision : state.previousDecisions()) {
395 if (!TrigCompositeUtils::isAnyIDPassing(decision, m_allowedIDs)) continue;
396
398 ATH_CHECK( viewLinkInfo.isValid() );
399 const SG::View* view = *viewLinkInfo.link;
400 if (views.find(view) != views.end()) continue; // tracks from this view have already been fetched
401
402 auto tracksHandle = ViewHelper::makeHandle(view, m_trackParticleContainerKey, state.context());
403 ATH_CHECK( tracksHandle.isValid() );
404 ATH_MSG_DEBUG( "tracks handle " << m_trackParticleContainerKey << " size: " << tracksHandle->size() );
405
406 std::vector<ElementLink<xAOD::TrackParticleContainer>> tracksFromView;
407 tracksFromView.reserve(tracksHandle->size());
408 for (size_t idx = 0; idx < tracksHandle->size(); ++idx) {
409 tracksFromView.emplace_back(ViewHelper::makeLink<xAOD::TrackParticleContainer>(view, tracksHandle, idx));
410 }
411
412 for (const auto& trackEL : tracksFromView) {
413 const xAOD::TrackParticle* track = *trackEL;
414 if (track->definingParametersCovMatrixVec().empty() || track->pt() < m_trkPtMin) continue;
415
416 if (views.empty() || !m_applyOverlapRemoval ||
417 std::find_if(tracks.begin(), tracks.end(),
418 [this, track](const auto& x){ return isIdenticalTracks(track, *x); }) == tracks.end()) {
419 tracks.emplace_back(trackEL);
420 }
421 }
422 views.insert(view);
423 }
424 std::sort(tracks.begin(), tracks.end(), [](const auto& lhs, const auto& rhs){ return ((*lhs)->pt() > (*rhs)->pt()); });
425
426 if (msgLvl(MSG::DEBUG)) {
427 ATH_MSG_DEBUG( "Dump found tracks before vertex fit: " << tracks.size() << " candidates" );
428 for (const auto& trackEL : tracks) {
429 const xAOD::TrackParticle* track = *trackEL;
430 ATH_MSG_DEBUG( " -- track pt/eta/phi/q: " << track->pt() << " / " << track->eta() << " / " << track->phi() << " / " << track->charge() );
431 }
432 }
433
434 return StatusCode::SUCCESS;
435}

◆ 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.

◆ passedDeltaRcut()

bool TrigMultiTrkComboHypo::passedDeltaRcut ( const std::vector< xAOD::TrackParticle::GenVecFourMom_t > & momenta) const
private

Definition at line 1040 of file TrigMultiTrkComboHypo.cxx.

1040 {
1041
1042 if (m_deltaRMax == std::numeric_limits<float>::max() && m_deltaRMin == std::numeric_limits<float>::lowest()) {
1043 return true;
1044 }
1045 for (size_t i = 0; i < p.size(); ++i) {
1046 for (size_t j = i + 1; j < p.size(); ++j) {
1047 double deltaR = ROOT::Math::VectorUtil::DeltaR(p[i], p[j]);
1048 if (deltaR > m_deltaRMax || deltaR < m_deltaRMin) return false;
1049 }
1050 }
1051 return true;
1052}
Scalar deltaR(const MatrixBase< Derived > &vec) const
Gaudi::Property< float > m_deltaRMin
Gaudi::Property< float > m_deltaRMax

◆ processMergedElectrons()

StatusCode TrigMultiTrkComboHypo::processMergedElectrons ( TrigMultiTrkState< xAOD::ElectronContainer > & state) const
private

Make all possible combinations from electrons originating from the same BPH-0DR3-EM7J15 cluster, only one track should pass e5_lhvloose requirements; initialRoI word will be saved into xAOD::TrigBphys objects.

Definition at line 706 of file TrigMultiTrkComboHypo.cxx.

706 {
707
708 ATH_MSG_DEBUG( "Try to find electrons originating from the same EM cluster" );
709
710 // some electrons can be already attached to the list, add electrons from BPH-0DR3-EM7J15 clusters to the end
711 auto& leptons = state.leptons();
712
713 if (m_mergedElectronIDs.empty()) {
714 ATH_MSG_DEBUG( "no chains similar to BPH-0DR3-EM7J15 have been requested, should not look for close-by electrons" );
715 return StatusCode::SUCCESS;
716 }
717
718 const std::vector<double> particleMasses(2, PDG::mElectron);
719
720 for (const Decision* decision : state.previousDecisions()) {
722
723 auto decisionEL = TrigCompositeUtils::decisionToElementLink(decision, state.context());
725 auto electronEL = decision->objectLink<xAOD::ElectronContainer>(TrigCompositeUtils::featureString());
726 const auto electron = *electronEL;
727 if (!electron->trackParticle()) continue;
728
729 // get all electrons from the same SG::View, i.e. from the same initialRoI
730 const auto electronContainer = electronEL.getStorableObjectPointer();
731 if (electronContainer->size() <= 1) continue;
732
733 // add electron from decision to state.leptons
736 leptons.push_back({electronEL, std::vector<ElementLink<DecisionContainer>>(1, decisionEL), decisionIDs});
737
738 // get initialRoI this electron originating from
740 ATH_CHECK( roiInfo.isValid() );
741 auto initialRoI = *roiInfo.link;
742
743 // try to build di-electron pairs: first electron is always from decision, second from the same SG::View
744 std::vector<ElementLink<xAOD::TrackParticleContainer>> tracklist(2);
745 tracklist[0] = electron->trackParticleLink();
746 for (size_t i = 0; i < electronContainer->size(); ++i) {
747 const auto electronFromView = electronContainer->at(i);
748 if (electronFromView == electron) continue;
749 if (!electronFromView->trackParticle()) continue;
750 if (!electronFromView->caloCluster() || electronFromView->caloCluster()->et() < m_caloClusterEtThreshold) continue;
751 tracklist[1] = electronFromView->trackParticleLink();
752
753 auto fitterState = m_vertexFitter->makeState(state.context());
754 auto vertex = fit(tracklist, particleMasses, *fitterState);
755 if (!vertex) continue;
756 xAOD::TrigBphys* trigBphys = makeTrigBPhys(*vertex, particleMasses, state.beamSpot(), *fitterState);
757 trigBphys->setRoiId(initialRoI->roiWord());
758 state.addTrigBphysObject(trigBphys, std::vector<size_t>(1, leptons.size() - 1));
759 }
760 }
761
762 return StatusCode::SUCCESS;
763}
xAOD::ElectronContainer * electronContainer
const T * at(size_type n) const
Access an element, as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
Gaudi::Property< double > m_caloClusterEtThreshold
void setRoiId(uint32_t id)
set method: roiId
ElectronContainer_v1 ElectronContainer
Definition of the current "electron container version".

◆ 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 }

◆ 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}
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.

◆ triggerMultiplicityMap()

const Combo::MultiplicityReqMap & ComboHypo::triggerMultiplicityMap ( ) const
inlineprotectedinherited

Definition at line 43 of file ComboHypo.h.

43{ return m_multiplicitiesReqMap.value(); }
Gaudi::Property< Combo::MultiplicityReqMap > m_multiplicitiesReqMap
Definition ComboHypo.h:56

◆ 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

◆ m_allowedIDs

TrigCompositeUtils::DecisionIDContainer TrigMultiTrkComboHypo::m_allowedIDs
private

Definition at line 309 of file TrigMultiTrkComboHypo.h.

◆ m_applyOverlapRemoval

Gaudi::Property<bool> TrigMultiTrkComboHypo::m_applyOverlapRemoval
private
Initial value:
{this, "applyOverlapRemoval", true,
"apply overlap removal for the close-by same-sign objects from different views"}

Definition at line 274 of file TrigMultiTrkComboHypo.h.

274 {this, "applyOverlapRemoval", true,
275 "apply overlap removal for the close-by same-sign objects from different views"};

◆ m_beamSpotKey

SG::ReadCondHandleKey<InDet::BeamSpotData> TrigMultiTrkComboHypo::m_beamSpotKey {this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot"}
private

Definition at line 262 of file TrigMultiTrkComboHypo.h.

262{this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot"};

◆ m_caloClusterEtThreshold

Gaudi::Property<double> TrigMultiTrkComboHypo::m_caloClusterEtThreshold
private
Initial value:
{this, "caloClusterEtThreshold", 5.,
"minimum transverse energy of the cluster, associated with close-by electron"}

Definition at line 300 of file TrigMultiTrkComboHypo.h.

300 {this, "caloClusterEtThreshold", 5.,
301 "minimum transverse energy of the cluster, associated with close-by electron"};

◆ m_checkMultiplicity

Gaudi::Property<bool> TrigMultiTrkComboHypo::m_checkMultiplicity
private
Initial value:
{this, "checkMultiplicity", false,
"check that we have enough leptons to fire the chain, needed for HLT_mu6_2mu4_bJpsi_L1MU5VF_3MU3VF"}

Definition at line 280 of file TrigMultiTrkComboHypo.h.

280 {this, "checkMultiplicity", false,
281 "check that we have enough leptons to fire the chain, needed for HLT_mu6_2mu4_bJpsi_L1MU5VF_3MU3VF"};

◆ m_checkMultiplicityMap

Gaudi::Property<bool> ComboHypo::m_checkMultiplicityMap
privateinherited
Initial value:
{ this, "CheckMultiplicityMap", true,
"Perform a consistency check of the MultiplicitiesMap"}

Definition at line 62 of file ComboHypo.h.

62 { this, "CheckMultiplicityMap", true,
63 "Perform a consistency check of the MultiplicitiesMap"};

◆ m_chi2

Gaudi::Property<float> TrigMultiTrkComboHypo::m_chi2
private
Initial value:
{this, "chi2", 150.,
"chi2 cut for the fitted vertex"}

Definition at line 288 of file TrigMultiTrkComboHypo.h.

288 {this, "chi2", 150.,
289 "chi2 cut for the fitted vertex"};

◆ m_combineInputDecisionCollections

Gaudi::Property<bool> TrigMultiTrkComboHypo::m_combineInputDecisionCollections
private
Initial value:
{this, "combineInputDecisionCollections", false,
"combine objects attached to decisions from different input collections, needed for HLT_mu4_ivarloose_mu4_b10invmAB120vtx20_L12MU3V chains"}

Definition at line 276 of file TrigMultiTrkComboHypo.h.

276 {this, "combineInputDecisionCollections", false,
277 "combine objects attached to decisions from different input collections, needed for HLT_mu4_ivarloose_mu4_b10invmAB120vtx20_L12MU3V chains"};

◆ m_deltaR

Gaudi::Property<float> TrigMultiTrkComboHypo::m_deltaR
private
Initial value:
{this, "deltaR", 0.01,
"minimum deltaR between same-sign tracks (overlap removal)"}

Definition at line 282 of file TrigMultiTrkComboHypo.h.

282 {this, "deltaR", 0.01,
283 "minimum deltaR between same-sign tracks (overlap removal)"};

◆ m_deltaRMax

Gaudi::Property<float> TrigMultiTrkComboHypo::m_deltaRMax
private
Initial value:
{this, "deltaRMax", std::numeric_limits<float>::max(),
"maximum deltaR between tracks in a candidate"}

Definition at line 284 of file TrigMultiTrkComboHypo.h.

284 {this, "deltaRMax", std::numeric_limits<float>::max(),
285 "maximum deltaR between tracks in a candidate"};

◆ m_deltaRMin

Gaudi::Property<float> TrigMultiTrkComboHypo::m_deltaRMin
private
Initial value:
{this, "deltaRMin", std::numeric_limits<float>::lowest(),
"maximum deltaR between tracks in a candidate"}

Definition at line 286 of file TrigMultiTrkComboHypo.h.

286 {this, "deltaRMin", std::numeric_limits<float>::lowest(),
287 "maximum deltaR between tracks in a candidate"};

◆ 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_doElectrons

Gaudi::Property<bool> TrigMultiTrkComboHypo::m_doElectrons
private
Initial value:
{this, "doElectrons", false,
"use electrons if true, otherwise use muons"}

Definition at line 294 of file TrigMultiTrkComboHypo.h.

294 {this, "doElectrons", false,
295 "use electrons if true, otherwise use muons"};

◆ 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_hypoTools

ToolHandleArray< ComboHypoToolBase > ComboHypo::m_hypoTools {this, "ComboHypoTools", {}, "Tools to perform selection"}
privateinherited

Definition at line 104 of file ComboHypo.h.

104{this, "ComboHypoTools", {}, "Tools to perform selection"};

◆ m_inputs

SG::ReadHandleKeyArray<TrigCompositeUtils::DecisionContainer> ComboHypo::m_inputs { this, "HypoInputDecisions", {}, "Input Decisions" }
privateinherited

Definition at line 50 of file ComboHypo.h.

50{ this, "HypoInputDecisions", {}, "Input Decisions" };

◆ m_isMuTrkMode

Gaudi::Property<bool> TrigMultiTrkComboHypo::m_isMuTrkMode
private
Initial value:
{this, "isMuTrkMode", false,
"make pairs between muon from SG::View and tracks from the same SG::View"}

Definition at line 292 of file TrigMultiTrkComboHypo.h.

292 {this, "isMuTrkMode", false,
293 "make pairs between muon from SG::View and tracks from the same SG::View"};

◆ m_isStreamer

Gaudi::Property<bool> TrigMultiTrkComboHypo::m_isStreamer
private
Initial value:
{this, "isStreamer", false,
"if true we do not create trigger objects, just copy all appropriate decisions to the next step or break the chain"}

Definition at line 290 of file TrigMultiTrkComboHypo.h.

290 {this, "isStreamer", false,
291 "if true we do not create trigger objects, just copy all appropriate decisions to the next step or break the chain"};

◆ m_legToInputCollectionMap

Gaudi::Property< Combo::LegMap > ComboHypo::m_legToInputCollectionMap
privateinherited
Initial value:
{this, "LegToInputCollectionMap", {},
"Map from the chain name to the per-leg index in this algorithm's ReadHandleKeyArray which should be used as the source of incoming Decision Objects on the leg."}

Definition at line 59 of file ComboHypo.h.

59 {this, "LegToInputCollectionMap", {},
60 "Map from the chain name to the per-leg index in this algorithm's ReadHandleKeyArray which should be used as the source of incoming Decision Objects on the leg."};

◆ m_massRange

Gaudi::Property<std::vector<std::pair<double, double> > > TrigMultiTrkComboHypo::m_massRange
private
Initial value:
{this, "massRange", { {0., 100000.} },
"range of the invariant mass of the track combinations"}

Definition at line 272 of file TrigMultiTrkComboHypo.h.

272 {this, "massRange", { {0., 100000.} },
273 "range of the invariant mass of the track combinations"};

◆ m_mergedElectronChains

Gaudi::Property<std::vector<std::string> > TrigMultiTrkComboHypo::m_mergedElectronChains
private
Initial value:
{this, "mergedElectronChains", {"BPH-0DR3-EM7J15"},
"patterns for BPH-0DR3-EM7J15 like chains"}

Definition at line 298 of file TrigMultiTrkComboHypo.h.

298 {this, "mergedElectronChains", {"BPH-0DR3-EM7J15"},
299 "patterns for BPH-0DR3-EM7J15 like chains"};

◆ m_mergedElectronIDs

TrigCompositeUtils::DecisionIDContainer TrigMultiTrkComboHypo::m_mergedElectronIDs
private

Definition at line 311 of file TrigMultiTrkComboHypo.h.

◆ m_monTool

ToolHandle<GenericMonitoringTool> TrigMultiTrkComboHypo::m_monTool {this, "MonTool", "", "monitoring tool"}
private

Definition at line 307 of file TrigMultiTrkComboHypo.h.

307{this, "MonTool", "", "monitoring tool"};

◆ m_multiplicitiesReqMap

Gaudi::Property< Combo::MultiplicityReqMap > ComboHypo::m_multiplicitiesReqMap
privateinherited
Initial value:
{this, "MultiplicitiesMap", {},
"Map from the chain name to multiplicities required at each input"}

Definition at line 56 of file ComboHypo.h.

56 {this, "MultiplicitiesMap", {},
57 "Map from the chain name to multiplicities required at each input"};

◆ m_nTrk

Gaudi::Property<std::vector<unsigned int> > TrigMultiTrkComboHypo::m_nTrk
private
Initial value:
{this, "nTracks", {2},
"number of tracks to be fitted into the common vertex"}

Definition at line 264 of file TrigMultiTrkComboHypo.h.

264 {this, "nTracks", {2},
265 "number of tracks to be fitted into the common vertex"};

◆ m_nTrkCharge

Gaudi::Property<std::vector<int> > TrigMultiTrkComboHypo::m_nTrkCharge
private
Initial value:
{this, "totalCharge", {},
"magnitude of the total charge to accept, negative is none"}

Definition at line 266 of file TrigMultiTrkComboHypo.h.

266 {this, "totalCharge", {},
267 "magnitude of the total charge to accept, negative is none"};

◆ m_outputs

SG::WriteHandleKeyArray<TrigCompositeUtils::DecisionContainer> ComboHypo::m_outputs { this, "HypoOutputDecisions", {}, "Output Decisions" }
privateinherited

Definition at line 51 of file ComboHypo.h.

51{ this, "HypoOutputDecisions", {}, "Output Decisions" };

◆ m_requireUniqueROI

Gaudi::Property<bool> ComboHypo::m_requireUniqueROI
privateinherited
Initial value:
{this, "RequireUniqueROI", false,
"Require each Feature in each leg of the combination to come from a unique L1 seeding ROI."}

Definition at line 53 of file ComboHypo.h.

53 {this, "RequireUniqueROI", false,
54 "Require each Feature in each leg of the combination to come from a unique L1 seeding ROI."};

◆ m_resolvedElectronIDs

TrigCompositeUtils::DecisionIDContainer TrigMultiTrkComboHypo::m_resolvedElectronIDs
private

Definition at line 310 of file TrigMultiTrkComboHypo.h.

◆ m_trackParticleContainerKey

SG::ReadHandleKey<xAOD::TrackParticleContainer> TrigMultiTrkComboHypo::m_trackParticleContainerKey {this, "TrackCollectionKey", "Tracks", "input TrackParticle container name"}
private

Definition at line 256 of file TrigMultiTrkComboHypo.h.

256{this, "TrackCollectionKey", "Tracks", "input TrackParticle container name"};

◆ m_trigBphysContainerKey

SG::WriteHandleKey<xAOD::TrigBphysContainer> TrigMultiTrkComboHypo::m_trigBphysContainerKey {this, "TrigBphysCollectionKey", "TrigBphysContainer", "output TrigBphysContainer name"}
private

Definition at line 259 of file TrigMultiTrkComboHypo.h.

259{this, "TrigBphysCollectionKey", "TrigBphysContainer", "output TrigBphysContainer name"};

◆ m_trigLevel

Gaudi::Property<std::string> TrigMultiTrkComboHypo::m_trigLevel
private
Initial value:
{this, "trigLevel", "EF",
"trigger Level to set for created TrigBphys objects: L2, L2IO, L2MT or EF"}

Definition at line 296 of file TrigMultiTrkComboHypo.h.

296 {this, "trigLevel", "EF",
297 "trigger Level to set for created TrigBphys objects: L2, L2IO, L2MT or EF"};

◆ m_trkMass

Gaudi::Property<std::vector<std::vector<double> > > TrigMultiTrkComboHypo::m_trkMass
private
Initial value:
{this, "trackMasses", {},
"track masses for vertex reco (one per track); muon mass is used by default"}

Definition at line 268 of file TrigMultiTrkComboHypo.h.

268 {this, "trackMasses", {},
269 "track masses for vertex reco (one per track); muon mass is used by default"};

◆ m_trkPt

Gaudi::Property<std::vector<std::vector<double> > > TrigMultiTrkComboHypo::m_trkPt
private
Initial value:
{this, "trackPtThresholds", { {3650., 3650.} },
"minimum track transverse momenta (one per track, sorted descending)"}

Definition at line 270 of file TrigMultiTrkComboHypo.h.

270 {this, "trackPtThresholds", { {3650., 3650.} },
271 "minimum track transverse momenta (one per track, sorted descending)"};

◆ m_trkPtMin

double TrigMultiTrkComboHypo::m_trkPtMin = 0.0
private

Definition at line 313 of file TrigMultiTrkComboHypo.h.

◆ m_useLeptonMomentum

Gaudi::Property<bool> TrigMultiTrkComboHypo::m_useLeptonMomentum
private
Initial value:
{this, "useLeptonMomentum", false,
"use 4-momentum of the xAOD::Muon to make fast calculation of the xAOD::TrigBphys mass, needed for consistency with TrigComboHypoTool::compute()"}

Definition at line 278 of file TrigMultiTrkComboHypo.h.

278 {this, "useLeptonMomentum", false,
279 "use 4-momentum of the xAOD::Muon to make fast calculation of the xAOD::TrigBphys mass, needed for consistency with TrigComboHypoTool::compute()"};

◆ m_v0Tools

ToolHandle<Trk::V0Tools> TrigMultiTrkComboHypo::m_v0Tools {this, "V0Tools", "", "tool to calculate Lxy/LxyError of dimuon candidate wrt beam spot"}
private

Definition at line 305 of file TrigMultiTrkComboHypo.h.

305{this, "V0Tools", "", "tool to calculate Lxy/LxyError of dimuon candidate wrt beam spot"};

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vertexFitter

ToolHandle<Trk::TrkVKalVrtFitter> TrigMultiTrkComboHypo::m_vertexFitter {this, "VertexFitter", "", "VKalVrtFitter tool to fit tracks into the common vertex"}
private

Definition at line 304 of file TrigMultiTrkComboHypo.h.

304{this, "VertexFitter", "", "VKalVrtFitter tool to fit tracks into the common vertex"};

◆ m_vertexPointEstimator

ToolHandle<InDet::VertexPointEstimator> TrigMultiTrkComboHypo::m_vertexPointEstimator {this, "VertexPointEstimator", "", "tool to find starting point for the vertex fitter"}
private

Definition at line 303 of file TrigMultiTrkComboHypo.h.

303{this, "VertexPointEstimator", "", "tool to find starting point for the vertex fitter"};

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.


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