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

Correction for time dependent HV effect in barrel presampler scale in 2012. More...

#include <corr_HV_EMBPS.h>

Collaboration diagram for corr_HV_EMBPS:

Public Member Functions

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

Static Private Member Functions

static float getRecoCorrection (float hv, float eta)
static float Respo (float e, float e_nominal, float tempe)
static float InvCharge (float e)
static float vdrift (float e, float tempe)
static float getDataCorrection (float hv, float eta)
static float getExtraScaleFactor ()

Private Attributes

TProfile2D * m_hHV [6][2] {}
TFile * m_file

Detailed Description

Correction for time dependent HV effect in barrel presampler scale in 2012.

Perform correction of E0_RAW (EMPBS) as function of run number, eta, phi to recover effective scale at 1600V at the beginning of 2012 Data are divided intwo 6 time snapshots

runs 200804-204910 (nominal 1600V for most EMBPS, 2000V EMECPS) runs 204932-205071 (20 EMBPS lines lowered to 1400V) runs 205112-211620 (22 EMBPS lines lowered from 1400 to 1200V) runs 211670-212272 (all EMBPS lowered to 1200V) (+changed in OFC / Mphys/Mcal) [main change for EMPBS] runs 212619-212742 (14 EMBPS lines lowered to 800V, all EMECPS lowered to 1600V) [main change for EMECPS] runs 212809-216432 ( 6 EMBPS lines lowered to 800V, 7 EMECPS lines lowered to 1200V)

The 1200->1600 and 800->1600 HV scale factors for the >211670 runs are derived from data For run<211670 an extra small correction is applied on top to account for the change in OFC / MphysMcal computation at run 211670 (using drift time for HV=1200V instead of old drift time for HV=2000V)

Definition at line 32 of file corr_HV_EMBPS.h.

Constructor & Destructor Documentation

◆ corr_HV_EMBPS()

corr_HV_EMBPS::corr_HV_EMBPS ( )

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

Definition at line 19 of file corr_HV_EMBPS.cxx.

19 {
20 std::string filename = PathResolverFindCalibFile("egammaLayerRecalibTool/v1/HVmaps_embps_2012.root");
21 m_file = TFile::Open(filename.c_str());
22 if (not m_file or m_file->IsZombie()) {
23 std::cerr << "FATAL: cannot open " << filename << std::endl;
24 }
25 else {
26 for (int iperiod=0;iperiod<6;iperiod++) {
27 for (int igap=0;igap<2;igap++) {
28 char name[12];
29 snprintf(name,sizeof(name),"h%d%d",iperiod+1,igap);
30 m_hHV[iperiod][igap]=(TProfile2D*) (m_file->Get(name));
31 if (not m_hHV[iperiod][igap]) {
32 std::cerr << "FATAL: cannot find " << name << std::endl;
33 }
34 }
35 }
36 }
37}
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
TProfile2D * m_hHV[6][2]

◆ ~corr_HV_EMBPS()

corr_HV_EMBPS::~corr_HV_EMBPS ( )

Definition at line 39 of file corr_HV_EMBPS.cxx.

40{
41 m_file->Close();
42}

Member Function Documentation

◆ getCorr()

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

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

inputs: run number, eta and phi

Definition at line 45 of file corr_HV_EMBPS.cxx.

46{
47 if (std::fabs(eta) > 1.5) return 1.;
48 if (run<200804) return 1.;
49 int iperiod;
50 if (run<204932) iperiod=0;
51 else if (run<205112) iperiod=1;
52 else if (run<211670) iperiod=2;
53 else if (run<212619) iperiod=3;
54 else if (run<212809) iperiod=4;
55 else iperiod=5;
56 //std::cout << "period " << iperiod << std::endl;
57
58 int ieta=(int)((eta+1.6)/0.4);
59 if (ieta<0 || ieta>7) return 1.;
60 if (phi<0.) phi=phi+2*M_PI;
61 int iphi=(int)(phi*(32./(2*M_PI)));
62 if (iphi<0) iphi=0;
63 if (iphi>31) iphi=31;
64 //std::cout << " ieta,iphi " << ieta << " " << iphi << std::endl;
65
66 float hv1 = m_hHV[iperiod][0]->GetBinContent(ieta+1,iphi+1);
67 float hv2 = m_hHV[iperiod][1]->GetBinContent(ieta+1,iphi+1);
68
69 float c1= getRecoCorrection(hv1,eta);
70 float c2= getRecoCorrection(hv2,eta);
71 float delta1 = getDataCorrection(hv1,eta);
72 float delta2 = getDataCorrection(hv2,eta);
73
74 float extra1=1.;
75 if (run<211670 && hv1<1300.) extra1 = getExtraScaleFactor();
76 float extra2=1.;
77 if (run<211670 && hv2<1300.) extra2 = getExtraScaleFactor();
78
79 if ((c1+c2)<0.01) return 1.;
80 float oldScale = (c1+c2);
81 float newScale = (c1*delta1*extra1 + c2*delta2*extra2);
82
83 float newCorr = oldScale/newScale;
84
85 //std::cout << " run,eta,phi " << run << " " << eta << " " << phi << " " << " period "<< iperiod << std::endl;
86 //std::cout << " HV values " << hv1 << " " << hv2 << std::endl;
87 //std::cout << " reco scales " << c1 << " " << c2 << std::endl;
88 //std::cout << " delta scale " << delta1 << " " << delta2 << std::endl;
89 //std::cout << " extra " << extra1 << " " << extra2 << std::endl;
90 //std::cout << " old/new scales " << oldScale << " " << newScale << std::endl;
91 //std::cout << " newCorr " << newCorr << std::endl;
92
93 return newCorr;
94
95}
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
static float getDataCorrection(float hv, float eta)
static float getRecoCorrection(float hv, float eta)
static float getExtraScaleFactor()

