BOSS 7.0.9
BESIII Offline Software System
Loading...
Searching...
No Matches
calib_barrel_common.cxx
Go to the documentation of this file.
2#include "TF1.h"
3
5
6 nKind = 4;
7 nBinPerCounter = nzbin;
8
12 nCanvas = 4;
13 CanvasName.push_back( static_cast<string>("Offset") );
14 CanvasName.push_back( static_cast<string>("Offset-PlusMinus") );
15 CanvasName.push_back( static_cast<string>("Sigma") );
16 CanvasName.push_back( static_cast<string>("Sigma-TCorrelation") );
17 nGraphPerCanvas.push_back(2);
18 nGraphPerCanvas.push_back(2);
19 nGraphPerCanvas.push_back(2);
20 nGraphPerCanvas.push_back(3);
21
22 int numGraphs = 0;
23 std::vector<unsigned int>::iterator iter = nGraphPerCanvas.begin();
24 for( ; iter!=nGraphPerCanvas.end(); iter++ ) {
25 numGraphs = numGraphs + (*iter);
26 }
27 if( numGraphs != nGraphTotalCommon ) {
28 cout << "tofcalgsec::calib_barrel_common: the number of Graphs is NOT reasonable!!!" << endl;
29 exit(0);
30 }
31
32 m_name = string("calib_barrel_common");
33
34 const int tbin = 150;
35 const double tbegin = -1.5;
36 const double tend = 1.5;
37
38 char hname[256];
39 // histograms
40 for( unsigned int j=0; j<nKind; j++ ) {
41 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
42 if( j==0 ) { sprintf( hname, "tleft-z%i", k); }
43 else if( j==1 ) { sprintf( hname, "tright-z%i", k); }
44 else if( j==2 ) { sprintf( hname, "tplus-z%i", k); }
45 else if( j==3 ) { sprintf( hname, "tminus-z%i", k); }
46 m_fitresult.push_back( HepVector( nParCommon,0 ) );
47 m_histograms.push_back( new TH1F( hname, hname, tbin, tbegin, tend ) );
48 }
49 }
50 sprintf( hname, "deltaT" );
51 m_fitresult.push_back( HepVector( nParCommon,0 ) );
52 m_histograms.push_back( new TH1F( hname, hname, tbin, tbegin, tend ) );
53
54 zpos.resize( nBinPerCounter );
55 zposerr.resize( nBinPerCounter );
56 zstep = ( zend - zbegin )/nBinPerCounter;
57 for( unsigned int i=0; i<nBinPerCounter; i++ ) {
58 zpos[i] = zbegin + ( i+0.5 )*zstep;
59 zposerr[i] = 0.5*zstep;
60 }
61
62}
63
64
66 m_fitresult.clear();
67 zpos.clear();
68 zposerr.clear();
69}
70
71
72void calib_barrel_common::calculate( RecordSet*& data, unsigned int icounter ) {
73
74 std::cout << setiosflags(ios::left) << setw(10) << icounter << setw(8) << data->size() << setw(30) << name() << std::endl;
75
76 if( data->size() > 0 ) {
77 std::vector<Record*>::iterator iter = data->begin();
78 for( ; iter!=data->end(); iter++ ) {
79 fillRecord( (*iter) );
80 }
81 }
82
83 if( icounter==(NBarrel-1) ) {
84 fitHistogram();
85 fillGraph();
86 fitGraph();
87 }
88
89 return;
90}
91
92
93
94void calib_barrel_common::fillRecord( const Record* r ) {
95 double zhit = r->zrhit();
96 if( zhit<zbegin || zhit>zend ) return;
97 int zbin = static_cast<int>((zhit-zbegin)/zstep);
98 if( ( zbin<0 ) || ( zbin>static_cast<int>(nBinPerCounter-1) ) ) {
99 cout << "tofcalgsec::calib_barrel_common: zhit is out of range, zhit=" << zhit << " zbin=" << zbin << endl;
100 return;
101 }
102
103 std::vector<TH1F*>::iterator iter = m_histograms.begin();
104 (*(iter+zbin))->Fill( r->tleft() );
105 (*(iter+nBinPerCounter+zbin))->Fill( r->tright() );
106 (*(iter+2*nBinPerCounter+zbin))->Fill( (r->tleft()+r->tright())/2.0 );
107 (*(iter+3*nBinPerCounter+zbin))->Fill( (r->tleft()-r->tright())/2.0 );
108 (*(iter+nKind*nBinPerCounter))->Fill( r->tleft() );
109 (*(iter+nKind*nBinPerCounter))->Fill( r->tright() );
110
111 return;
112}
113
114
115void calib_barrel_common::fitHistogram() {
116 TF1* g = new TF1("g", "gaus");
117 g->SetLineColor(2);
118 g->SetLineWidth(1);
119
120 std::vector<TH1F*>::iterator iter1 = m_histograms.begin();
121 std::vector<HepVector>::iterator iter2 = m_fitresult.begin();
122 for( ; iter1!=m_histograms.end(); iter1++, iter2++ ) {
123 (*iter1)->Fit( g, "Q");
124 (*iter2)[0] = g->GetParameter(1);
125 (*iter2)[1] = g->GetParError(1);
126 (*iter2)[2] = g->GetParameter(2);
127 (*iter2)[3] = g->GetParError(2);
128 }
129
130 iter2 = m_fitresult.end() - 1;
131 X[2] = (*iter2)[0];
132 X[3] = (*iter2)[1];
133
134 return;
135
136}
137
138
139void calib_barrel_common::fillGraph() {
140
141 char gname1[256], gname2[256];
142 TH1F *gra[nGraphTotalCommon];
143
144 // 4 canvas all counter,
145 // 1. offset of tleft and tright vs z --- gra[0] and gra[1]
146 // 2. common of tleft and tright vs z --- gra[2] and gra[3]
147 // 3. offset of tplus and tminus vs z --- gra[4] and gra[5]
148 // 4. common of tplus, tminus and T_Correlation vs z
149 // --- gra[6], gra[7] and gra[8]
150
151 std::vector<double> toffset, toffseterr;
152 std::vector<double> tsigma, tsigmaerr;
153 toffset.resize( nBinPerCounter );
154 toffseterr.resize( nBinPerCounter );
155 tsigma.resize( nBinPerCounter );
156 tsigmaerr.resize( nBinPerCounter );
157
158 unsigned int number = 0;
159 std::vector<HepVector>::iterator iter = m_fitresult.begin();
160 for( unsigned int j=0; j<nKind; j++ ) {
161 if( j==0 ) {
162 sprintf( gname1, "tlefttright-offset" );
163 sprintf( gname2, "tlefttright-sigma" );
164 }
165 else if( j==1 ) {
166 sprintf( gname1, "tcommon-offset" );
167 sprintf( gname2, "tcommon-sigma" );
168 }
169 else if( j==2 ) {
170 sprintf( gname1, "tplusminus-offset" );
171 sprintf( gname2, "tplusminus-sigma" );
172 }
173 else if( j==3 ) {
174 sprintf( gname1, "tcorrelation-offset" );
175 sprintf( gname2, "tcorrelation-sigma" );
176 }
177
178 gra[j] = new TH1F( gname1, gname1, nBinPerCounter, zbegin, zend );
179 gra[j+4] = new TH1F( gname2, gname2, nBinPerCounter, zbegin, zend );
180
181 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
182 number = j*nBinPerCounter + k;
183 toffset[k] = (*(iter+number))[0];
184 toffseterr[k] = (*(iter+number))[1];
185 tsigma[k] = (*(iter+number))[2];
186 tsigmaerr[k] = (*(iter+number))[3];
187 gra[j]->SetBinContent( k+1, toffset[k] );
188 gra[j]->SetBinError( k+1, toffseterr[k] );
189 gra[j+4]->SetBinContent( k+1, tsigma[k] );
190 gra[j+4]->SetBinError( k+1, tsigmaerr[k] );
191 }
192 }
193
194 gra[nGraphTotalCommon-1] = new TH1F( "Sigma", "Sigma", nBinPerCounter, zbegin, zend );
195 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
196 number = (nKind-1)*nBinPerCounter + k;
197 double sigPlus = (*(iter+number-nBinPerCounter))[2];
198 double sigMinus = (*(iter+number))[2];
199 double sigErrPlus = (*(iter+number-nBinPerCounter))[3];
200 double sigErrMinus = (*(iter+number))[3];
201 if( sigPlus>sigMinus ) {
202 tsigma[k] = sqrt( sigPlus*sigPlus - sigMinus*sigMinus );
203 }
204 else {
205 tsigma[k] = 0.0 - sqrt( sigMinus*sigMinus - sigPlus*sigPlus );
206 }
207 tsigmaerr[k] = sqrt( sigErrPlus*sigErrPlus + sigErrMinus*sigErrMinus );
208
209 gra[nGraphTotalCommon-1]->SetBinContent( k+1, tsigma[k] );
210 gra[nGraphTotalCommon-1]->SetBinError( k+1, tsigmaerr[k] );
211 }
212
213 for( int j=0; j<nGraphTotalCommon; j++, j++ ) {
214 gra[j]->SetMarkerSize(1.5);
215 gra[j]->SetMarkerStyle(20);
216 gra[j]->SetMarkerColor(2);
217 if( j==4 ) {
218 gra[j]->SetMaximum( 0.22 );
219 gra[j]->SetMinimum( 0.07 );
220 }
221 else if( j==6 ) {
222 gra[j]->SetMaximum( 0.20 );
223 gra[j]->SetMinimum( -0.02 );
224 }
225 }
226 for( int j=1; j<nGraphTotalCommon; j++, j++ ) {
227 gra[j]->SetMarkerSize(1.5);
228 gra[j]->SetMarkerStyle(21);
229 gra[j]->SetMarkerColor(4);
230 }
231 gra[nGraphTotalCommon-1]->SetMarkerStyle(4);
232 gra[nGraphTotalCommon-1]->SetMarkerColor(6);
233
234 for( int j=0; j<nGraphTotalCommon; j++ ) {
235 m_graphs.push_back( gra[j] );
236 }
237
238 return;
239}
240
241
242void calib_barrel_common::fitGraph() {
243 TF1* p0 = new TF1("p0", "pol0");
244 p0->SetLineColor(1);
245 p0->SetLineWidth(1);
246
247 std::vector<TH1F*>::iterator iter=m_graphs.end()-1;
248 (*iter)->Fit( "p0", "Q" );
249 X[0] = p0->GetParameter(0);
250 X[1] = p0->GetParError(0);
251
252 m_result.push_back( X );
253 return;
254}
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)
TTree * data
const double zend
Definition: TofCalibFit.h:13
const double zbegin
Definition: TofCalibFit.h:12
const int nBarrelCommon
Definition: TofCalibFit.h:19
std::vector< Record * > RecordSet
Definition: TofDataSet.h:98
const unsigned int NBarrel
Definition: TofDataSet.h:12
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
const int nParCommon
const int nGraphTotalCommon
double tleft() const
Definition: TofDataSet.h:59
double tright() const
Definition: TofDataSet.h:60
double zrhit() const
Definition: TofDataSet.h:61
string m_name
Definition: TofCalibFit.h:50
std::vector< HepVector > m_result
Definition: TofCalibFit.h:55
unsigned int nCanvas
Definition: TofCalibFit.h:47
std::vector< TH1F * > m_histograms
Definition: TofCalibFit.h:53
unsigned int nBinPerCounter
Definition: TofCalibFit.h:41
unsigned int nKind
Definition: TofCalibFit.h:40
unsigned int nCanvasPerCounter
Definition: TofCalibFit.h:44
unsigned int nHistPerCounter
Definition: TofCalibFit.h:43
const string & name() const
Definition: TofCalibFit.h:27
HepVector X
Definition: TofCalibFit.h:51
std::vector< unsigned int > nGraphPerCanvas
Definition: TofCalibFit.h:48
std::vector< string > CanvasName
Definition: TofCalibFit.h:58
unsigned int nHistogram
Definition: TofCalibFit.h:46
std::vector< TH1F * > m_graphs
Definition: TofCalibFit.h:54
calib_barrel_common(const unsigned int nzbin)
void calculate(RecordSet *&data, unsigned int icounter)