ATLAS Offline Software
Loading...
Searching...
No Matches
CaloHadDMCoeffFit.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef CALOLOCALHADCALIB_CALOHADDMCOEFFFIT_H
6#define CALOLOCALHADCALIB_CALOHADDMCOEFFFIT_H
7
21
22#include <vector>
23#include <string>
24#include <math.h>
26
30class TProfile;
31class TH2F;
32class TH1F;
33class TProfile2D;
34
35
37 public :
38
39 // average and rms on-the-flight calculation
40 class PrepData {
41 public:
42 unsigned int n_entries;
43 double m_sum;
44 double m_aver;
45 double m_rms;
46 double m_sw;
47 PrepData(): n_entries(0),m_sum(0.0),m_aver(0.0),m_rms(0.0),m_sw(0.0) { }
49 void add(double xx, double w=1.0)
50 {
51 n_entries++;
52 m_sum += xx;
53 m_rms = (m_sw/(m_sw+w))*(m_rms+(w/(m_sw+w))*(xx-m_aver)*(xx-m_aver));
54 m_aver = m_aver+(xx-m_aver)*w/(m_sw+w);
55 m_sw += w;
56 }
57 unsigned int size() { return n_entries;}
58 };
59
60 // fit data
61 class FitData {
62 public:
63 bool isOK;
64 float p0, p1, s0, s1;
65 std::string descr;
66 FitData() : isOK(true), p0(0), p1(0), s0(0), s1(0), descr("def") {}
67 FitData(float the_p0, float the_s0, float the_p1, float the_s1)
68 : isOK(true), p0(the_p0), p1(the_p1), s0(the_s0), s1(the_s1), descr("def") {}
69 void getInverted(float &p0inv, float &s0inv, float &p1inv, float &s1inv)
70 {
71 p0inv=0.0; s0inv=0.0; p1inv=0.0; s1inv=0.0;
72 if(p1!=0) {
73 p0inv = -1.0*p0/p1;
74 p1inv = 1./p1;
75 s0inv = sqrt((s0*s0)/(p1*p1) + (s1*s1)*(p0*p0)/(p1*p1*p1*p1));
76 s1inv = s1/(p1*p1);
77 }
78 }
79 };
80
83
84 CaloLocalHadCoeff * process(CaloHadDMCoeffData *myData, CaloLocalHadCoeff *myHadDMCoeff, bool isSingleParticle=true, bool tbflag=false);
85 void make_report(std::string &sfname);
86 void SetNormalizationType(std::string &stype) {m_NormalizationType = stype;}
87
88 private:
91
93
98
101
102 std::vector<PrepData *> m_engClus;
103 std::vector<PrepData *> m_engPrep;
104 std::vector<PrepData *> m_engDm;
105 std::vector<PrepData *> m_engDmOverClus;
106 std::vector<TProfile *> m_hp_DmVsPrep;
107 std::vector<TH2F *> m_h2_DmVsPrep;
108 std::vector<TH1F *> m_h1_engDmOverClus;
109 std::vector<FitData *> m_FitData;
110
111 std::vector<TProfile2D *> m_hp2_DmWeight;
112
113 void clear();
114 double ProfileRefiner(TProfile *pH, double ev_ratio=0.92);
115 double GetAverageWithoutRightTail(TH1F *pH, double ev_ratio=0.92);
116 int getFirstEnerLambdaBin(int ibin);
117
120};
121
122#endif
123
Define macros for attributes used to control the static checker.
#define ATLAS_NOT_THREAD_SAFE
getNoisyStrip() Find noisy strips from hitmaps and write out into xml/db formats
Data to read from special DeadMaterialTree.
FitData(float the_p0, float the_s0, float the_p1, float the_s1)
void getInverted(float &p0inv, float &s0inv, float &p1inv, float &s1inv)
void add(double xx, double w=1.0)
To fill and fit TProfile histograms using special dead material tree.
CaloLocalHadCoeffHelper * m_HadDMHelper
std::vector< PrepData * > m_engClus
std::string m_NormalizationType
CaloHadDMCoeffData * m_data
double ProfileRefiner(TProfile *pH, double ev_ratio=0.92)
CaloLocalHadCoeff * m_HadDMCoeff
std::vector< PrepData * > m_engPrep
std::vector< FitData * > m_FitData
std::vector< PrepData * > m_engDmOverClus
double GetAverageWithoutRightTail(TH1F *pH, double ev_ratio=0.92)
CaloHadDMCoeffFit(const CaloHadDMCoeffFit &)
std::vector< PrepData * > m_engDm
std::vector< TProfile * > m_hp_DmVsPrep
void SetNormalizationType(std::string &stype)
std::vector< TProfile2D * > m_hp2_DmWeight
std::vector< TH1F * > m_h1_engDmOverClus
std::vector< TH2F * > m_h2_DmVsPrep
int getFirstEnerLambdaBin(int ibin)
Hold binned correction data for local hadronic calibration procedure.
const std::string process
TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)