BOSS 7.0.4
BESIII Offline Software System
Loading...
Searching...
No Matches
TofCalibFit.cxx
Go to the documentation of this file.
1//**************************************************************************
2// several value:
3// npar m_npar: the parameter numbers of result
4// nCounter: counter number of TOF, 176 for barrel, 96 for endcap
5// nKind: number of values, barrel_q0 1, barrel_sigma 4, endcap_q 1
6// nBinPerCounter: number of bins divide one counter ( z(r) bins number)
7// total: nCounter*nKind*nBinPerCounter
8//**************************************************************************
9#include "tofcalgsec/TofCalibFit.h"
10#include "TDirectory.h"
11#include "TFile.h"
12#include "TStyle.h"
13#include <iostream>
14#include <ios>
15#include <fstream>
16
17TofCalibFit::TofCalibFit( bool isbarrel, const int npar ) {
18 if( isbarrel ) { nCounter = NBarrel; }
19 else { nCounter = NEndcap; }
20
21 nKind = 0;
23
26 nHistogram = 0;
27 nCanvas = 0;
28
29 m_npar = npar;
30 X = HepVector(m_npar,0);
31 m_name = string("calibfit");
32
33 m_tcorrelation = HepVector( nBarrelCommon,0);
34}
35
36
38 std::vector<TH1F*>::iterator iter1 = m_histograms.begin();
39 for( ; iter1 != m_histograms.end(); iter1++ ) {
40 delete (*iter1);
41 }
42 m_histograms.clear();
43 std::vector<TH1F*>::iterator iter2 = m_graphs.begin();
44 for( ; iter2 != m_graphs.end(); iter2++ ) {
45 delete (*iter2);
46 }
47 m_graphs.clear();
48 m_result.clear();
49}
50
51
52void TofCalibFit::fillTxt( const char* file ) {
53 std::ofstream out(file,ios::out);
54 if( out ) {
55 std::vector<HepVector>::iterator it;
56 for( it=m_result.begin(); it!=m_result.end(); it++ ) {
57 for( int i=0; i<(*it).num_row(); i++ ) {
58 out << (*it)[i] << " ";
59 }
60 out << std::endl;
61 }
62 out.close();
63 }
64 else{
65 cerr << "error when open file " << file << " for write in " << name() << "::fillTxt()" << std::endl;
66 cout << "print all parameters to srceen: in total " << m_result.size() << " items" << std::endl;
67 std::vector<HepVector>::iterator it;
68 for( it=m_result.begin(); it!=m_result.end(); it++ ) {
69 for( int i=0; i<(*it).num_row(); i++ ) {
70 cout << (*it)[i] << " ";
71 }
72 cout << std::endl;
73 }
74 }
75
76 return;
77}
78
79
80void TofCalibFit::fillRoot( const char* file ) {
81
82 unsigned int nhist = m_histograms.size();
83 if( nhist != (nCounter*nHistPerCounter + nHistogram) ) {
84 std::cout<<" tofcalgsec::TofCalibFit:" << m_name << ": the number of histograms is NOT same as the number of histograms saved!" << " nhist=" << nhist << " calculated=" << (nCounter*nHistPerCounter + nHistogram) << " nCounter=" << nCounter << " nHistPerCounter=" << nHistPerCounter << " nHistogram=" << nHistogram << std::endl;
85 exit(0);
86 }
87
88 unsigned int numgraph1 = 0;
89 unsigned int numgraph2 = 0;
90 if( nCanvasPerCounter!=0 ) {
91 std::vector<unsigned int>::iterator iter = nGraphPerCanvasPerCounter.begin();
92 for( ; iter!=nGraphPerCanvasPerCounter.end(); iter++ ) {
93 numgraph1 = numgraph1 + (*iter);
94 }
95 }
96 if( nCanvas!=0 ) {
97 std::vector<unsigned int>::iterator iter = nGraphPerCanvas.begin();
98 for( ; iter!=nGraphPerCanvas.end(); iter++ ) {
99 numgraph2 = numgraph2 + (*iter);
100 }
101 }
102 unsigned int ngraph = m_graphs.size();
103 if( ngraph != ( nCounter*numgraph1+numgraph2) ) {
104 std::cout<<" tofcalgsec::TofCalibFit:"<< m_name << ": the number of graphs is NOT same as the number of graphs saved!" << " ngraph=" << ngraph << " nCounter=" << nCounter << " numgraph1=" << numgraph1 << " numgraph2=" << numgraph2 << std::endl;
105 exit(0);
106 }
107
108 TFile f(file,"RECREATE");
109
110 gStyle->SetOptStat(2211);
111 gStyle->SetOptFit(1111);
112 gStyle->SetLabelSize(0.03,"x");
113 gStyle->SetLabelSize(0.03,"y");
114
115 char dirname[256];
116 char canvasname[256];
118 std::vector<TH1F*>::iterator iter1 = m_histograms.begin();
119 std::vector<TH1F*>::iterator iter2 = m_graphs.begin();
120 for( unsigned int i=0; i<nCounter; i++ ) {
121 sprintf( dirname, "tofid%i", i );
122 TDirectory* cdresult = f.mkdir( dirname );
123 cdresult->cd();
124
125 for( unsigned int j=0; j<nHistPerCounter; j++ ) {
126 (*(iter1+j))->Write();
127 }
128 iter1 = iter1 + nHistPerCounter;
129
130 for( unsigned int j=0; j<nCanvasPerCounter; j++ ) {
131 std::vector<string>::iterator it1 = CanvasPerCounterName.begin() + j;
132 std::vector<unsigned int>::iterator it2 = nGraphPerCanvasPerCounter.begin() + j;
133 sprintf( canvasname, "%s-tofid-%i", (*it1).c_str(), i );
134 TCanvas* c1 = new TCanvas( canvasname, canvasname, 1);
135 c1->SetFillColor(10);
136 for( unsigned int k=0; k<(*it2); k++ ) {
137 if( k==0 ) {
138 (*(iter2+k))->Draw("E");
139 }
140 else {
141 (*(iter2+k))->Draw("Esame");
142 }
143 }
144 iter2 = iter2 + (*it2);
145 c1->Write();
146 }
147 }
148 }
149
150 if( nHistogram>0 || nCanvas>0 ) {
151 sprintf( dirname, "summary" );
152 TDirectory* cdresult = f.mkdir( dirname );
153 cdresult->cd();
154
155 std::vector<TH1F*>::iterator iter1 = m_histograms.begin() + nCounter*nHistPerCounter;
156 std::vector<TH1F*>::iterator iter2 = m_graphs.begin() + nCounter*numgraph1;
157 for( ; iter1 != m_histograms.end(); iter1++ ) {
158 (*iter1)->Write();
159 }
160
161 for( unsigned int j=0; j<nCanvas; j++ ) {
162 std::vector<string>::iterator it1 = CanvasName.begin() + j;
163 std::vector<unsigned int>::iterator it2 = nGraphPerCanvas.begin() + j;
164 sprintf( canvasname, (*it1).c_str() );
165 TCanvas* c1 = new TCanvas( canvasname, canvasname, 1);
166 c1->SetFillColor(10);
167 for( unsigned int k=0; k<(*it2); k++ ) {
168 if( k==0 ) {
169 (*(iter2+k))->Draw("E");
170 }
171 else {
172 (*(iter2+k))->Draw("Esame");
173 }
174 }
175 iter2 = iter2 + (*it2);
176 c1->Write();
177 }
178
179 }
180
181 f.Close();
182
183 return;
184}
char * file
Definition: DQA_TO_DB.cxx:15
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
void fillRoot(const char *file)
Definition: TofCalibFit.cxx:80
TofCalibFit(bool isbarrel, const int npar)
Definition: TofCalibFit.cxx:17
void fillTxt(const char *file)
Definition: TofCalibFit.cxx:52
std::vector< unsigned int > nGraphPerCanvasPerCounter
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)
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")