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

An internal EDM object which stores information about systems of associated tracks and calorimeter clusters. More...

#include <eflowCaloObject.h>

Collaboration diagram for eflowCaloObject:

Public Member Functions

 eflowCaloObject ()=default
 ~eflowCaloObject ()
void addTrack (eflowRecTrack *track)
void addCluster (eflowRecCluster *cluster)
void addTrackClusterLinks (const std::vector< eflowTrackClusterLink * > &trackClusterLink)
void addTracks (const std::vector< eflowRecTrack * > &tracks)
void addClusters (const std::vector< eflowRecCluster * > &clusters)
void setTrackClusterLinkSubtractionStatus (unsigned int index, const std::pair< float, float > &energyRatio_energyValPair)
const eflowRecTrackefRecTrack (int i) const
eflowRecTrackefRecTrack (int i)
unsigned nTracks () const
void clearTracks ()
const eflowRecClusterefRecCluster (int i) const
eflowRecClusterefRecCluster (int i)
unsigned nClusters () const
void clearClusters ()
const std::vector< std::pair< eflowTrackClusterLink *, std::pair< float, float > > > & efRecLink () const
void clearLinks ()
double getExpectedEnergy () const
double getExpectedVariance () const
double getClusterEnergy () const
void simulateShower (eflowLayerIntegrator *integrator, const eflowEEtaBinnedParameters *binnedParameters, const PFEnergyPredictorTool *energyP, bool useLegacyEnergyBinIndexing)

Private Member Functions

void addTrackClusterLink (eflowTrackClusterLink *trackClusterLink)

Private Attributes

std::vector< eflowRecCluster * > m_eflowRecClusters
std::vector< std::pair< eflowTrackClusterLink *, std::pair< float, float > > > m_trackClusterLinks
std::vector< eflowRecTrack * > m_eflowRecTracks

Detailed Description

An internal EDM object which stores information about systems of associated tracks and calorimeter clusters.

Specifically it stores vectors of pointers to eflowRecTracks, eflowRecClusters and eflowTrackClusterLinks. In addition it stores links to an xAOD::CaloClusterContainer and its associated aux container. This class also calculates the expected energy deposit in the calorimeter from a track in the system, and stores that information so that clients can retrieve it. It also calculates the calorimeter cell ordering to be used in the subtraction. Both of these things are done in the simulateShower method which uses the data stored in an eflowEEtaBinnedParameters object, which is filled by e.g the eflowCellEOverP_mc12_JetETMiss tool.

Definition at line 34 of file eflowCaloObject.h.

Constructor & Destructor Documentation

◆ eflowCaloObject()

eflowCaloObject::eflowCaloObject ( )
default

◆ ~eflowCaloObject()

eflowCaloObject::~eflowCaloObject ( )
default

Member Function Documentation

◆ addCluster()

void eflowCaloObject::addCluster ( eflowRecCluster * cluster)
inline

Definition at line 41 of file eflowCaloObject.h.

41{ m_eflowRecClusters.push_back(cluster); }
std::vector< eflowRecCluster * > m_eflowRecClusters

◆ addClusters()

void eflowCaloObject::addClusters ( const std::vector< eflowRecCluster * > & clusters)

Definition at line 36 of file eflowCaloObject.cxx.

36 {
37 m_eflowRecClusters.insert(m_eflowRecClusters.end(), clusters.begin(), clusters.end());
38 }

◆ addTrack()

void eflowCaloObject::addTrack ( eflowRecTrack * track)
inline

Definition at line 40 of file eflowCaloObject.h.

40{ m_eflowRecTracks.push_back(track); }
std::vector< eflowRecTrack * > m_eflowRecTracks

◆ addTrackClusterLink()

void eflowCaloObject::addTrackClusterLink ( eflowTrackClusterLink * trackClusterLink)
inlineprivate

Definition at line 78 of file eflowCaloObject.h.

78{ m_trackClusterLinks.push_back(std::pair(trackClusterLink,std::pair(NAN,NAN))); }
std::vector< std::pair< eflowTrackClusterLink *, std::pair< float, float > > > m_trackClusterLinks

◆ addTrackClusterLinks()

