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