ATLAS Offline Software
Loading...
Searching...
No Matches
TrigEgammaTLAPhotonReAlgo.cxx
Go to the documentation of this file.
1 /*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4
11
12
13class ISvcLocator;
14 // helpers functions, if needed
16
17 bool operator()(const xAOD::Photon* l, const xAOD::Photon* r) const {
18 return l->p4().Et() > r->p4().Et();
19 }
20
21};
22
24
27
28 bool operator()(const xAOD::Photon* myphoton) const {
29 return myphoton->p4().Pt() > thresholdPt;
30 }
31
32 };
33
34
35TrigEgammaTLAPhotonReAlgo::TrigEgammaTLAPhotonReAlgo(const std::string & name, ISvcLocator* pSvcLocator)
36 : AthReentrantAlgorithm (name, pSvcLocator)
37 {
38
39 }
40
42{
43
44 // here we simply check that the class properties have been correctly initialized at the config stage
45 ATH_MSG_DEBUG( "Initializing " << name() << " ...") ;
46 ATH_MSG_DEBUG( "pT threshold for saving Photons to TLA Collection: " << m_photonPtThreshold) ;
47 ATH_MSG_DEBUG( "Maximum number of Photons to be saved: " << m_maxNPhotons );
48
49 // check if the data handles have been initialized
50 ATH_CHECK(m_inputPhotonsKeys.initialize());
51 ATH_CHECK(m_TLAOutPhotonsKey.initialize());
52
53
54 if (static_cast<int>(m_maxNPhotons) == 0){ // allow m_maxNPhotons == -1 as flag to save all Photons
55 ATH_MSG_ERROR("This algorithm is initialized to save NO (0) Photons, so it shouldn't be running at all!");
56 return StatusCode::FAILURE;
57 }
58
59
60 return StatusCode::SUCCESS;
61}
62
63StatusCode TrigEgammaTLAPhotonReAlgo::execute(const EventContext& ctx) const
64{
65 using namespace xAOD;
66
67 ATH_MSG_DEBUG("Executing " << name() << " ...");
68
71
72 // effectively make the TLA Photon Container
73
74 ATH_MSG_DEBUG("Retrieving Photons from " << h_inputPhotons.key() );
75 ATH_MSG_DEBUG("Placing <selected> Photons in " << h_TLAPhotons.key() );
76
77 //ATH_CHECK(h_TLAPhotons.record (std::make_unique<xAOD::PhotonContainer>(),
78 // std::make_unique<PhotonAuxContainer>()) );
79
80
81 const xAOD::PhotonContainer* inputPhotons = h_inputPhotons.get();
82 std::vector<const xAOD::Photon*> originalPhotons(inputPhotons->begin(), inputPhotons->end());
83
84 // define the maximum number of photons we care about: either equivalent to m_maxNPhotons if smaller than size of vector, or keep all photons (in case of negative value)
85 std::vector<const xAOD::Photon*>::iterator it_maxPhotonBound;
86
87 int maxNPhotons = static_cast<int>(m_maxNPhotons);
88 int sizeOfOriginalPhotonContainer = static_cast<int>(originalPhotons.size());
89
90 // take the largest between the number of photons requested and the size of the vector, in iterator form
91 if (m_maxNPhotons > 0) it_maxPhotonBound = maxNPhotons < sizeOfOriginalPhotonContainer ? originalPhotons.begin() + maxNPhotons : originalPhotons.end();
92 else it_maxPhotonBound = originalPhotons.end();
93
94 // check the sort order of the input container is ok
95 // use a partial sort to save some time
96 std::partial_sort (originalPhotons.begin(), it_maxPhotonBound, originalPhotons.end(), DescendingEt());
97
98 // get an iterator to the last element above the pT threshold (because we ordered the photons, this is the last one we want)
99 std::vector<const xAOD::Photon*>::iterator it_ptThresholdBound;
100 it_ptThresholdBound = std::partition(originalPhotons.begin(), it_maxPhotonBound, HasPtAboveThreshold(static_cast<float>(m_photonPtThreshold)));
101
102 //make the output photon container
103 ATH_CHECK(h_TLAPhotons.record (std::make_unique<xAOD::PhotonContainer>(),
104 std::make_unique<PhotonAuxContainer>()) ); // in FastPhotonFex this was TrigEMClusterAuxContainer
105
106
107 //loop on all the photons from the beginning to the last photon we want, and put them the output photon collection
108 //also, set a link to their parent photon
109
110
111 for( auto it_ph=originalPhotons.begin(); it_ph!=it_ptThresholdBound; ++it_ph ) {
112 xAOD::Photon* copiedPhoton = new xAOD::Photon(*(*it_ph));
113
114 ATH_CHECK(setOriginalObjectLink(*(*it_ph),*copiedPhoton));
115
116 h_TLAPhotons->push_back(copiedPhoton);
117 ATH_MSG_DEBUG("Selected photon pT: " << copiedPhoton->pt());
118
119
120
121 }
122
123 return StatusCode::SUCCESS;
124}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
An algorithm that can be simultaneously executed in multiple threads.
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.
const_pointer_type get() const
Dereference the pointer, but don't cache anything.
virtual const std::string & key() const override final
Return the StoreGate ID for the referenced object.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
SG::WriteHandleKey< xAOD::PhotonContainer > m_TLAOutPhotonsKey
virtual StatusCode initialize() override
Gaudi::Property< float > m_photonPtThreshold
virtual StatusCode execute(const EventContext &ctx) const override
TrigEgammaTLAPhotonReAlgo(const std::string &name, ISvcLocator *pSvcLocator)
SG::ReadHandleKey< xAOD::PhotonContainer > m_inputPhotonsKeys
Gaudi::Property< int > m_maxNPhotons
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Definition Egamma_v1.cxx:66
virtual FourMom_t p4() const override final
The full 4-momentum of the particle as a TLoretzVector.
Definition Egamma_v1.cxx:94
int r
Definition globals.cxx:22
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
DataModel_detail::iterator< DVL > partition(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of partition for DataVector/List.
void partial_sort(DataModel_detail::iterator< DVL > beg, DataModel_detail::iterator< DVL > mid, DataModel_detail::iterator< DVL > end)
Specialization of partial_sort for DataVector/List.
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
PhotonContainer_v1 PhotonContainer
Definition of the current "photon container version".
bool setOriginalObjectLink(const IParticle &original, IParticle &copy)
This function should be used by CP tools when they make a deep copy of an object in their correctedCo...
Photon_v1 Photon
Definition of the current "egamma version".
bool operator()(const xAOD::Photon *l, const xAOD::Photon *r) const
HasPtAboveThreshold(double thresholdPt)
bool operator()(const xAOD::Photon *myphoton) const