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

Correction for wrong HV EMEC presampler correction for 20.7 processed data (2015,2016) More...

#include <corr_HV_EMECPS.h>

Collaboration diagram for corr_HV_EMECPS:

Public Member Functions

 corr_HV_EMECPS ()
 constructor (initialization done there reading a root file for the HV maps per period
 ~corr_HV_EMECPS ()
float getCorr (int run, float eta, float phi) const
 get correction factor to apply to raw EMECPS energy : corrected raw EMECPS energy = correction factor * raw uncorrected EMECPS energy

Private Attributes

TH2D * m_HV [2] {}
TFile * m_file

Detailed Description

Correction for wrong HV EMEC presampler correction for 20.7 processed data (2015,2016)

Perform correction of E0_RAW (EMEC PS) as function of run number and phi

Definition at line 14 of file corr_HV_EMECPS.h.

Constructor & Destructor Documentation

◆ corr_HV_EMECPS()

corr_HV_EMECPS::corr_HV_EMECPS ( )

constructor (initialization done there reading a root file for the HV maps per period

Definition at line 16 of file corr_HV_EMECPS.cxx.

16 {
17 std::string filename = PathResolverFindCalibFile("egammaLayerRecalibTool/v2/emecps_hv_corr_2015-2016_207X.root");
18
19 m_file = TFile::Open(filename.c_str());
20 if (not m_file or m_file->IsZombie()) {
21 std::cerr << "FATAL: cannot open " << filename << std::endl;
22 }
23 else {
24 m_HV[0]=(TH2D*) (m_file->Get("hc"));
25 if (!m_HV[0]) {
26 std::cerr << " cannot read histogram hc " << std::endl;
27 }
28 m_HV[1]=(TH2D*) (m_file->Get("ha"));
29 if (!m_HV[1]) {
30 std::cerr << " cannot read histogram hc " << std::endl;
31 }
32 }
33}
std::string PathResolverFindCalibFile(const std::string &logical_file_name)

◆ ~corr_HV_EMECPS()

corr_HV_EMECPS::~corr_HV_EMECPS ( )

Definition at line 35 of file corr_HV_EMECPS.cxx.

36{
37 m_file->Close();
38}

Member Function Documentation

◆ getCorr()

float corr_HV_EMECPS::getCorr ( int run,
float eta,
float phi ) const

get correction factor to apply to raw EMECPS energy : corrected raw EMECPS energy = correction factor * raw uncorrected EMECPS energy

inputs: run number, eta and phi

Definition at line 41 of file corr_HV_EMECPS.cxx.

42{
43 float newCorr=1.;
44 if (eta>-1.5 && eta<1.5) return newCorr;
45 if (eta<-1.85 || eta > 1.85) return newCorr;
46
47 int iside=0;
48 if (eta>0.) iside=1;
49
50 if (!m_HV[iside]) {
51 std::cerr << " cannot find histogram to apply correction " << std::endl;
52 } else {
53 int ibin = std::as_const(*m_HV[iside]).GetXaxis()->FindFixBin(((Double_t)(run+0.1)));
54 if (ibin<1 || ibin>7) return newCorr;
55 if (phi>3.14159) phi=phi-6.283185;
56 int jbin = std::as_const(*m_HV[iside]).GetYaxis()->FindFixBin(phi);
57 if (jbin<1 || jbin>2048) return newCorr;
58 newCorr = m_HV[iside]->GetBinContent(ibin,jbin);
59 }
60 //std::cout << " in corr_HV_EMECPS run,eta,phi,corr " << run << " " << eta << " " << phi << " " << newCorr << std::endl;
61 return newCorr;
62
63}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method

Member Data Documentation

◆ m_file

TFile* corr_HV_EMECPS::m_file
private

Definition at line 29 of file corr_HV_EMECPS.h.

◆ m_HV

TH2D* corr_HV_EMECPS::m_HV[2] {}
private

Definition at line 28 of file corr_HV_EMECPS.h.

28{};

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