BOSS 7.0.6
BESIII Offline Software System
Loading...
Searching...
No Matches
calib_etf_bunch.cxx
Go to the documentation of this file.
2#include "TF1.h"
3
5
6 nKind = nbunch;
8
12 nCanvas = 8;
13 CanvasName.push_back( static_cast<string>("Offset-bunch0") );
14 CanvasName.push_back( static_cast<string>("Offset-bunch1") );
15 CanvasName.push_back( static_cast<string>("Offset-bunch2") );
16 CanvasName.push_back( static_cast<string>("Offset-bunch3") );
17 CanvasName.push_back( static_cast<string>("Sigma-bunch0") );
18 CanvasName.push_back( static_cast<string>("Sigma-bunch1") );
19 CanvasName.push_back( static_cast<string>("Sigma-bunch2") );
20 CanvasName.push_back( static_cast<string>("Sigma-bunch3") );
21
22 nGraphPerCanvas.push_back(1);
23 nGraphPerCanvas.push_back(1);
24 nGraphPerCanvas.push_back(1);
25 nGraphPerCanvas.push_back(1);
26 nGraphPerCanvas.push_back(1);
27 nGraphPerCanvas.push_back(1);
28 nGraphPerCanvas.push_back(1);
29 nGraphPerCanvas.push_back(1);
30
31 int numGraphs = 0;
32 std::vector<unsigned int>::iterator iter = nGraphPerCanvas.begin();
33 for( ; iter!=nGraphPerCanvas.end(); iter++ ) {
34 numGraphs = numGraphs + (*iter);
35 }
36 if( numGraphs != nGraphTotalBunch ) {
37 cout << "tofcalgsec::calib_barrel_common: the number of Graphs is NOT reasonable!!!" << endl;
38 exit(0);
39 }
40
41 m_name = string("calib_etf_bunch");
42
43 const int tbin = 160;
44 const double tbegin = -0.8;
45 const double tend = 0.8;
46
47 char hname[256];
48 // histograms
49 for( unsigned int j=0; j<4; j++ ) {
50 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
51 if( k==0 ) { sprintf( hname, "tleft-bunch%i", j); }
52 else if( k==1 ) { sprintf( hname, "tright-bunch%i", j); }
53 else if( k==2 ) { sprintf( hname, "tcombine-bunch%i", j); }
54 m_fitresult.push_back( HepVector( nParEtfBunch,0 ) );
55 m_histograms.push_back( new TH1F( hname, hname, tbin, tbegin, tend ) );
56 }
57 }
58 m_fitresult.push_back( HepVector( nParEtfBunch,0 ) );
59
60 modpos.resize( nBinPerCounter );
61 modposerr.resize( nBinPerCounter );
62 for( unsigned int i=0; i<nBinPerCounter; i++ ) {
63 modpos[i] = 1.0*i;
64 modposerr[i] = 0.5;
65 }
66
67 return;
68}
69
70
72 m_fitresult.clear();
73 modpos.clear();
74 modposerr.clear();
75}
76
77
78void calib_etf_bunch::calculate( RecordSet*& data, unsigned int icounter ) {
79
80 std::cout << setiosflags(ios::left) << setw(10) << icounter << setw(8) << data->size() << setw(30) << name() << std::endl;
81
82 if( data->size() > 0 ) {
83 std::vector<Record*>::iterator iter = data->begin();
84 for( ; iter!=data->end(); iter++ ) {
85 fillRecord( (*iter) );
86 }
87 }
88
89 if( icounter==(NEtf*NStrip-1) ) {
90 fitHistogram();
91 fillGraph();
92 fitGraph();
93 }
94
95 return;
96}
97
98
99void calib_etf_bunch::fillRecord( const Record* r ) {
100 double t0 = r->phi();
101 int ibunch = static_cast<int>(t0/(12000./499.8/nKind)+0.1)%nKind;
102 if( ibunch<0 || ibunch>int(nKind) ) return;
103
104 std::vector<TH1F*>::iterator iter = m_histograms.begin();
105 (*(iter+nBinPerCounter*ibunch+0))->Fill( r->tleft() );
106 (*(iter+nBinPerCounter*ibunch+1))->Fill( r->tright() );
107 (*(iter+nBinPerCounter*ibunch+2))->Fill( r->t0() );
108
109 return;
110}
111
112
113void calib_etf_bunch::fitHistogram() {
114 TF1* g = new TF1("g", "gaus");
115 g->SetLineColor(2);
116 g->SetLineWidth(1);
117
118 std::vector<TH1F*>::iterator iter1 = m_histograms.begin();
119 std::vector<HepVector>::iterator iter2 = m_fitresult.begin();
120 for( ; iter1!=m_histograms.end(); iter1++, iter2++ ) {
121 (*iter1)->Fit( g, "Q");
122 (*iter2)[0] = g->GetParameter(1);
123 (*iter2)[1] = g->GetParError(1);
124 (*iter2)[2] = g->GetParameter(2);
125 (*iter2)[3] = g->GetParError(2);
126 }
127
128 return;
129}
130
131
132void calib_etf_bunch::fillGraph() {
133
134 char gname1[256], gname2[256];
135 TH1F *gra[nGraphTotalBunch];
136
137 // 8canvas all counter,
138 // 1. offset of tleft, tright and tcombine of bunch0 --- gra[0]
139 // 2. offset of tleft, tright and tcombine of bunch1 --- gra[1]
140 // 3. offset of tleft, tright and tcombine of bunch2 --- gra[2]
141 // 4. offset of tleft, tright and tcombine of bunch3 --- gra[3]
142 // 5. sigma of tleft, tright and tcombine of bunch0 --- gra[4]
143 // 6. sigma of tleft, tright and tcombine of bunch1 --- gra[5]
144 // 7. sigma of tleft, tright and tcombine of bunch2 --- gra[6]
145 // 8. sigma of tleft, tright and tcombine of bunch3 --- gra[7]
146
147
148 std::vector<double> toffset, toffseterr;
149 std::vector<double> tsigma, tsigmaerr;
150 toffset.resize( nBinPerCounter );
151 toffseterr.resize( nBinPerCounter );
152 tsigma.resize( nBinPerCounter );
153 tsigmaerr.resize( nBinPerCounter );
154
155 unsigned int number = 0;
156 std::vector<HepVector>::iterator iter = m_fitresult.begin();
157 for( unsigned int j=0; j<4; j++ ) {
158 sprintf( gname1, "bunch%i-offset", j );
159 sprintf( gname2, "bunch%i-sigma", j );
160
161 gra[j] = new TH1F( gname1, gname1, nBinPerCounter, -0.5, 2.5 );
162 gra[j+4] = new TH1F( gname2, gname2, nBinPerCounter, -0.5, 2.5 );
163
164 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
165 number = j*nBinPerCounter + k;
166 toffset[k] = (*(iter+number))[0];
167 toffseterr[k] = (*(iter+number))[1];
168 tsigma[k] = (*(iter+number))[2];
169 tsigmaerr[k] = (*(iter+number))[3];
170 gra[j]->SetBinContent( k+1, toffset[k] );
171 gra[j]->SetBinError( k+1, toffseterr[k] );
172 gra[j+4]->SetBinContent( k+1, tsigma[k] );
173 gra[j+4]->SetBinError( k+1, tsigmaerr[k] );
174 }
175 }
176
177 for( int j=0; j<nGraphTotalBunch; j++, j++ ) {
178 gra[j]->SetMarkerSize(1.5);
179 gra[j]->SetMarkerStyle(20);
180 gra[j]->SetMarkerColor(2);
181 m_graphs.push_back( gra[j] );
182
183 gra[j+1]->SetMarkerSize(1.5);
184 gra[j+1]->SetMarkerStyle(21);
185 gra[j+1]->SetMarkerColor(4);
186 m_graphs.push_back( gra[j+1] );
187 }
188
189 return;
190}
191
192
193void calib_etf_bunch::fitGraph() {
194 TF1* p0 = new TF1("p0", "pol0");
195 p0->SetLineColor(1);
196 p0->SetLineWidth(1);
197
198 X = HepVector( 2, 0 );
199 std::vector<TH1F*>::iterator iter=m_graphs.begin();
200 for( int i=0; i<nParEtfBunch; i++ ) {
201 (*(iter+i))->Fit( "p0", "Q" );
202 X[0] = p0->GetParameter(0);
203 X[1] = p0->GetParError(0);
204 m_result.push_back( X );
205 }
206
207 return;
208}
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
std::vector< Record * > RecordSet
Definition: TofDataSet.h:98
const unsigned int NStrip
Definition: TofDataSet.h:15
const unsigned int NEtf
Definition: TofDataSet.h:14
const int nParEtfBunch
const int nGraphTotalBunch
double tleft() const
Definition: TofDataSet.h:59
double tright() const
Definition: TofDataSet.h:60
double t0() const
Definition: TofDataSet.h:68
double phi() const
Definition: TofDataSet.h:65
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
void calculate(RecordSet *&data, unsigned int ibunch)
calib_etf_bunch(const unsigned int nbunch)
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)