ATLAS Offline Software
egPhotonWrtPoint.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 
8 #include "xAODEgamma/Egamma.h"
9 
11  double z) {
12  std::pair<double, float> RZ1 = {0., 0.};
13  float etaBE1 = ph.caloCluster()->etaBE(1);
14  if(std::abs(etaBE1) < 10) RZ1 = CP::ShowerDepthTool::getRZ(etaBE1, 1);
15 
16  double rCalo = RZ1.first;
17  double zCalo = RZ1.second;
18  double correctedZ = zCalo - z;
19  double eta = rCalo>0. ? std::asinh(correctedZ / rCalo) : ph.eta();
20  return {ph.e() / std::cosh(eta), eta, ph.phi()};
21 }
22 
24  auto corr = photonWrtPoint::PtEtaPhiWrtZ(ph, z);
25  ph.setP4(corr.pt, corr.eta, corr.phi, ph.m());
26 }
27 
egPhotonWrtPoint.h
ShowerDepthTool.h
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:83
xAOD::Egamma_v1::e
virtual double e() const override final
The total energy of the particle.
Definition: Egamma_v1.cxx:90
xAOD::Egamma_v1
Definition: Egamma_v1.h:56
photonWrtPoint::PtEtaPhiWrtZ
PtEtaPhi PtEtaPhiWrtZ(const xAOD::Egamma &ph, double z)
Function to get the kinematics of a photon cluster wrt (0,0,z0)
Definition: egPhotonWrtPoint.cxx:10
xAOD::CaloCluster_v1::etaBE
float etaBE(const unsigned layer) const
Get the eta in one layer of the EM Calo.
Definition: CaloCluster_v1.cxx:644
Egamma.h
photonWrtPoint::PtEtaPhi
egamma clusters kinematics are always wrt the ATLAS frame (0,0,0).
Definition: egPhotonWrtPoint.h:20
xAOD::Egamma_v1::setP4
void setP4(float pt, float eta, float phi, float m)
set the 4-vec
Definition: Egamma_v1.cxx:104
z
#define z
CP::ShowerDepthTool::getRZ
static std::pair< float, float > getRZ(float eta, int sampling)
Shower depth in R,Z for the given sampling.
Definition: ShowerDepthTool.cxx:107
xAOD::Egamma_v1::caloCluster
const xAOD::CaloCluster * caloCluster(size_t index=0) const
Pointer to the xAOD::CaloCluster/s that define the electron candidate.
Definition: Egamma_v1.cxx:388
python.Dumpers.asinh
def asinh(x)
helper methods ---------------------------------------------------------—
Definition: Dumpers.py:89
xAOD::Egamma_v1::phi
virtual double phi() const override final
The azimuthal angle ( ) of the particle.
Definition: Egamma_v1.cxx:75
photonWrtPoint::correctForZ
void correctForZ(xAOD::Egamma &ph, double z)
Function to modify in place the kinematics of a photon wrt (0,0,z0)
Definition: egPhotonWrtPoint.cxx:23
xAOD::Egamma_v1::eta
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Definition: Egamma_v1.cxx:70
xAOD::Egamma_v1::m
virtual double m() const override final
The invariant mass of the particle.
Definition: Egamma_v1.cxx:80