CGEM BOSS 6.6.5.h
BESIII Offline Software System
Loading...
Searching...
No Matches
calib_endcap_atten.cxx
Go to the documentation of this file.
2#include "TF1.h"
3
4//Q fit function
5static double endcapQFunc(double *x, double *par) {
6 double xx = x[0];
7 double f = par[0] + par[1]*(xx-44.5) + par[2]*(xx-44.5)*(xx-44.5);
8
9 return f;
10}
11
12
14
15 nKind = 1; // pulse height
16 nBinPerCounter = nrbin + 1;
17
20 CanvasPerCounterName.push_back( static_cast<string>("Pulse Height Most Probable Value") );
21 CanvasPerCounterName.push_back( static_cast<string>("Pulse Height Sigma") );
22 nGraphPerCanvasPerCounter.push_back(1);
23 nGraphPerCanvasPerCounter.push_back(1);
24
25 nHistogram = 0;
26 nCanvas = 2;
27 CanvasName.push_back( static_cast<string>("Pulse Height Most Probable Value vs TOF Counter Number") );
28 CanvasName.push_back( static_cast<string>("Pulse Height Sigma vs TOF Counter Number") );
29 nGraphPerCanvas.push_back(1);
30 nGraphPerCanvas.push_back(1);
31
32 int numGraphs1 = 0;
33 std::vector<unsigned int>::iterator iter = nGraphPerCanvasPerCounter.begin();
34 for( ; iter!=nGraphPerCanvasPerCounter.end(); iter++ ) {
35 numGraphs1 = numGraphs1 + (*iter);
36 }
37 if( numGraphs1 != nGraphEcAtten ) {
38 cout << "tofcalgsec::calib_endcap_atten: the number of Graphs is NOT reasonable!!!" << endl;
39 exit(0);
40 }
41 int numGraphs2 = 0;
42 iter = nGraphPerCanvas.begin();
43 for( ; iter!=nGraphPerCanvas.end(); iter++ ) {
44 numGraphs2 = numGraphs2 + (*iter);
45 }
46 if( numGraphs2 != nGraphEcAtten ) {
47 cout << "tofcalgsec::calib_endcap_atten: the number of Graphs is NOT reasonable!!!" << endl;
48 exit(0);
49 }
50
51 m_name = string("calib_endcap_atten");
52
53 const int qbin = 100;
54 const double qbegin = 0.0;
55 const double qend = 5000.0;
56
57 // histograms per counter
58 char hname[256];
59 for( unsigned int i=0; i<NEndcap; i++ ) {
60 m_result.push_back( HepVector(nEndcapAtten,0) );
61 for( unsigned int k=0; k<nrbin; k++ ) {
62 sprintf( hname, "Q-id%i-r%i", i, k);
63 m_histograms.push_back( new TH1F( hname, hname, qbin, qbegin, qend ) );
64 m_fitresult.push_back( HepVector(nParEcAtten,0) );
65 }
66 sprintf( hname, "Q0-id%i", i );
67 m_histograms.push_back( new TH1F( hname, hname, qbin, qbegin, qend ) );
68 m_fitresult.push_back( HepVector(nParEcAtten,0) );
69 }
70
71 rpos.resize( nrbin );
72 rposerr.resize( nrbin );
73 rstep = ( rend - rbegin )/nrbin;
74 for( unsigned int i=0; i<nrbin; i++ ) {
75 rpos[i] = rbegin + ( i+0.5 )*rstep;
76 rposerr[i] = 0.5*rstep;
77 }
78 itofid.resize( NEndcap );
79 itofiderr.resize( NEndcap );
80 itofidstep = 1.0;
81 for( unsigned int i=0; i<NEndcap; i++ ) {
82 itofid[i] = i*1.0;
83 itofiderr[i] = 0.5;
84 }
85
86}
87
88
90 m_fitresult.clear();
91 rpos.clear();
92 rposerr.clear();
93 itofid.clear();
94 itofiderr.clear();
95}
96
97
98void calib_endcap_atten::calculate( RecordSet*& data, unsigned int icounter ) {
99
100 std::cout << setiosflags(ios::left) << setw(10) << icounter << setw(8) << data->size() << setw(30) << name() << std::endl;
101
102 if( data->size() > 0 ) {
103 std::vector<Record*>::iterator iter = data->begin();
104 for( ; iter!=data->end(); iter++ ) {
105 fillRecord( (*iter), icounter );
106 }
107 fitHistogram( icounter );
108 fillGraph( icounter );
109 fitGraph( icounter );
110 }
111 else {
112 fillGraph( icounter ); // keep the m_graphs is not empty()
113 }
114
115 if( data->size() > 0 ) {
116 std::vector<Record*>::iterator iter = data->begin();
117 for( ; iter!=data->end(); iter++ ) {
118 updateData( (*iter), icounter );
119 fillRecordQ0( (*iter), icounter );
120 }
121 fitHistogramQ0( icounter );
122 }
123
124 if( icounter==(NEndcap-1) ) {
125 fillGraphQ0();
126 }
127
128 return;
129}
130
131
132void calib_endcap_atten::fillRecord( const Record* r, unsigned int icounter ) {
133
134 double rhit = r->zrhit();
135 if( rhit<rbegin || rhit>rend ) return;
136 int rbin = static_cast<int>((rhit-rbegin)/rstep);
137 if( rbin<0 ) { rbin = 0; }
138 else if( rbin>static_cast<int>(nBinPerCounter-1-1) ) {
139 cout << "tofcalgsec::calib_endcap_atten:fillRecord: rhit is out of range, rhit=" << rhit << " rbin=" << rbin << endl;
140 return;
141 }
142
143 std::vector<TH1F*>::iterator iter = m_histograms.begin() + icounter*nKind*nBinPerCounter + rbin;
144 (*iter)->Fill( r->qleft()*abs(r->theta()) );
145
146 return;
147}
148
149
150void calib_endcap_atten::fitHistogram( unsigned int icounter ) {
151 TF1* ld = new TF1("ld", "landau");
152 ld->SetLineColor(2);
153 ld->SetLineWidth(1);
154
155 std::vector<TH1F*>::iterator iter1 = m_histograms.begin() + icounter*nKind*nBinPerCounter;
156 std::vector<HepVector>::iterator iter2 = m_fitresult.begin() + icounter*nKind*nBinPerCounter;
157 for( unsigned int j=0; j<nBinPerCounter-1; j++, iter1++, iter2++ ) {
158 (*iter1)->Fit( ld, "Q");
159 (*iter2)[0] = ld->GetParameter(1);
160 (*iter2)[1] = ld->GetParError(1);
161 (*iter2)[2] = ld->GetParameter(2);
162 (*iter2)[3] = ld->GetParError(2);
163 }
164
165 return;
166
167}
168
169
170void calib_endcap_atten::fillGraph( unsigned int icounter ) {
171
172 char gname1[256], gname2[256];
173
174 // fill graphs
175 // 2 canvas per counter,
176 // 1. Q MPV vs z
177 // 2. Q sigma vs z
178 std::vector<double> toffset, toffseterr;
179 std::vector<double> tsigma, tsigmaerr;
180 toffset.resize( nBinPerCounter-1 );
181 toffseterr.resize( nBinPerCounter-1 );
182 tsigma.resize( nBinPerCounter-1 );
183 tsigmaerr.resize( nBinPerCounter-1 );
184
185 std::vector<HepVector>::iterator iter = m_fitresult.begin() + icounter*nKind*nBinPerCounter;
186 for( unsigned int k=0; k<nBinPerCounter-1; k++ ) {
187 toffset[k] = log((*(iter+k))[0]);
188 toffseterr[k] = log((*(iter+k))[0])*((*(iter+k))[1])/((*(iter+k))[0]);
189 tsigma[k] = (*(iter+k))[2];
190 tsigmaerr[k] = (*(iter+k))[3];
191 }
192
193 sprintf( gname1, "Q MPV-tofid-%i", icounter );
194 m_graphs.push_back( new TGraphErrors( nBinPerCounter-1, &rpos[0], &toffset[0], &rposerr[0], &toffseterr[0]) );
195 std::vector<TGraphErrors*>::iterator itgraph = m_graphs.end() - 1;
196 (*itgraph)->SetTitle(gname1);
197 (*itgraph)->SetMarkerSize(1.5);
198 (*itgraph)->SetMarkerStyle(20);
199 (*itgraph)->SetMarkerColor(2);
200
201 sprintf( gname2, "Q sigma-tofid-%i", icounter );
202 m_graphs.push_back( new TGraphErrors( nBinPerCounter-1, &rpos[0], &tsigma[0], &rposerr[0], &tsigmaerr[0]) );
203 itgraph = m_graphs.end() - 1;
204 (*itgraph)->SetTitle(gname2);
205 (*itgraph)->SetMarkerSize(1.5);
206 (*itgraph)->SetMarkerStyle(21);
207 (*itgraph)->SetMarkerColor(4);
208
209 return;
210}
211
212
213void calib_endcap_atten::fitGraph( unsigned int icounter ) {
214
215 TF1 *fsingleq = new TF1( "fsingleq", endcapQFunc, rbegin, rend, 3 );
216 fsingleq->SetLineColor(1);
217 fsingleq->SetLineWidth(1);
218 fsingleq->SetParameters( 6.5, 0.0, 0.0 );
219
220 std::vector<TGraphErrors*>::iterator itgraph = m_graphs.begin() + icounter*nGraphEcAtten;
221
222 (*itgraph)->Fit( "fsingleq", "Q", "", rbegin, rend );
223 X = HepVector( m_npar, 0 );
224 X[0] = fsingleq->GetParameter(0);
225 X[1] = fsingleq->GetParameter(1);
226 X[2] = fsingleq->GetParameter(2);
227 X[3] = 0.;
228 X[4] = 0.;
229 X[5] = 0.;
230
231 std::vector<HepVector>::iterator iter = m_result.begin() + icounter;
232 (*iter) = X;
233
234 return;
235}
236
237
238void calib_endcap_atten::updateData( Record* r, unsigned int icounter ) {
239 double rhit = r->zrhit();
240 double q = r->qleft();
241 double costheta = abs(r->theta());
242
243 double par[3];
244 std::vector<HepVector>::iterator iter = m_result.begin() + icounter;
245 for( unsigned int i=0; i<3; i++ ) {
246 par[i] = (*iter)[i];
247 }
248
249 par[0] = 0.;
250 double q0 = q*costheta/exp(endcapQFunc(&rhit,par));
251 r->setQ0( q0 );
252
253 return;
254}
255
256
257void calib_endcap_atten::fillRecordQ0( const Record* r, unsigned int icounter ) {
258 std::vector<TH1F*>::iterator iter = m_histograms.begin() + icounter*nKind*nBinPerCounter + nBinPerCounter - 1;
259 (*iter)->Fill( r->q0() );
260
261 return;
262}
263
264
265void calib_endcap_atten::fitHistogramQ0( unsigned int icounter ) {
266 TF1* ld = new TF1("ld", "landau");
267 ld->SetLineColor(2);
268 ld->SetLineWidth(1);
269
270 std::vector<TH1F*>::iterator iter1 = m_histograms.begin() + icounter*nKind*nBinPerCounter + nBinPerCounter - 1;
271 std::vector<HepVector>::iterator iter2 = m_fitresult.begin() + icounter*nKind*nBinPerCounter + nBinPerCounter - 1;
272 (*iter1)->Fit( ld, "Q");
273 (*iter2)[0] = ld->GetParameter(1);
274 (*iter2)[1] = ld->GetParError(1);
275 (*iter2)[2] = ld->GetParameter(2);
276 (*iter2)[3] = ld->GetParError(2);
277
278 return;
279}
280
281
282void calib_endcap_atten::fillGraphQ0() {
283 char gname1[256], gname2[256];
284
285 // fill graphs
286 // 2 canvas per counter,
287 // 1. Q0 MPV vs z
288 // 2. Q0 sigma vs z
289 std::vector<double> toffset, toffseterr;
290 std::vector<double> tsigma, tsigmaerr;
291 toffset.resize( NEndcap );
292 toffseterr.resize( NEndcap );
293 tsigma.resize( NEndcap );
294 tsigmaerr.resize( NEndcap );
295
296 unsigned int number = 0;
297 std::vector<HepVector>::iterator iter1 = m_fitresult.begin() + nBinPerCounter - 1;
298 std::vector<HepVector>::iterator iter2 = m_result.begin();
299 for( unsigned int i=0; i<NEndcap; i++ ) {
300 number = i*nBinPerCounter;
301 toffset[i] = (*(iter1+number))[0];
302 toffseterr[i] = (*(iter1+number))[1];
303 tsigma[i] = (*(iter1+number))[2];
304 tsigmaerr[i] = (*(iter1+number))[3];
305
306
307 (*(iter2+i))[3] = toffset[i]/toffset[0];
308 (*(iter2+i))[4] = toffset[i];
309 }
310
311 sprintf( gname1, "Q0 MPV vs TOF Counter Number" );
312 m_graphs.push_back( new TGraphErrors( NEndcap, &itofid[0], &toffset[0], &itofiderr[0], &toffseterr[0]) );
313 std::vector<TGraphErrors*>::iterator itgraph = m_graphs.end() - 1;
314 (*itgraph)->SetTitle(gname1);
315 (*itgraph)->SetMarkerSize(1.5);
316 (*itgraph)->SetMarkerStyle(20);
317 (*itgraph)->SetMarkerColor(2);
318
319 sprintf( gname2, "Q0 Sigma vs TOF Counter Number" );
320 m_graphs.push_back( new TGraphErrors( NEndcap, &itofid[0], &tsigma[0], &itofiderr[0], &tsigmaerr[0]) );
321 itgraph = m_graphs.end() - 1;
322 (*itgraph)->SetTitle(gname2);
323 (*itgraph)->SetMarkerSize(1.5);
324 (*itgraph)->SetMarkerStyle(21);
325 (*itgraph)->SetMarkerColor(4);
326
327 return;
328}
TTree * data
Double_t x[10]
EvtComplex exp(const EvtComplex &c)
double abs(const EvtComplex &c)
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
****INTEGER imax DOUBLE PRECISION m_pi *DOUBLE PRECISION m_amfin DOUBLE PRECISION m_Chfin DOUBLE PRECISION m_Xenph DOUBLE PRECISION m_sinw2 DOUBLE PRECISION m_GFermi DOUBLE PRECISION m_MfinMin DOUBLE PRECISION m_ta2 INTEGER m_out INTEGER m_KeyFSR INTEGER m_KeyQCD *COMMON c_Semalib $ !copy of input $ !CMS energy $ !beam mass $ !final mass $ !beam charge $ !final charge $ !smallest final mass $ !Z mass $ !Z width $ !EW mixing angle $ !Gmu Fermi $ alphaQED at q
Definition KKsem.h:33
const double rend
Definition TofCalibFit.h:19
const double rbegin
Definition TofCalibFit.h:18
std::vector< Record * > RecordSet
Definition TofDataSet.h:89
const unsigned int NEndcap
Definition TofDataSet.h:13
const int nEndcapAtten
const int nGraphEcAtten
const int nParEcAtten
double qleft() const
Definition TofDataSet.h:51
double q0() const
Definition TofDataSet.h:62
void setQ0(double q0)
Definition TofDataSet.h:69
double theta() const
Definition TofDataSet.h:59
double zrhit() const
Definition TofDataSet.h:55
string m_name
Definition TofCalibFit.h:52
std::vector< HepVector > m_result
Definition TofCalibFit.h:57
unsigned int nCanvas
Definition TofCalibFit.h:49
std::vector< TH1F * > m_histograms
Definition TofCalibFit.h:55
unsigned int nBinPerCounter
Definition TofCalibFit.h:43
std::vector< TGraphErrors * > m_graphs
Definition TofCalibFit.h:56
unsigned int nKind
Definition TofCalibFit.h:42
std::vector< unsigned int > nGraphPerCanvasPerCounter
Definition TofCalibFit.h:47
unsigned int nCanvasPerCounter
Definition TofCalibFit.h:46
unsigned int nHistPerCounter
Definition TofCalibFit.h:45
const string & name() const
Definition TofCalibFit.h:29
HepVector X
Definition TofCalibFit.h:53
std::vector< unsigned int > nGraphPerCanvas
Definition TofCalibFit.h:50
std::vector< string > CanvasName
Definition TofCalibFit.h:61
unsigned int nHistogram
Definition TofCalibFit.h:48
std::vector< string > CanvasPerCounterName
Definition TofCalibFit.h:60
calib_endcap_atten(const unsigned int nrbin)
void calculate(RecordSet *&data, unsigned int icounter)