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