BOSS 6.6.4.p01
BESIII Offline Software System
Loading...
Searching...
No Matches
calib_barrel_q0.cxx
Go to the documentation of this file.
2#include "TF1.h"
3
4using namespace CLHEP;
5
7
8 nKind = 1;
10
13
16 nHistogram = NBarrel*nKind*nBinPerCounter; // NO directory
17 nCanvas = 2;
18 CanvasName.push_back( static_cast<string>("Most Probable Value of Q0 vs TOF counter Number (Barrel part)") );
19 CanvasName.push_back( static_cast<string>("Sigma of Q0 vs TOF Counter Number (Barrel part)") );
20 nGraphPerCanvas.push_back(1);
21 nGraphPerCanvas.push_back(1);
22
23 m_name = string("calib_barrel_q0");
24
25 int numGraphs = 0;
26 std::vector<unsigned int>::iterator iter = nGraphPerCanvas.begin();
27 for( ; iter!=nGraphPerCanvas.end(); iter++ ) {
28 numGraphs = numGraphs + (*iter);
29 }
30 if( numGraphs != nGraphTotalQ0 ) {
31 cout << "tofcalgsec::calib_barrel_q0: the number of Graphs is NOT reasonable!!!" << endl;
32 exit(0);
33 }
34
35 const int qbin = 100;
36 const double qbegin = 0.0;
37 const double qend = 5000.0;
38
39 // histograms
40 char hname[256];
41 for( unsigned int i=0; i<NBarrel; i++ ) {
42 m_result.push_back( HepVector(nBarrelQ0,0) );
43
44 sprintf( hname, "Q0-tofid-%i", i );
45 m_histograms.push_back( new TH1F( hname, hname, qbin, qbegin, qend ) );
46
47 m_fitresult.push_back( HepVector(nParQ0,0) );
48 }
49
50 itofid.resize( NBarrel );
51 itofiderr.resize( NBarrel );
52 itofidstep = 1.0;
53 for( unsigned int i=0; i<NBarrel; i++ ) {
54 itofid[i] = i*1.0;
55 itofiderr[i] = 0.5;
56 }
57}
58
59
61 m_fitresult.clear();
62 itofid.clear();
63 itofiderr.clear();
64}
65
66
67void calib_barrel_q0::calculate( RecordSet*& data, unsigned int icounter ) {
68
69 std::cout << setiosflags(ios::left) << setw(10) << icounter << setw(8) << data->size() << setw(30) << name() << std::endl;
70
71 if( data->size() > 0 ) {
72 std::vector<Record*>::iterator iter = data->begin();
73 for( ; iter!=data->end(); iter++ ) {
74 fillRecord( (*iter), icounter );
75 }
76 }
77 fitHistogram( icounter );
78
79 if( icounter==(NBarrel-1) ) {
80 fillGraph();
81 fitGraph();
82 }
83
84 return;
85}
86
87
88void calib_barrel_q0::fillRecord( const Record* r, unsigned int icounter ) {
89 std::vector<TH1F*>::iterator iter = m_histograms.begin() + icounter;
90 (*iter)->Fill( r->q0() );
91 return;
92}
93
94
95void calib_barrel_q0::fitHistogram( unsigned int icounter ) {
96 TF1* ld = new TF1("ld", "landau");
97 ld->SetLineColor(2);
98 ld->SetLineWidth(1);
99
100 std::vector<TH1F*>::iterator iter1 = m_histograms.begin() + icounter;
101 std::vector<HepVector>::iterator iter2 = m_fitresult.begin() + icounter;
102 (*iter1)->Fit( ld, "Q");
103 (*iter2)[0] = ld->GetParameter(1);
104 (*iter2)[1] = ld->GetParError(1);
105 (*iter2)[2] = ld->GetParameter(2);
106 (*iter2)[3] = ld->GetParError(2);
107
108 return;
109}
110
111
112void calib_barrel_q0::fillGraph() {
113
114 std::vector<double> qmean, qmeanerr;
115 std::vector<double> qsig, qsigerr;
116 qmean.resize( NBarrel );
117 qmeanerr.resize( NBarrel );
118 qsig.resize( NBarrel );
119 qsigerr.resize( NBarrel );
120
121 std::vector<HepVector>::iterator iter = m_fitresult.begin();
122 for( unsigned int i=0; i<NBarrel; i++, iter++ ) {
123 qmean[i] = (*iter)[0];
124 qmeanerr[i] = (*iter)[1];
125 qsig[i] = (*iter)[2];
126 qsigerr[i] = (*iter)[3];
127 }
128
129 TGraphErrors* graph1 = new TGraphErrors( NBarrel, &itofid[0], &qmean[0], &itofiderr[0], &qmeanerr[0]);
130 graph1->SetTitle( CanvasName[0].c_str() );
131 graph1->SetMarkerSize(1.5);
132 graph1->SetMarkerStyle(20);
133 graph1->SetMarkerColor(2);
134 m_graphs.push_back( graph1 );
135
136 TGraphErrors* graph2 = new TGraphErrors( NBarrel, &itofid[0], &qsig[0], &itofiderr[0], &qsigerr[0]);
137 graph2->SetTitle( CanvasName[1].c_str() );
138 graph2->SetMarkerSize(1.5);
139 graph2->SetMarkerStyle(20);
140 graph2->SetMarkerColor(4);
141 m_graphs.push_back( graph2 );
142
143 return;
144}
145
146
147void calib_barrel_q0::fitGraph() {
148 unsigned int number = 0;
149 std::vector<HepVector>::iterator iter1 = m_result.begin();
150 std::vector<HepVector>::iterator iter2 = m_fitresult.begin();
151 for( ; iter1!=m_result.end(); iter1++, number++ ) {
152 (*iter1)[0] = (*(iter2+number))[0]/ (*iter2)[0];
153 (*iter1)[1] = (*(iter2+number))[0];
154 (*iter1)[2] = (*(iter2+number))[2];
155 }
156
157 return;
158}
TTree * data
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
std::vector< Record * > RecordSet
Definition: TofDataSet.h:89
const unsigned int NBarrel
Definition: TofDataSet.h:12
const int nBarrelQ0
const int nParQ0
const int nGraphTotalQ0
double q0() const
Definition: TofDataSet.h:62
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
unsigned int nCanvasPerCounter
Definition: TofCalibFit.h:46
unsigned int nHistPerCounter
Definition: TofCalibFit.h:45
const string & name() const
Definition: TofCalibFit.h:29
std::vector< unsigned int > nGraphPerCanvas
Definition: TofCalibFit.h:50
std::vector< string > CanvasName
Definition: TofCalibFit.h:61
unsigned int nHistogram
Definition: TofCalibFit.h:48
void calculate(RecordSet *&data, unsigned int icounter)
char * c_str(Index i)
Definition: EvtCyclic3.cc:252