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

Update the input RoI's z-coordinate and z-width to the HitZ regression estimation, if the regression sigma is below the indicated threshold. More...

#include <TrigTauHitZRoiUpdater.h>

Inheritance diagram for TrigTauHitZRoiUpdater:
Collaboration diagram for TrigTauHitZRoiUpdater:

Public Member Functions

 TrigTauHitZRoiUpdater (const std::string &, ISvcLocator *)
virtual StatusCode initialize () override
virtual StatusCode execute (const EventContext &) const override
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual bool isClonable () const override
 Specify if the algorithm is clonable.
virtual unsigned int cardinality () const override
 Cardinality (Maximum number of clones that can exist) special value 0 means that algorithm is reentrant.
virtual StatusCode sysExecute (const EventContext &ctx) override
 Execute an algorithm.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
virtual bool filterPassed (const EventContext &ctx) const
virtual void setFilterPassed (bool state, const EventContext &ctx) const
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

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

Private Attributes

Gaudi::Property< float > m_z0HalfWidth {this, "z0HalfWidth", 10.0, "z0 half width for the output RoI, if the HitZ regression is successful, in mm"}
Gaudi::Property< float > m_etaHalfWidth {this, "etaHalfWidth", 0.1, "eta half width for the output RoI, if the HitZ regression is successful"}
Gaudi::Property< float > m_phiHalfWidth {this, "phiHalfWidth", 0.1, "phi half width for the output RoI, if the HitZ regression is successful"}
Gaudi::Property< float > m_maxPt {this, "maxPt", 1e5, "Maximum Pt of the taus to apply the HitZ regression RoI update, in GeV"}
Gaudi::Property< float > m_maxSigma {this, "maxSigma", 5.0, "Maximum HitZ regression sigma to update the RoI, in mm"}
SG::ReadHandleKey< TrigRoiDescriptorCollectionm_roIInputKey {this, "RoIInputKey", "", "Input RoI key"}
SG::WriteHandleKey< TrigRoiDescriptorCollectionm_roIOutputKey {this, "RoIOutputKey", "", "Output RoI key"}
SG::ReadHandleKey< xAOD::TauJetContainerm_tauKey {this, "TauKey", "", "Input Tau container key"}
SG::ReadDecorHandleKey< xAOD::TauJetContainerm_zDecorKey {this, "zDecorKey", "", "HitZ regression z decoration key"}
SG::ReadDecorHandleKey< xAOD::TauJetContainerm_sigmaDecorKey {this, "sigmaDecorKey", "", "HitZ regression sigma decoration key"}
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

Update the input RoI's z-coordinate and z-width to the HitZ regression estimation, if the regression sigma is below the indicated threshold.

Definition at line 21 of file TrigTauHitZRoiUpdater.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

◆ TrigTauHitZRoiUpdater()

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

Definition at line 18 of file TrigTauHitZRoiUpdater.cxx.

19 : AthReentrantAlgorithm(name, pSvcLocator)
20{
21
22}

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}

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

Definition at line 53 of file TrigTauHitZRoiUpdater.cxx.