void eflowCaloObject::addTrackClusterLinks ( const std::vector< eflowTrackClusterLink * > & trackClusterLink)

Definition at line 26 of file eflowCaloObject.cxx.

26 {
27 for (auto *ptr : trackClusterLink) {
29 }
30}
void addTrackClusterLink(eflowTrackClusterLink *trackClusterLink)

◆ addTracks()

void eflowCaloObject::addTracks ( const std::vector< eflowRecTrack * > & tracks)

Definition at line 32 of file eflowCaloObject.cxx.

32 {
33 m_eflowRecTracks.insert(m_eflowRecTracks.end(), tracks.begin(), tracks.end());
34 }

◆ clearClusters()

void eflowCaloObject::clearClusters ( )
inline

Definition at line 61 of file eflowCaloObject.h.

61{ m_eflowRecClusters.clear(); }

◆ clearLinks()

void eflowCaloObject::clearLinks ( )
inline

Definition at line 66 of file eflowCaloObject.h.

66{ m_trackClusterLinks.clear(); }

◆ clearTracks()

void eflowCaloObject::clearTracks ( )
inline

Definition at line 55 of file eflowCaloObject.h.

55{ m_eflowRecTracks.clear(); }

◆ efRecCluster() [1/2]

eflowRecCluster * eflowCaloObject::efRecCluster ( int i)
inline

Definition at line 59 of file eflowCaloObject.h.

59{ return m_eflowRecClusters[i]; }

◆ efRecCluster() [2/2]

const eflowRecCluster * eflowCaloObject::efRecCluster ( int i) const
inline

Definition at line 58 of file eflowCaloObject.h.

58{ return m_eflowRecClusters[i]; }

◆ efRecLink()

const std::vector< std::pair< eflowTrackClusterLink *, std::pair< float, float > > > & eflowCaloObject::efRecLink ( ) const
inline

Definition at line 65 of file eflowCaloObject.h.

65{ return m_trackClusterLinks; }

◆ efRecTrack() [1/2]

eflowRecTrack * eflowCaloObject::efRecTrack ( int i)
inline

Definition at line 53 of file eflowCaloObject.h.

53{ return m_eflowRecTracks[i]; }

◆ efRecTrack() [2/2]

const eflowRecTrack * eflowCaloObject::efRecTrack ( int i) const
inline

Definition at line 52 of file eflowCaloObject.h.

52{ return m_eflowRecTracks[i]; }

◆ getClusterEnergy()

double eflowCaloObject::getClusterEnergy ( ) const

Definition at line 56 of file eflowCaloObject.cxx.

56 {
57 double clusterEnergy(0.0);
58 for(unsigned iCluster=0; iCluster < nClusters(); ++iCluster) {
59 clusterEnergy += efRecCluster(iCluster)->getCluster()->e();
60 }
61 return clusterEnergy;
62}
unsigned nClusters() const
const eflowRecCluster * efRecCluster(int i) const
xAOD::CaloCluster * getCluster()
virtual double e() const
The total energy of the particle.

◆ getExpectedEnergy()

double eflowCaloObject::getExpectedEnergy ( ) const

Definition at line 40 of file eflowCaloObject.cxx.

40 {
41 double expectedEnergy(0.0);
42 for(unsigned iTrack=0; iTrack < nTracks(); ++iTrack) {
43 if (!efRecTrack(iTrack)->isInDenseEnvironment()) expectedEnergy += efRecTrack(iTrack)->getEExpect();
44 }
45 return expectedEnergy;
46}
const eflowRecTrack * efRecTrack(int i) const
unsigned nTracks() const
double getEExpect() const

◆ getExpectedVariance()

double eflowCaloObject::getExpectedVariance ( ) const

Definition at line 48 of file eflowCaloObject.cxx.

48 {
49 double expectedVariance(0.0);
50 for(unsigned iTrack=0; iTrack < nTracks(); ++iTrack) {
51 if (!efRecTrack(iTrack)->isInDenseEnvironment()) expectedVariance += efRecTrack(iTrack)->getVarEExpect();
52 }
53 return expectedVariance;
54}
double getVarEExpect() const

