BOSS 7.0.2
BESIII Offline Software System
Loading...
Searching...
No Matches
calib_endcap_sigma.cxx
Go to the documentation of this file.
1#include "tofcalgsec/calib_endcap_sigma.h"
2#include "TF1.h"
3
5
6 nKind = 1; // 0:tleft
7 nBinPerCounter = nrbin;
8
11 CanvasPerCounterName.push_back( static_cast<string>("Endcap-offset") );
12 CanvasPerCounterName.push_back( static_cast<string>("Endcap-sigma") );
13
14 nGraphPerCanvasPerCounter.push_back(1);
15 nGraphPerCanvasPerCounter.push_back(1);
16
17 nHistogram = 0;
18 nCanvas = 0;
19
20 int numGraphs = 0;
21 std::vector<unsigned int>::iterator iter = nGraphPerCanvasPerCounter.begin();
22 for( ; iter!=nGraphPerCanvasPerCounter.end(); iter++ ) {
23 numGraphs = numGraphs + (*iter);
24 }
25 if( numGraphs != nGraphEcSigma ) {
26 cout << "tofcalgsec::calib_endcap_sigma: the number of Graphs is NOT reasonable!!!" << endl;
27 exit(0);
28 }
29
30 m_name = string("calib_endcap_sigma");
31
32 const int tbin = 150;
33 const double tbegin = -1.5;
34 const double tend = 1.5;
35
36 // histograms per counter
37 char hname[256];
38 for( unsigned int i=0; i<NEndcap; i++ ) {
39 m_result.push_back( HepVector(nEndcapSigma,0) );
40 for( unsigned int j=0; j<nKind; j++ ) {
41 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
42 sprintf( hname, "tleft-tofid%i-r%i", i, k);
43 m_histograms.push_back( new TH1F( hname, hname, tbin, tbegin, tend ) );
44
45 m_fitresult.push_back( HepVector(nParEcSigma,0) );
46 }
47 }
48 }
49
50 rpos.resize( nBinPerCounter );
51 rposerr.resize( nBinPerCounter );
52 rstep = ( rend - rbegin )/nBinPerCounter;
53 for( unsigned int i=0; i<nBinPerCounter; i++ ) {
54 rpos[i] = rbegin + ( i+0.5 )*rstep;
55 rposerr[i] = 0.5*rstep;
56 }
57
58}
59
60
62 m_fitresult.clear();
63 rpos.clear();
64 rposerr.clear();
65}
66
67
68void calib_endcap_sigma::calculate( RecordSet*& data, unsigned int icounter ) {
69
70 std::cout << setiosflags(ios::left) << setw(10) << icounter << setw(8) << data->size() << setw(30) << name() << std::endl;
71
72 if( data->size() > 0 ) {
73 std::vector<Record*>::iterator iter = data->begin();
74 for( ; iter!=data->end(); iter++ ) {
75 fillRecord( (*iter), icounter );
76 }
77 }
78 fitHistogram( icounter );
79 fillGraph( icounter );
80 fitGraph( icounter );
81
82 return;
83}
84
85
86void calib_endcap_sigma::fillRecord( const Record* r, unsigned int icounter ) {
87
88 double rhit = r->zrhit();
89 if( rhit<rbegin || rhit>rend ) return;
90 int rbin = static_cast<int>((rhit-rbegin)/rstep);
91 if( rbin<0 ) { rbin = 0; }
92 else if( rbin>static_cast<int>(nBinPerCounter-1) ) {
93 cout << "tofcalgsec::calib_endcap_sigma:fillRecord: rhit is out of range, rhit=" << rhit << " rbin=" << rbin << endl;
94 return;
95 }
96
97 std::vector<TH1F*>::iterator iter = m_histograms.begin() + icounter*nKind*nBinPerCounter + rbin;
98 (*iter)->Fill( r->tleft() );
99
100 return;
101}
102
103
104void calib_endcap_sigma::fitHistogram( unsigned int icounter ) {
105 TF1* g = new TF1("g", "gaus");
106 g->SetLineColor(2);
107 g->SetLineWidth(1);
108
109 std::vector<TH1F*>::iterator iter1 = m_histograms.begin() + icounter*nKind*nBinPerCounter;
110 std::vector<HepVector>::iterator iter2 = m_fitresult.begin() + icounter*nKind*nBinPerCounter;
111 for( unsigned int j=0; j<nBinPerCounter; j++, iter1++, iter2++ ) {
112 (*iter1)->Fit( g, "Q");
113 (*iter2)[0] = g->GetParameter(1);
114 (*iter2)[1] = g->GetParError(1);
115 (*iter2)[2] = g->GetParameter(2);
116 (*iter2)[3] = g->GetParError(2);
117 }
118
119 return;
120
121}
122
123
124void calib_endcap_sigma::fillGraph( unsigned int icounter ) {
125
126 char gname1[256], gname2[256];
127
128 // fill graphs per counter
129 // 4 canvas per counter,
130 // 1. offset of tleft vs z
131 // 2. sigma of tleft vs z
132 std::vector<double> toffset, toffseterr;
133 std::vector<double> tsigma, tsigmaerr;
134 toffset.resize( nBinPerCounter );
135 toffseterr.resize( nBinPerCounter );
136 tsigma.resize( nBinPerCounter );
137 tsigmaerr.resize( nBinPerCounter );
138
139 std::vector<HepVector>::iterator iter = m_fitresult.begin() + icounter*nKind*nBinPerCounter;
140 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
141 toffset[k] = (*(iter+k))[0];
142 toffseterr[k] = (*(iter+k))[1];
143 tsigma[k] = (*(iter+k))[2];
144 tsigmaerr[k] = (*(iter+k))[3];
145 }
146
147 sprintf( gname1, "endcap-offset-tofid-%i", icounter );
148 m_graphs.push_back( new TH1F( gname1, gname1, nBinPerCounter, rbegin, rend ) );
149 std::vector<TH1F*>::iterator itgraph = m_graphs.end() - 1;
150 (*itgraph)->SetMarkerSize(1.5);
151 (*itgraph)->SetMarkerStyle(20);
152 (*itgraph)->SetMarkerColor(2);
153 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
154 (*itgraph)->SetBinContent( k+1, toffset[k] );
155 (*itgraph)->SetBinError( k+1, toffseterr[k] );
156 }
157
158 sprintf( gname2, "endcap-sigma-tofid-%i", icounter );
159 m_graphs.push_back( new TH1F( gname2, gname2, nBinPerCounter, rbegin, rend ) );
160 itgraph = m_graphs.end() - 1;
161 (*itgraph)->SetMarkerSize(1.5);
162 (*itgraph)->SetMarkerStyle(21);
163 (*itgraph)->SetMarkerColor(4);
164 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
165 (*itgraph)->SetBinContent( k+1, tsigma[k] );
166 (*itgraph)->SetBinError( k+1, tsigmaerr[k] );
167 }
168
169 return;
170}
171
172
173void calib_endcap_sigma::fitGraph( unsigned int icounter ) {
174
175 TF1* p2 = new TF1("p2", "pol2", rbegin, rend );
176
177 // std::vector<unsigned int>::iterator itnumber = nGraphPerCanvasPerCounter.begin();
178 std::vector<TH1F*>::iterator itgraph = m_graphs.begin() + icounter*nGraphEcSigma + 1;
179
180 (*itgraph)->Fit( "p2", "Q" );
181 X = HepVector( m_npar, 0 );
182 X[0] = p2->GetParameter(0);
183 X[1] = p2->GetParameter(1);
184 X[2] = p2->GetParameter(2);
185
186 std::vector<HepVector>::iterator iter = m_result.begin() + icounter;
187 (*iter) = X;
188
189 return;
190}
TTree * data
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
std::vector< Record * > RecordSet
std::vector< unsigned int > nGraphPerCanvasPerCounter
calib_endcap_sigma(const unsigned int nrbin)
void calculate(RecordSet *&data, unsigned int icounter)