54{
55 ATH_MSG_DEBUG("Running " << name());
56
57 //---------------------------------------------------------------
58 // Prepare I/O
59 //---------------------------------------------------------------
60
61 // Prepare output RoI container
62 std::unique_ptr<TrigRoiDescriptorCollection> roiCollection = std::make_unique<TrigRoiDescriptorCollection>();
63 SG::WriteHandle<TrigRoiDescriptorCollection> outputRoIHandle(m_roIOutputKey, ctx);
64 ATH_CHECK(outputRoIHandle.record(std::move(roiCollection)));
65
66
67 // Retrieve Input Tau container
68 SG::ReadHandle<xAOD::TauJetContainer> tauHandle(m_tauKey, ctx);
69 ATH_CHECK(tauHandle.isValid());
70 const xAOD::TauJetContainer* tauContainer = tauHandle.get();
71 if(!tauContainer) {
72 ATH_MSG_ERROR("No tau container found, the Tau RoI updater should not be scheduled");
73 return StatusCode::FAILURE;
74 }
75 ATH_MSG_DEBUG("Found " << tauContainer->size() << " taus, updating the RoI");
76 const xAOD::TauJet *tau = tauContainer->at(0); // We only have one tau in the container
77
78
79 // Retrieve input RoI descriptor
80 SG::ReadHandle<TrigRoiDescriptorCollection> roisHandle(m_roIInputKey, ctx);
81 ATH_CHECK(roisHandle.isValid());
82 ATH_MSG_DEBUG("Size of roisHandle: " << roisHandle->size());
83 if(roisHandle->size() != 1) {
84 ATH_MSG_ERROR("Expected exactly one RoI!");
85 return StatusCode::FAILURE;
86 }
87 const TrigRoiDescriptor* roiDescriptor = roisHandle->at(0); // We only have one RoI in the handle
88
89
90 // Fill local variables for RoI reference position
91 const float eta = roiDescriptor->eta();
92 const float phi = roiDescriptor->phi();
93 float zed = roiDescriptor->zed();
94 float zedMinus = roiDescriptor->zedMinus();
95 float zedPlus = roiDescriptor->zedPlus();
96
97
98 //---------------------------------------------------------------
99 // Get HitZ regression results
100 //---------------------------------------------------------------
101
102 SG::ReadDecorHandle<xAOD::TauJetContainer, float> zDecor(m_zDecorKey, ctx);
103 ATH_CHECK(zDecor.isAvailable());
104 const float hitz_z0 = zDecor(*tau);
105
106 SG::ReadDecorHandle<xAOD::TauJetContainer, float> sigmaDecor(m_sigmaDecorKey, ctx);
107 ATH_CHECK(sigmaDecor.isAvailable());
108 const float hitz_sigma = sigmaDecor(*tau);
109
110 ATH_MSG_DEBUG("HitZ Tau pt: " << tau->pt() << ", eta: " << tau->eta() << ", phi: " << tau->phi() << ", HitZ z0: " << hitz_z0 << ", HitZ sigma: " << hitz_sigma);
111
112 //---------------------------------------------------------------
113 // Update the RoI
114 //---------------------------------------------------------------
115 // If the HitZ regression is valid (z0 != -1111) and the regression sigma is below the threshold,
116 // update the RoI z position and width:
117 constexpr float invalid_z0 = -1111;
118 if(hitz_z0 != invalid_z0 && tau->pt() < m_maxPt * Gaudi::Units::GeV && hitz_sigma < m_maxSigma) {
119 zed = hitz_z0;
120 zedMinus = zed - m_z0HalfWidth;
121 zedPlus = zed + m_z0HalfWidth;
122 }
123
124 const float etaMinus = eta - m_etaHalfWidth;
125 const float etaPlus = eta + m_etaHalfWidth;
126 const float phiMinus = CxxUtils::wrapToPi(phi - m_phiHalfWidth);
127 const float phiPlus = CxxUtils::wrapToPi(phi + m_phiHalfWidth);
128
129 // Create the new RoI
130 outputRoIHandle->push_back(std::make_unique<TrigRoiDescriptor>(
131 roiDescriptor->roiWord(), roiDescriptor->l1Id(), roiDescriptor->roiId(),
132 eta, etaMinus, etaPlus,
133 phi, phiMinus, phiPlus,
134 zed, zedMinus, zedPlus
135 ));
136
137
138 ATH_MSG_DEBUG("Input RoI: " << *roiDescriptor);
139 ATH_MSG_DEBUG("Output RoI: " << *outputRoIHandle->back());
140
141 return StatusCode::SUCCESS;
142}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
Athena::TPCnvVers::Current TrigRoiDescriptor
const T * get(size_type n) const
Access an element, as an rvalue.
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.
virtual double zed() const override final
virtual double phi() const override final
Methods to retrieve data members.
virtual double zedPlus() const override final
z at the most forward end of the RoI
virtual double zedMinus() const override final
z at the most backward end of the RoI
virtual double eta() const override final
virtual unsigned int roiWord() const override final
virtual unsigned int roiId() const override final
these quantities probably don't need to be used any more
virtual unsigned int l1Id() const override final
Gaudi::Property< float > m_maxSigma
SG::ReadDecorHandleKey< xAOD::TauJetContainer > m_zDecorKey
Gaudi::Property< float > m_phiHalfWidth
SG::ReadHandleKey< xAOD::TauJetContainer > m_tauKey
Gaudi::Property< float > m_etaHalfWidth
SG::WriteHandleKey< TrigRoiDescriptorCollection > m_roIOutputKey
Gaudi::Property< float > m_maxPt
Gaudi::Property< float > m_z0HalfWidth
SG::ReadHandleKey< TrigRoiDescriptorCollection > m_roIInputKey
SG::ReadDecorHandleKey< xAOD::TauJetContainer > m_sigmaDecorKey
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double eta() const
The pseudorapidity ( ) of the particle.
T wrapToPi(T phi)
Wrap angle in radians to [-pi, pi].
Definition phihelper.h:24
TauJet_v3 TauJet
Definition of the current "tau version".
TauJetContainer_v3 TauJetContainer
Definition of the current "taujet container version".

