ATLAS Offline Software
Loading...
Searching...
No Matches
PFChargedFlowElementCreatorAlgorithm.cxx
Go to the documentation of this file.
1//eflowRec includes
6
7//EDM includes
11#include "xAODPFlow/PFODefs.h"
13
22
23StatusCode PFChargedFlowElementCreatorAlgorithm::execute(const EventContext& ctx) const {
24
25 ATH_MSG_DEBUG("Starting PFOChargedCreatorAlgorithm::execute");
26
28 ATH_CHECK(chargedFlowElementContainerWriteHandle.record(std::make_unique<xAOD::FlowElementContainer>(),std::make_unique<xAOD::FlowElementAuxContainer>()));
29
30 /* Create Charged FlowElements from all eflowCaloObjects */
32 for (const auto *const thisEflowCaloObject : *eflowCaloObjectContainerReadHandle) createChargedFlowElements(*thisEflowCaloObject,true,chargedFlowElementContainerWriteHandle);
33
34 std::sort(chargedFlowElementContainerWriteHandle->begin(), chargedFlowElementContainerWriteHandle->end(), [] (const xAOD::FlowElement* flowElement1, const xAOD::FlowElement* flowElement2) {return flowElement1->pt()>flowElement2->pt();});
35
36 return StatusCode::SUCCESS;
37
38}
39
40void PFChargedFlowElementCreatorAlgorithm::createChargedFlowElements(const eflowCaloObject& energyFlowCaloObject, bool addClusters, SG::WriteHandle<xAOD::FlowElementContainer>& chargedFlowElementContainerWriteHandle) const {
41
42 /* Loop over all tracks in the eflowCaloObject */
43 for (unsigned int iTrack = 0; iTrack < energyFlowCaloObject.nTracks(); ++iTrack) {
44
45 const eflowRecTrack* efRecTrack = energyFlowCaloObject.efRecTrack(iTrack);
46
47 /* Skip tracks that haven't been subtracted */
48 if (!m_eOverPMode && !efRecTrack->isSubtracted()){continue;}
49
50 /* Create new xAOD::FlowElement and set the type to charged PFlow */
52 chargedFlowElementContainerWriteHandle->push_back(thisFE);
54
55 /* Get the track elementLink and add it to the xAOD:FE.
56 Note we first have to convert it to an IParticle ElementLink. */
58 ElementLink< xAOD::IParticleContainer > theIParticleTrackLink (theTrackLink);
59 std::vector<ElementLink<xAOD::IParticleContainer> > vecIParticleTrackLinkContainer;
60 vecIParticleTrackLinkContainer.push_back(theIParticleTrackLink);
61 thisFE->setChargedObjectLinks(vecIParticleTrackLinkContainer);
62
63 //Now set the charge
64 thisFE->setCharge(efRecTrack->getTrack()->charge());
66
67 std::pair<double,double> etaPhi(0.0,0.0);
68
69 if (m_eOverPMode){
70 /* In EOverPMode want charged eflowObjects to have extrapolated eta,phi as coordinates
71 * (needed for analysis of EOverP Data) */
72 etaPhi = efRecTrack->getTrackCaloPoints().getEM2etaPhi();
73
74 /*add information to xAOD*/
75 const static SG::AuxElement::Accessor<int> accLHED("LayerHED");
76 accLHED(*thisFE) = efRecTrack->getLayerHED();
77
78 const static SG::AuxElement::Accessor<std::vector<int> > accCellOrderVector("LayerVectorCellOrdering");
79 accCellOrderVector(*thisFE) = efRecTrack->getLayerCellOrderVector();
80
81 const static SG::AuxElement::Accessor<std::vector<float> > accRadiusCellOrderVector("RadiusVectorCellOrdering");
82 accRadiusCellOrderVector(*thisFE) = efRecTrack->getRadiusCellOrderVector();
83
84 const static SG::AuxElement::Accessor<std::vector<float> > accAvgEDensityCellOrderVector("AvgEdensityVectorCellOrdering");
85 accAvgEDensityCellOrderVector(*thisFE) = efRecTrack->getAvgEDensityCellOrderVector();
86 } else {
87 /* In normal mode we want the track eta,phi at the perigee */
88 etaPhi.first = efRecTrack->getTrack()->eta();
89 etaPhi.second = efRecTrack->getTrack()->phi();
90 }
91
92 /* Set the 4-vector of the xAOD::PFO */
93 thisFE->setP4(efRecTrack->getTrack()->pt(), etaPhi.first, etaPhi.second, efRecTrack->getTrack()->m());
94
95 ATH_MSG_DEBUG("Created charged PFO with E, pt, eta and phi of " << thisFE->e() << ", " << thisFE->pt() << ", " << thisFE->eta() << " and " << thisFE->phi());
96
97 /* Add the amount of energy the track was expected to deposit in the calorimeter - this is needed to calculate the charged weight in the jet finding */
98 const static SG::AuxElement::Accessor<float> accTracksExpectedEnergyDeposit("TracksExpectedEnergyDeposit");
99 accTracksExpectedEnergyDeposit(*thisFE) = efRecTrack->getEExpect();
100 ATH_MSG_DEBUG("Have set that PFO's expected energy deposit to be " << efRecTrack->getEExpect());
101
102 /* Flag if this track was in a dense environment for later checking */
103 //There is an issue using bools - when written to disk they convert to chars. So lets store the bool as an int.
104 const static SG::AuxElement::Accessor<int> accIsInDenseEnvironment("IsInDenseEnvironment");
105 accIsInDenseEnvironment(*thisFE) = efRecTrack->isInDenseEnvironment();
106
107 /* Optionally we add the links to clusters to the xAOD::PFO */
108 if (addClusters){
109
110 std::vector<std::pair<eflowTrackClusterLink*,std::pair<float,float> > > trackClusterLinkPairs = energyFlowCaloObject.efRecLink();
111
112 ATH_MSG_DEBUG("Have got " << trackClusterLinkPairs.size() << " trackClusterLinkPairs");
113
114 std::vector<eflowTrackClusterLink*> thisTracks_trackClusterLinks = efRecTrack->getClusterMatches();
115
116 ATH_MSG_DEBUG("Have got " << thisTracks_trackClusterLinks.size() << " cluster matches");
117
121
122 std::vector<eflowTrackClusterLink*> thisTracks_trackClusterLinksSubtracted;
123
124 //Create vector of pairs which map each CaloCluster to the asbolute value of its unstracted energy
125 std::vector<std::pair<ElementLink<xAOD::CaloClusterContainer>, double> > vectorClusterToSubtractedEnergies;
126
127 for (auto& trackClusterLink : thisTracks_trackClusterLinks){
128 for (auto& trackClusterLinkPair : trackClusterLinkPairs){
129 //If either of trackClusterLinkPair.second.first or trackClusterLinkPair.second.second is NAN, then no charegd shower
130 //subraction occurred and hence we don't want to save information about that.
131 if (!m_eOverPMode && trackClusterLinkPair.first == trackClusterLink && !std::isnan(trackClusterLinkPair.second.first)) {
132 thisTracks_trackClusterLinksSubtracted.push_back(trackClusterLink);
133 eflowRecCluster* efRecCluster = trackClusterLinkPair.first->getCluster();
134 ElementLink<xAOD::CaloClusterContainer> theOriginalClusterLink = efRecCluster->getOriginalClusElementLink();
135 ElementLink<xAOD::CaloClusterContainer> theSisterClusterLink = (*theOriginalClusterLink)->getSisterClusterLink();
136 ATH_MSG_DEBUG("Will add cluster with E, ratio and absolute subtracted energy " << (*theOriginalClusterLink)->e() << ", " << trackClusterLinkPair.second.first << ", " << trackClusterLinkPair.second.second);
137 if (theSisterClusterLink.isValid()) vectorClusterToSubtractedEnergies.emplace_back(std::pair(theSisterClusterLink,trackClusterLinkPair.second.second));
138 else vectorClusterToSubtractedEnergies.emplace_back(std::pair(theOriginalClusterLink,trackClusterLinkPair.second.second));
139 }
140 else if (m_eOverPMode && trackClusterLinkPair.first == trackClusterLink){
141 thisTracks_trackClusterLinksSubtracted.push_back(trackClusterLink);
142 eflowRecCluster* efRecCluster = trackClusterLinkPair.first->getCluster();
143 ElementLink<xAOD::CaloClusterContainer> theOriginalClusterLink = efRecCluster->getOriginalClusElementLink();
144 ElementLink<xAOD::CaloClusterContainer> theSisterClusterLink = (*theOriginalClusterLink)->getSisterClusterLink();
145 ATH_MSG_DEBUG("Will add cluster with E, ratio and absolute subtracted energy " << (*theOriginalClusterLink)->e() << ", " << 1.0 << ", " << 0.0);
146 if (theSisterClusterLink.isValid()) vectorClusterToSubtractedEnergies.emplace_back(theSisterClusterLink,0.0);
147 else vectorClusterToSubtractedEnergies.emplace_back(theOriginalClusterLink,0.0);
148 }
149 }
150 }
151
152 //sort the vectorClusterToSubtractedEnergies in order of subtracted energy from high (most subtracted) to low (least subtracted)
153 std::sort(vectorClusterToSubtractedEnergies.begin(),vectorClusterToSubtractedEnergies.end(), [](auto const& a, auto const&b){return a.second > b.second;});
154 //now split this into two vectors, ready to be used by the FlowElement
155 std::vector<ElementLink<xAOD::IParticleContainer> > theClusters;
156 std::vector<float> theClusterWeights;
157 for (auto thePair : vectorClusterToSubtractedEnergies){
158 ElementLink< xAOD::IParticleContainer > theIParticleTrackLink(thePair.first);
159 theClusters.push_back(theIParticleTrackLink);
160 theClusterWeights.push_back(thePair.second);
161 }
162
163 //Set the vector of clusters and vector of corresponding energies.
164 thisFE->setOtherObjectLinks(theClusters,theClusterWeights);
165
166 }//if we add clusters
167
168 //Add detailed CP data
169 if (m_addCPData){
170 const static SG::AuxElement::Accessor<float> accEtaEM2("EtaEM2");
171 accEtaEM2(*thisFE) = efRecTrack->getTrackCaloPoints().getEM2etaPhi().first;
172
173 const static SG::AuxElement::Accessor<float> accPhiEM2("PhiEM2");
174 accPhiEM2(*thisFE) = efRecTrack->getTrackCaloPoints().getEM2etaPhi().second;
175
176 const static SG::AuxElement::Accessor<char> accIsRecovered("isRecovered");
177 accIsRecovered(*thisFE) = efRecTrack->isRecovered();
178
179 const static SG::AuxElement::Accessor<unsigned int> accNumMatchedClusters("numMatchedClusters");
180 accNumMatchedClusters(*thisFE) = efRecTrack->getClusterMatches().size();
181
182 const static SG::AuxElement::Accessor<std::vector<float> > accDRPrimes("dRPrimes");
183 accDRPrimes(*thisFE) = efRecTrack->getDRPrimes();
184
185 const static SG::AuxElement::Accessor<float > accPull15("Pull15");
186 accPull15(*thisFE) = efRecTrack->getpull15();
187
188 const static SG::AuxElement::Accessor<std::vector<std::pair<ElementLink<CaloCellContainer>, double> > > accSubtractedCaloCells("SubtractedCaloCells");
189 accSubtractedCaloCells(*thisFE) = efRecTrack->getSubtractedCaloCells();
190 }
191
192 }//loop over eflowRecTracks
193
194}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
static Double_t a
SG::WriteHandleKey< xAOD::FlowElementContainer > m_chargedFlowElementContainerWriteHandleKey
WriteHandleKey for charged PFO.
SG::ReadHandleKey< eflowCaloObjectContainer > m_eflowCaloObjectContainerReadHandleKey
ReadHandleKey for eflowCaloObjectContainer.
Gaudi::Property< bool > m_eOverPMode
Toggle EOverP algorithm mode, whereby no charged shower subtraction is performed.
void createChargedFlowElements(const eflowCaloObject &energyFlowCaloObject, bool addClusters, SG::WriteHandle< xAOD::FlowElementContainer > &chargedFlowElementContainerWriteHandle) const
Create the charged PFO.
Gaudi::Property< bool > m_addCPData
Toggle whether to decorate FlowElements with addutional data for Combined Performance studies.
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
An internal EDM object which stores information about systems of associated tracks and calorimeter cl...
const eflowRecTrack * efRecTrack(int i) const
const std::vector< std::pair< eflowTrackClusterLink *, std::pair< float, float > > > & efRecLink() const
unsigned nTracks() const
This class extends the information about a xAOD::CaloCluster.
ElementLink< xAOD::CaloClusterContainer > getOriginalClusElementLink() const
xAOD::CaloCluster * getCluster()
This class extends the information about a xAOD::Track.
double getEExpect() const
const std::vector< float > & getDRPrimes() const
const std::vector< eflowTrackClusterLink * > & getClusterMatches() const
const std::vector< float > & getRadiusCellOrderVector() const
bool isSubtracted() const
double getpull15() const
const xAOD::TrackParticle * getTrack() const
const eflowTrackCaloPoints & getTrackCaloPoints() const
const std::vector< int > & getLayerCellOrderVector() const
const std::vector< std::pair< ElementLink< CaloCellContainer >, double > > & getSubtractedCaloCells() const
bool isInDenseEnvironment() const
int getLayerHED() const
bool isRecovered() const
const std::vector< float > & getAvgEDensityCellOrderVector() const
ElementLink< xAOD::TrackParticleContainer > getTrackElemLink() const
const std::pair< float, float > getEM2etaPhi() const
virtual double pt() const override
void setP4(float pt, float eta, float phi, float m)
virtual double phi() const override
The azimuthal angle ( ) of the particle.
void setCharge(float c)
virtual double eta() const override
The pseudorapidity ( ) of the particle.
void setChargedObjectLinks(const std::vector< ElementLink< IParticleContainer > > &elV)
void setSignalType(signal_t t)
void setOtherObjectLinks(const std::vector< ElementLink< IParticleContainer > > &elV)
virtual double e() const override
The total energy of the particle.
virtual double m() const override final
The invariant mass of the particle..
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
float charge() const
Returns the charge.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition FlowElement.h:16