◆ getDataCorrection()

float corr_HV_EMBPS::getDataCorrection ( float hv,
float eta )
staticprivate

Definition at line 151 of file corr_HV_EMBPS.cxx.

156{
157 if (hv>1300.) return 1.;
158 if (hv>1000.) {
159 if (eta<-1.2) return 0.9925;
160 if (eta<-0.8) return 0.9918;
161 if (eta<-0.4) return 0.9889;
162 if (eta<0.) return 0.9935;
163 if (eta<0.4) return 0.9908;
164 if (eta<0.8) return 0.99197;
165 if (eta<1.2) return 0.9974;
166 return 0.9971;
167 }
168 return 1.015;
169}

◆ getExtraScaleFactor()

float corr_HV_EMBPS::getExtraScaleFactor ( )
staticprivate

Definition at line 172 of file corr_HV_EMBPS.cxx.

173{
174 return 1.013;
175}

◆ getRecoCorrection()

float corr_HV_EMBPS::getRecoCorrection ( float hv,
float eta )
staticprivate

Definition at line 99 of file corr_HV_EMBPS.cxx.

100{
101 float nominal = 2000.;
102 float T = 88.37;
103
104 int ieta=(std::fabs(eta)/0.025);
105 float d;
106 if (ieta>=0 && ieta<16) d = 0.196; //cm
107 else if (ieta>=16 && ieta<32) d = 0.193; //cm
108 else if (ieta>=32 && ieta<48) d = 0.2; //cm
109 else d = 0.190; //cm
110
111 if (hv>(nominal-2.)) return 1.;
112 double efield=std::fabs(hv)/(d*1e+3);
113 double enominal=nominal/(d*1e+3);
114 double scale=Respo(efield, enominal,T);
115
116 return scale;
117}
static float Respo(float e, float e_nominal, float tempe)
unsigned long long T

◆ InvCharge()

float corr_HV_EMBPS::InvCharge ( float e)
staticprivate

Definition at line 131 of file corr_HV_EMBPS.cxx.

132{
133 float q = 1.;
134 if ( e > 2.) q=(1.+0.36/e);
135 else if (e>1e-6) q =(1.04+0.25/e);
136 return q;
137}

◆ Respo()

float corr_HV_EMBPS::Respo ( float e,
float e_nominal,
float tempe )
staticprivate

Definition at line 120 of file corr_HV_EMBPS.cxx.

121{
122 if (e < -999.) return 1.;
123 if (e < 0.01) return 0;
124 if ( e > e_nominal ) return 1;
125 float den = InvCharge(e)*vdrift(e_nominal,tempe);
126 float resp = den == 0 ? 0 : (InvCharge(e_nominal)*vdrift(e,tempe))/den;
127 return resp;
128}
static float InvCharge(float e)
static float vdrift(float e, float tempe)

◆ vdrift()

float corr_HV_EMBPS::vdrift ( float e,
float tempe )
staticprivate

Definition at line 140 of file corr_HV_EMBPS.cxx.

141{
142 const float T = tempe;
143 static const float P[6] = {-0.01481,-0.0075,0.141,12.4,1.627,0.317};
144 if ( e < -999.) return 0.;
145 float vd = (P[0]*T+1)*( P[2]*e*std::log(1+ (P[3]/e)) + P[4]*std::pow(e,P[5])) + P[1]*T; // vdrift formula walcowialk mm/micro_s
146 return vd;
147}
static Double_t P(Double_t *tt, Double_t *par)

Member Data Documentation

◆ m_file

TFile* corr_HV_EMBPS::m_file
private

Definition at line 64 of file corr_HV_EMBPS.h.

◆ m_hHV

TProfile2D* corr_HV_EMBPS::m_hHV[6][2] {}
private

Definition at line 63 of file corr_HV_EMBPS.h.

63{};

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