◆ extraDeps_update_handler()

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

Add StoreName to extra input/output deps as needed.

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

◆ extraOutputDeps()

const DataObjIDColl & AthCommonReentrantAlgorithm< Gaudi::Algorithm >::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

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

Definition at line 94 of file AthCommonReentrantAlgorithm.cxx.

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

◆ filterPassed()

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

Definition at line 96 of file AthCommonReentrantAlgorithm.h.

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

◆ initialize()

StatusCode TrigTauHitZRoiUpdater::initialize ( )
overridevirtual

Definition at line 25 of file TrigTauHitZRoiUpdater.cxx.

26{
27 ATH_MSG_DEBUG("Initializing " << name() );
28 ATH_MSG_DEBUG("z0HalfWidth: " << m_z0HalfWidth.value() );
29 ATH_MSG_DEBUG("etaHalfWidth: " << m_etaHalfWidth.value() );
30 ATH_MSG_DEBUG("phiHalfWidth: " << m_phiHalfWidth.value() );
31 ATH_MSG_DEBUG("maxSigma: " << m_maxSigma.value() );
32 ATH_MSG_DEBUG("maxPt: " << m_maxPt.value() );
33
34 if(m_z0HalfWidth < 0 || m_etaHalfWidth < 0 || m_phiHalfWidth < 0) {
35 ATH_MSG_ERROR("Incorrect parameters");
36 return StatusCode::FAILURE;
37 }
38
39 ATH_MSG_DEBUG("Initialising HandleKeys");
40 ATH_CHECK(m_roIInputKey.initialize());
41 ATH_CHECK(m_roIOutputKey.initialize());
42
43 ATH_CHECK(m_tauKey.initialize());
44 m_zDecorKey = m_tauKey.key() + "." + m_zDecorKey.key();
45 ATH_CHECK(m_zDecorKey.initialize());
46 m_sigmaDecorKey = m_tauKey.key() + "." + m_sigmaDecorKey.key();
47 ATH_CHECK(m_sigmaDecorKey.initialize());
48
49 return StatusCode::SUCCESS;
50}

◆ inputHandles()

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

Return this algorithm's input handles.

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

◆ isClonable()

◆ msg()

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

Definition at line 24 of file AthCommonMsg.h.

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

◆ msgLvl()

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

Definition at line 30 of file AthCommonMsg.h.

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

◆ outputHandles()

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

Return this algorithm's output handles.

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

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

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

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ setFilterPassed()

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

Definition at line 100 of file AthCommonReentrantAlgorithm.h.

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

◆ sysExecute()

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

Execute an algorithm.

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

Definition at line 85 of file AthCommonReentrantAlgorithm.cxx.

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

◆ sysInitialize()

StatusCode AthCommonReentrantAlgorithm< Gaudi::Algorithm >::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

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

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

Reimplemented from AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >.

Reimplemented in HypoBase, and InputMakerBase.

Definition at line 61 of file AthCommonReentrantAlgorithm.cxx.

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