◆ nClusters()

unsigned eflowCaloObject::nClusters ( ) const
inline

Definition at line 60 of file eflowCaloObject.h.

60{ return m_eflowRecClusters.size(); }

◆ nTracks()

unsigned eflowCaloObject::nTracks ( ) const
inline

Definition at line 54 of file eflowCaloObject.h.

54{ return m_eflowRecTracks.size(); }

◆ setTrackClusterLinkSubtractionStatus()

void eflowCaloObject::setTrackClusterLinkSubtractionStatus ( unsigned int index,
const std::pair< float, float > & energyRatio_energyValPair )
inline

Definition at line 49 of file eflowCaloObject.h.

49{ m_trackClusterLinks[index].second = energyRatio_energyValPair; }
str index
Definition DeMoScan.py:362

◆ simulateShower()

void eflowCaloObject::simulateShower ( eflowLayerIntegrator * integrator,
const eflowEEtaBinnedParameters * binnedParameters,
const PFEnergyPredictorTool * energyP,
bool useLegacyEnergyBinIndexing )

Definition at line 64 of file eflowCaloObject.cxx.

65 {
66
67 for (auto *thisEfRecTrack : m_eflowRecTracks) {
68
69 std::vector<eflowRecCluster*> matchedClusters;
70 const std::vector<eflowTrackClusterLink*>& links = thisEfRecTrack->getClusterMatches();
71 for (auto* itLink : links) {
72 matchedClusters.push_back(itLink->getCluster());
73 }
74
75 double trackEM1eta = thisEfRecTrack->getTrackCaloPoints().getEM1eta();
76 /* if a track is in the forward EM (2.5 < eta < 3.2) then there is no EM1 -> need to use EM2 */
77 if(trackEM1eta<-998.) trackEM1eta = thisEfRecTrack->getTrackCaloPoints().getEM2eta();
78 /* if a track is not in the EM region (3.2 < eta < 4.0) then should use FCAL0 */
79 if(trackEM1eta<-998.) trackEM1eta = thisEfRecTrack->getTrackCaloPoints().getFCAL0eta();
80
81 double trackE = thisEfRecTrack->getTrack()->e();
82
83 if (!binnedParameters->binExists(trackE, trackEM1eta)) {
84 thisEfRecTrack->setHasBin(false);
85 return;
86 }
87
88 /* Determine the LFI */
89 integrator->measureNewClus(matchedClusters, thisEfRecTrack);
90 eflowFirstIntENUM j1st = integrator->getFirstIntLayer();
91
92 /*Save j1st info */
93 thisEfRecTrack->setLayerHED(j1st);
94
95 /* Get parameters for j1st */
96 eflowRingSubtractionManager& cellSubtractionManager = thisEfRecTrack->getCellSubtractionManager();
97 cellSubtractionManager.getOrdering(binnedParameters, trackE, trackEM1eta, j1st,useLegacyEnergyBinIndexing);
98
99 /* Set expected energy in the eflowRecTrack object */
100 const double expectedEnergy = energyP ? energyP->nnEnergyPrediction(thisEfRecTrack) : cellSubtractionManager.fudgeMean() * thisEfRecTrack->getTrack()->e();
101 const double expectedEnergySigma = std::fabs(cellSubtractionManager.fudgeStdDev() * thisEfRecTrack->getTrack()->e());
102
103 const std::vector<eflowTrackClusterLink*>* bestClusters_015 = thisEfRecTrack->getAlternativeClusterMatches("cone_015");
104 const std::vector<eflowTrackClusterLink*>* bestClusters_02 = thisEfRecTrack->getAlternativeClusterMatches("cone_02");
105
106 //First we calculate how much energy is found in a cone of 0.15
107 float totalE_015 = 0.0;
108
109 //This pointer can be a nullptr if no clusters were matched to a track in dR < 0.15
110 if (bestClusters_015){
111 for (eflowTrackClusterLink* thisLink : *bestClusters_015){
112 eflowRecCluster* thisBestCluster = thisLink->getCluster();
113 if (thisBestCluster){
114 const xAOD::CaloCluster* theCluster = thisBestCluster->getCluster();
115 if (theCluster) {
116 if (theCluster->e()>0.0){
117 totalE_015 += theCluster->e();
118 }
119 }
120 }
121 }
122 }//if vector of 0.15 clusters exists
123
124 double pull_015 = NAN;
125 if (expectedEnergySigma > 1e-6 ) pull_015 = (totalE_015-expectedEnergy)/expectedEnergySigma;
126 thisEfRecTrack->setpull15(pull_015);
127
128 double trackPt = thisEfRecTrack->getTrack()->pt();
129 //If the looked up expected energy deposit was 0.0, then pull_015 is NAN. In that case we should not try to apply the 2D cut described below.
130 if (!std::isnan(pull_015)){
131 //We use a 2D cut in the pull_015 and log10 of track pt plane to define a dense environment - if too dense then we disable the charged shower subtraction
132 //The cut values are such that above a track pt of 40 GeV we are very agrressive about disabling the cell by cell charged shower subtraction.
133 //The specific cut values were found by optimising the jet resolutions in the run 2 data taking environment.
134 //Further details can be found in these slides from C. Young:
135 //https://indico.cern.ch/event/396923/contributions/949333/attachments/795878/1090871/Presentation.pdf
136 if (pull_015 > 0.0 + (log10(40000)-log10(trackPt))*33.2 && 0.0 != expectedEnergySigma && bestClusters_015 && bestClusters_02){
137 thisEfRecTrack->setSubtracted(); //this tricks eflowRec into thinking this track was subtracted, and hence no further subtraction will be done
138 thisEfRecTrack->setIsInDenseEnvironment();
139 //recalculate the LHED and the ordering and find the new expected E + sigma of expected E (the new LHED can change the latter two values we find in the look up tables)
140 //we use a larger cone of 0.2 for this
141 std::vector<eflowRecCluster*> theBestEfRecClusters_02;
142 for (eflowTrackClusterLink* thisLink : *bestClusters_02) if (thisLink->getCluster()->getCluster()->e() > 0.0) theBestEfRecClusters_02.push_back(thisLink->getCluster());
143 integrator->measureNewClus(theBestEfRecClusters_02, thisEfRecTrack);
144 j1st = integrator->getFirstIntLayer();
145 cellSubtractionManager.getOrdering(binnedParameters, trackE, trackEM1eta, j1st,useLegacyEnergyBinIndexing);
146 thisEfRecTrack->setEExpect(cellSubtractionManager.fudgeMean() * trackE, std::fabs(cellSubtractionManager.fudgeStdDev()*trackE)*std::fabs(cellSubtractionManager.fudgeStdDev()*trackE));
147 }
148 else {
149 thisEfRecTrack->setEExpect(expectedEnergy, expectedEnergySigma*expectedEnergySigma);
150 }//ok to do subtraction, and so we just set the usual expected E + sigma of expected E needed for subtraction
151 }
152 }
153}
float nnEnergyPrediction(const eflowRecTrack *ptr) const
bool binExists(double e, double eta) const
void measureNewClus(const xAOD::CaloCluster *clus, const eflowTrackCaloPoints &trackCalo)
eflowFirstIntENUM getFirstIntLayer() const
bool getOrdering(const eflowEEtaBinnedParameters *binnedParameters, double e, double eta, eflowFirstIntENUM j1st, bool useLegacyEnergyBinIndexing)
eflowFirstIntRegions::J1STLAYER eflowFirstIntENUM
bool trackPt(const xAOD::TauJet &, const xAOD::TauTrack &track, float &out)
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.

Member Data Documentation

◆ m_eflowRecClusters

std::vector<eflowRecCluster*> eflowCaloObject::m_eflowRecClusters
private

Definition at line 82 of file eflowCaloObject.h.

◆ m_eflowRecTracks

std::vector<eflowRecTrack*> eflowCaloObject::m_eflowRecTracks
private

Definition at line 92 of file eflowCaloObject.h.

◆ m_trackClusterLinks

std::vector<std::pair<eflowTrackClusterLink*,std::pair<float,float> > > eflowCaloObject::m_trackClusterLinks
private

Definition at line 89 of file eflowCaloObject.h.


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