BOSS 7.0.9
BESIII Offline Software System
Loading...
Searching...
No Matches
ResiAlign.cpp
Go to the documentation of this file.
1#include "include/ResiAlign.h"
2#include "include/fun.h"
3
4#include <cmath>
5
6#include "TF1.h"
7#include "TStyle.h"
8#include "TCanvas.h"
9
11 cout << "Alignment type: ResiAlign" << endl;
12}
13
15}
16
17void ResiAlign::init(TObjArray* hlist, MdcCosGeom* pGeom){
18 m_pGeom = pGeom;
19 char hname[200];
20 m_hnTrk = new TH1F("mHNtrack", "", 10, -0.5, 9.5);
21 hlist->Add(m_hnTrk);
22
23 m_hnHit = new TH1F("mHNhit", "", 100, -0.5, 99.5);
24 hlist->Add(m_hnHit);
25
26 m_hlayHitmap = new TH1F("mHitmap", "", 43, -0.5, 42.5);
27 hlist->Add(m_hnHit);
28
29 m_hresAll = new TH1F("mHResAllInc", "", 200, -1.0, 1.0);
30 hlist->Add(m_hresAll);
31
32 m_hresInn = new TH1F("mHResInnInc", "", 200, -1.0, 1.0);
33 hlist->Add(m_hresInn);
34
35 m_hresStp = new TH1F("mHResStpInc", "", 200, -1.0, 1.0);
36 hlist->Add(m_hresStp);
37
38 m_hresOut = new TH1F("mHResOutInc", "", 200, -1.0, 1.0);
39 hlist->Add(m_hresOut);
40
41 for(int lay=0; lay<LAYERNMAX; lay++){
42 sprintf(hname, "mRes_Layer%02d", lay);
43 m_hresLay[lay] = new TH1F(hname, "", 200, -1.0, 1.0);
44 hlist->Add(m_hresLay[lay]);
45 }
46
47 for(int iEP=0; iEP<NEP; iEP++){
48 m_gr[iEP] = new TGraph();
49 sprintf(hname, "mgrResi%02d", iEP);
50 m_gr[iEP]->SetName(hname);
51 hlist->Add(m_gr[iEP]);
52
53 m_grSinPhi[iEP] = new TGraph();
54 sprintf(hname, "mgrResi_sinPhi%02d", iEP);
55 m_grSinPhi[iEP]->SetName(hname);
56 hlist->Add(m_grSinPhi[iEP]);
57
58 m_grCosPhi[iEP] = new TGraph();
59 sprintf(hname, "mgrResi_cosPhi%02d", iEP);
60 m_grCosPhi[iEP]->SetName(hname);
61 hlist->Add(m_grCosPhi[iEP]);
62
63 m_npoint[iEP] = 0;
64 }
65}
66
67void ResiAlign::mergeHist(TFile* fhist){
68 char hname[200];
69 TH1F* hist;
70 hist = (TH1F*)fhist->Get("HNtrack");
71 m_hnTrk->Add(hist);
72
73 hist = (TH1F*)fhist->Get("HNhit");
74 m_hnHit->Add(hist);
75
76 hist = (TH1F*)fhist->Get("Hitmap");
77 m_hlayHitmap->Add(hist);
78
79 hist = (TH1F*)fhist->Get("HResAllInc");
80 m_hresAll->Add(hist);
81
82 hist = (TH1F*)fhist->Get("HResInnInc");
83 m_hresInn->Add(hist);
84
85 hist = (TH1F*)fhist->Get("HResStpInc");
86 m_hresStp->Add(hist);
87
88 hist = (TH1F*)fhist->Get("HResOutInc");
89 m_hresOut->Add(hist);
90
91 for(int lay=0; lay<LAYERNMAX; lay++){
92 sprintf(hname, "Res_Layer%02d", lay);
93 hist = (TH1F*)fhist->Get(hname);
94 m_hresLay[lay]->Add(hist);
95 }
96
97 for(int iEP=0; iEP<NEP; iEP++){
98 sprintf(hname, "grResi%02d", iEP);
99 TGraph* gr = (TGraph*)fhist->Get(hname);
100 int np = gr->GetN();
101 double xx;
102 double yy;
103 for(int i=0; i<np; i++){
104 gr->GetPoint(i, xx, yy);
105 m_gr[iEP]->SetPoint(m_npoint[iEP], xx, yy);
106 m_grSinPhi[iEP]->SetPoint(m_npoint[iEP], sin(xx), yy);
107 m_grCosPhi[iEP]->SetPoint(m_npoint[iEP], cos(xx), yy);
108 m_npoint[iEP]++;
109 }
110 }
111}
112
114 int iEP;
115 double par[3];
116 double err[3];
117 double dx;
118 double dy;
119 double rz;
120 double rLayer[] = { 120.225, 205.0, 237.55, 270.175, 302.625, 334.775, 366.65, 500.0,
121 120.225, 205.0, 237.55, 270.175, 302.625, 334.775, 366.65, 500.0 };
122
123 TCanvas c1("c1", "c1", 10, 10, 700, 500);
124 c1.SetFillColor(10);
125
126 TF1* fResPhi = new TF1("fResPhi", funResi, 0, PI2, 3);
127 fResPhi->SetParameter(0, 0.0);
128 fResPhi->SetParameter(1, 0.0);
129 fResPhi->SetParameter(2, 0.0);
130
131 for(iEP=0; iEP<NEP; iEP++){
132 if((m_gr[iEP]->GetN()) > 500){
133 // align dx, dy, rz
134 m_gr[iEP]->Fit("fResPhi", "V");
135 par[0] = fResPhi->GetParameter(0);
136 par[1] = fResPhi->GetParameter(1);
137 par[2] = fResPhi->GetParameter(2);
138 err[0] = fResPhi->GetParError(0);
139 err[1] = fResPhi->GetParError(1);
140 err[2] = fResPhi->GetParError(2);
141
142 // align dx and rz
143// m_grSinPhi[iEP]->Fit("pol1");
144// par[0] = m_grSinPhi[iEP]->GetFunction("pol1")->GetParameter(0);
145// par[1] = m_grSinPhi[iEP]->GetFunction("pol1")->GetParameter(1);
146// par[2] = 0.0;
147// err[0] = m_grSinPhi[iEP]->GetFunction("pol1")->GetParError(0);
148// err[1] = m_grSinPhi[iEP]->GetFunction("pol1")->GetParError(1);
149// err[2] = 0.0;
150
151 // align dy
152// m_grCosPhi[iEP]->Fit("pol1");
153// par[0] = 0.0;
154// par[1] = 0.0;
155// par[2] = m_grCosPhi[iEP]->GetFunction("pol1")->GetParameter(1);
156// err[0] = 0.0;
157// err[1] = 0.0;
158// err[2] = m_grCosPhi[iEP]->GetFunction("pol1")->GetParError(1);
159
160 dx = -1.0 * par[1];
161 dy = par[2];
162 rz = par[0] / rLayer[iEP];
163
164 if (7==iEP || 15==iEP) {
165 dx = 0.0;
166 dy = 0.0;
167 rz = 0.0;
168 par[0] = 0.0;
169 par[1] = 0.0;
170 par[2] = 0.0;
171 }
172 alignPar->setDelDx(iEP, dx);
173 alignPar->setDelDy(iEP, dy);
174 alignPar->setDelRz(iEP, rz);
175
176 alignPar->setErrDx(iEP, err[1]);
177 alignPar->setErrDy(iEP, err[2]);
178 alignPar->setErrRz(iEP, err[0]/rLayer[iEP]);
179 }
180 }
181 renameHist();
182 delete fResPhi;
183}
184
185void ResiAlign::renameHist(){
186 char hname[200];
187 m_hnTrk->SetName("HNtrack");
188 m_hnHit->SetName("HNhit");
189 m_hlayHitmap->SetName("Hitmap");
190 m_hresAll->SetName("HResAllInc");
191 m_hresInn->SetName("HResInnInc");
192 m_hresStp->SetName("HResStpInc");
193 m_hresOut->SetName("HResOutInc");
194 for(int lay=0; lay<LAYERNMAX; lay++){
195 sprintf(hname, "Res_Layer%02d", lay);
196 m_hresLay[lay]->SetName(hname);
197 }
198 for(int iEP=0; iEP<NEP; iEP++){
199 sprintf(hname, "grResi%02d", iEP);
200 m_gr[iEP]->SetName(hname);
201
202 sprintf(hname, "grResi_sinPhi%02d", iEP);
203 m_grSinPhi[iEP]->SetName(hname);
204
205 sprintf(hname, "grResi_cosPhi%02d", iEP);
206 m_grCosPhi[iEP]->SetName(hname);
207 }
208}
209
210Double_t ResiAlign::funResi(Double_t* x, Double_t* par){
211 Double_t val;
212 val = par[0] + par[1]*sin(x[0]) + par[2]*cos(x[0]);
213 return val;
214}
double sin(const BesAngle a)
Definition: BesAngle.h:210
double cos(const BesAngle a)
Definition: BesAngle.h:213
sprintf(cut,"kal_costheta0_em>-0.93&&kal_costheta0_em<0.93&&kal_pxy0_em>=0.05+%d*0.1&&kal_pxy0_em<0.15+%d*0.1&&NGch>=2", j, j)
TGraph * gr
void setErrRz(int iEP, double val)
void setDelDy(int iEP, double val)
void setDelRz(int iEP, double val)
void setErrDx(int iEP, double val)
void setDelDx(int iEP, double val)
void setErrDy(int iEP, double val)
static Double_t funResi(double *x, double *par)
void mergeHist(TFile *fhist)
Definition: ResiAlign.cpp:67
void align(MdcAlignPar *alignPar)
Definition: ResiAlign.cpp:113
void init(TObjArray *hlist, MdcCosGeom *pGeom)
Definition: ResiAlign.cpp:17