◆ sysStart()

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

Handle START transition.

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

◆ 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_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_etaHalfWidth

Gaudi::Property<float> TrigTauHitZRoiUpdater::m_etaHalfWidth {this, "etaHalfWidth", 0.1, "eta half width for the output RoI, if the HitZ regression is successful"}
private

Definition at line 31 of file TrigTauHitZRoiUpdater.h.

31{this, "etaHalfWidth", 0.1, "eta half width for the output RoI, if the HitZ regression is successful"};

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

Gaudi::Property<float> TrigTauHitZRoiUpdater::m_maxPt {this, "maxPt", 1e5, "Maximum Pt of the taus to apply the HitZ regression RoI update, in GeV"}
private

Definition at line 34 of file TrigTauHitZRoiUpdater.h.

34{this, "maxPt", 1e5, "Maximum Pt of the taus to apply the HitZ regression RoI update, in GeV"};

◆ m_maxSigma

Gaudi::Property<float> TrigTauHitZRoiUpdater::m_maxSigma {this, "maxSigma", 5.0, "Maximum HitZ regression sigma to update the RoI, in mm"}
private

Definition at line 35 of file TrigTauHitZRoiUpdater.h.

35{this, "maxSigma", 5.0, "Maximum HitZ regression sigma to update the RoI, in mm"};

◆ m_phiHalfWidth

Gaudi::Property<float> TrigTauHitZRoiUpdater::m_phiHalfWidth {this, "phiHalfWidth", 0.1, "phi half width for the output RoI, if the HitZ regression is successful"}
private

Definition at line 32 of file TrigTauHitZRoiUpdater.h.

32{this, "phiHalfWidth", 0.1, "phi half width for the output RoI, if the HitZ regression is successful"};

◆ m_roIInputKey

SG::ReadHandleKey<TrigRoiDescriptorCollection> TrigTauHitZRoiUpdater::m_roIInputKey {this, "RoIInputKey", "", "Input RoI key"}
private

Definition at line 37 of file TrigTauHitZRoiUpdater.h.

37{this, "RoIInputKey", "", "Input RoI key"};

◆ m_roIOutputKey

SG::WriteHandleKey<TrigRoiDescriptorCollection> TrigTauHitZRoiUpdater::m_roIOutputKey {this, "RoIOutputKey", "", "Output RoI key"}
private

Definition at line 38 of file TrigTauHitZRoiUpdater.h.

38{this, "RoIOutputKey", "", "Output RoI key"};

◆ m_sigmaDecorKey

SG::ReadDecorHandleKey<xAOD::TauJetContainer> TrigTauHitZRoiUpdater::m_sigmaDecorKey {this, "sigmaDecorKey", "", "HitZ regression sigma decoration key"}
private

Definition at line 42 of file TrigTauHitZRoiUpdater.h.

42{this, "sigmaDecorKey", "", "HitZ regression sigma decoration key"};

◆ m_tauKey

SG::ReadHandleKey<xAOD::TauJetContainer> TrigTauHitZRoiUpdater::m_tauKey {this, "TauKey", "", "Input Tau container key"}
private

Definition at line 40 of file TrigTauHitZRoiUpdater.h.

40{this, "TauKey", "", "Input Tau container key"};

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.

◆ m_z0HalfWidth

Gaudi::Property<float> TrigTauHitZRoiUpdater::m_z0HalfWidth {this, "z0HalfWidth", 10.0, "z0 half width for the output RoI, if the HitZ regression is successful, in mm"}
private

Definition at line 30 of file TrigTauHitZRoiUpdater.h.

30{this, "z0HalfWidth", 10.0, "z0 half width for the output RoI, if the HitZ regression is successful, in mm"};

◆ m_zDecorKey

SG::ReadDecorHandleKey<xAOD::TauJetContainer> TrigTauHitZRoiUpdater::m_zDecorKey {this, "zDecorKey", "", "HitZ regression z decoration key"}
private

Definition at line 41 of file TrigTauHitZRoiUpdater.h.

41{this, "zDecorKey", "", "HitZ regression z decoration key"};

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