BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxCalibCostheta Class Reference

#include <DedxCalibCostheta.h>

Inheritance diagram for DedxCalibCostheta:

Public Member Functions

 DedxCalibCostheta (const std::string &name, ISvcLocator *pSvcLocator)
 ~DedxCalibCostheta ()
void initializing ()
void BookHists ()
void genNtuple ()
void FillHists ()
void AnalyseHists ()
void WriteHists ()
Public Member Functions inherited from DedxCalib
 DedxCalib (const std::string &name, ISvcLocator *pSvcLocator)
 ~DedxCalib ()
StatusCode initialize ()
StatusCode execute ()
StatusCode finalize ()

Additional Inherited Members

Protected Member Functions inherited from DedxCalib
double cal_dedx_bitrunc (float truncate, std::vector< double > phlist)
double cal_dedx (float truncate, std::vector< double > phlist)
void getCurvePar ()
void set_dEdx (int landau, float dEdx, int trkalg, int runflag, int dedxhit_use, float ptrk, float theta, float pl_rp, int vflag[3], double t0)
void ReadRecFileList ()
Protected Attributes inherited from DedxCalib
IMdcGeomSvcgeosvc
IDedxCorrecSvcexsvc
float truncate
vector< double > Curve_Para
vector< double > Sigma_Para
int vFlag [3]
double m_dedx_exp [5]
double m_ex_sigma [5]
double m_pid_prob [5]
double m_chi_ex [5]
vector< string > m_recFileVector
int ParticleFlag
int m_calibflag
int m_phShapeFlag
std::string m_eventType
std::string m_recFileList
std::string m_rootfile
std::string m_curvefile

Detailed Description

Definition at line 12 of file DedxCalibCostheta.h.

Constructor & Destructor Documentation

◆ DedxCalibCostheta()

DedxCalibCostheta::DedxCalibCostheta ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 21 of file DedxCalibCostheta.cxx.

22 : DedxCalib( name, pSvcLocator ) {
23 declareProperty( "Debug", m_debug = false );
24 declareProperty( "Ave", m_ave = true ); // unweighted average of positive and negative
25 declareProperty( "DebugMin", m_debug_min = 2 );
26 declareProperty( "DebugMax", m_debug_max = 2 );
27 declareProperty( "PMin", pMin = 0.2 );
28 declareProperty( "PMax", pMax = 2.3 );
29}
DedxCalib(const std::string &name, ISvcLocator *pSvcLocator)
Definition DedxCalib.cxx:12

◆ ~DedxCalibCostheta()

DedxCalibCostheta::~DedxCalibCostheta ( )
inline

Definition at line 15 of file DedxCalibCostheta.h.

15{};

Member Function Documentation

◆ AnalyseHists()

void DedxCalibCostheta::AnalyseHists ( )
virtual

Implements DedxCalib.

Definition at line 175 of file DedxCalibCostheta.cxx.

175 {
176 MsgStream log( msgSvc(), name() );
177 log << MSG::INFO << "DedxCalibCostheta::AnalyseHists()" << endmsg;
178
179 gStyle->SetOptFit( 1111 );
180 for ( int i = 0; i < NumBinCostheta; i++ )
181 {
182 if ( m_debug ) cout << "num of bin " << i << endl;
183 if ( m_costheta[i]->GetEntries() > 100 ) m_costheta[i]->Fit( "gaus", "Q" );
184 if ( m_pos_costheta[i]->GetEntries() > 100 ) m_pos_costheta[i]->Fit( "gaus", "Q" );
185 if ( m_neg_costheta[i]->GetEntries() > 100 ) m_neg_costheta[i]->Fit( "gaus", "Q" );
186 if ( m_chi[i]->GetEntries() > 100 ) m_chi[i]->Fit( "gaus", "Q" );
187 if ( m_pos_chi[i]->GetEntries() > 100 ) m_pos_chi[i]->Fit( "gaus", "Q" );
188 if ( m_neg_chi[i]->GetEntries() > 100 ) m_neg_chi[i]->Fit( "gaus", "Q" );
189 }
190}
IMessageSvc * msgSvc()

◆ BookHists()

void DedxCalibCostheta::BookHists ( )
virtual

Implements DedxCalib.

Definition at line 31 of file DedxCalibCostheta.cxx.

31 {
32 MsgStream log( msgSvc(), name() );
33 log << MSG::INFO << "DedxCalibCostheta::BookHists()" << endmsg;
34
36
37 m_costheta = new TH1F*[NumBinCostheta];
38 m_pos_costheta = new TH1F*[NumBinCostheta];
39 m_neg_costheta = new TH1F*[NumBinCostheta];
40 m_chi = new TH1F*[NumBinCostheta];
41 m_pos_chi = new TH1F*[NumBinCostheta];
42 m_neg_chi = new TH1F*[NumBinCostheta];
43 stringstream hlabel;
44 for ( int i = 0; i < NumBinCostheta; i++ )
45 {
46 hlabel.str( "" );
47 hlabel << "dEdx_costheta" << i;
48 m_costheta[i] = new TH1F( hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins,
50 hlabel.str( "" );
51 hlabel << "pos_dEdx_costheta" << i;
52 m_pos_costheta[i] = new TH1F( hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins,
54 hlabel.str( "" );
55 hlabel << "neg_dEdx_costheta" << i;
56 m_neg_costheta[i] = new TH1F( hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins,
58 hlabel.str( "" );
59 hlabel << "chi_costheta" << i;
60 m_chi[i] = new TH1F( hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinChiValue,
62 hlabel.str( "" );
63 hlabel << "pos_chi_costheta" << i;
64 m_pos_chi[i] = new TH1F( hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins,
66 hlabel.str( "" );
67 hlabel << "neg_chi_costheta" << i;
68 m_neg_chi[i] = new TH1F( hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins,
70 }
71 hlabel.str( "" );
72 hlabel << "dEdxVsCostheta";
73 m_dEdxVsCostheta = new TH1F( hlabel.str().c_str(), "dEdx vs costheta", NumBinCostheta,
75 m_dEdxVsCostheta->GetXaxis()->SetTitle( "cos#theta" );
76 m_dEdxVsCostheta->GetYaxis()->SetTitle( "dE/dx" );
77 hlabel.str( "" );
78 hlabel << "pos_dEdxVsCostheta";
79 m_pos_dEdxVsCostheta = new TH1F( hlabel.str().c_str(), "positron dEdx vs costheta",
81 m_pos_dEdxVsCostheta->GetXaxis()->SetTitle( "cos#theta" );
82 m_pos_dEdxVsCostheta->GetYaxis()->SetTitle( "dE/dx" );
83 hlabel.str( "" );
84 hlabel << "neg_dEdxVsCostheta";
85 m_neg_dEdxVsCostheta = new TH1F( hlabel.str().c_str(), "electron dEdx vs costheta",
87 m_neg_dEdxVsCostheta->GetXaxis()->SetTitle( "cos#theta" );
88 m_neg_dEdxVsCostheta->GetYaxis()->SetTitle( "dE/dx" );
89
90 Vec_dedx.clear();
91 Vec_costheta.clear();
92}
#define MinHistValue1
#define MaxChiValue
#define MaxHistValue1
#define MinChiValue
void ReadRecFileList()
Definition DedxCalib.cxx:84

◆ FillHists()

void DedxCalibCostheta::FillHists ( )
virtual

Implements DedxCalib.

Definition at line 94 of file DedxCalibCostheta.cxx.

94 {
95 MsgStream log( msgSvc(), name() );
96 log << MSG::INFO << "DedxCalibCostheta::FillHists()" << endmsg;
97
98 TFile* f;
99 TTree* n103;
100 string runlist;
101
102 int ndedxhit = 0;
103 double dEdx[100] = { 0 }, pathlength[100] = { 0 }, wid[100] = { 0 }, layid[100] = { 0 },
104 dd_in[100] = { 0 }, eangle[100] = { 0 }, zhit[100] = { 0 };
105 double dedx = 0;
106 float runNO = 0, evtNO = 0, runFlag = 0, costheta = 0, tes = 0, charge = 0, recalg = 0,
107 ptrk = 0, chidedx = 0;
108 int usedhit = 0;
109 vector<double> phlist;
110 cout << "m_recFileVector.size()= " << m_recFileVector.size() << endl;
111 for ( unsigned int i = 0; i < m_recFileVector.size(); i++ )
112 {
113 int m_current_trks( 0 );
114 runlist = m_recFileVector[i];
115 f = new TFile( runlist.c_str() );
116 n103 = (TTree*)f->Get( "n103" );
117 n103->SetBranchAddress( "ndedxhit", &ndedxhit );
118 n103->SetBranchAddress( "dEdx_hit", dEdx );
119 n103->SetBranchAddress( "pathlength_hit", pathlength );
120 n103->SetBranchAddress( "wid_hit", wid );
121 n103->SetBranchAddress( "layid_hit", layid );
122 n103->SetBranchAddress( "dd_in_hit", dd_in );
123 n103->SetBranchAddress( "eangle_hit", eangle );
124 n103->SetBranchAddress( "zhit_hit", zhit );
125 n103->SetBranchAddress( "runNO", &runNO );
126 n103->SetBranchAddress( "evtNO", &evtNO );
127 n103->SetBranchAddress( "runFlag", &runFlag );
128 n103->SetBranchAddress( "costheta", &costheta );
129 n103->SetBranchAddress( "t0", &tes );
130 n103->SetBranchAddress( "charge", &charge );
131 n103->SetBranchAddress( "recalg", &recalg );
132 n103->SetBranchAddress( "ndedxhit", &ndedxhit );
133 n103->SetBranchAddress( "ptrk", &ptrk );
134 n103->SetBranchAddress( "chidedx_e", &chidedx );
135 for ( int j = 0; j < n103->GetEntries(); j++ )
136 {
137 m_current_trks++;
138 if ( m_current_trks > Size ) break;
139 phlist.clear();
140 n103->GetEntry( j );
141 if ( costheta > CosthetaMax || costheta < CosthetaMin ) continue;
142 if ( tes > 1400 ) continue;
143 if ( ptrk > pMax || ptrk < pMin ) continue;
144 for ( int k = 0; k < ndedxhit; k++ )
145 {
146 dEdx[k] = exsvc->StandardHitCorrec( 0, runFlag, 2, runNO, evtNO, pathlength[k], wid[k],
147 layid[k], dEdx[k], dd_in[k], eangle[k], costheta );
148 phlist.push_back( dEdx[k] );
149 }
150 dedx = cal_dedx_bitrunc( truncate, phlist );
151 int iCos = (Int_t)floor( ( costheta - CosthetaMin ) / StepCostheta );
152 double pre_dedx = dedx;
153 if ( m_debug && iCos >= m_debug_min && iCos <= m_debug_max )
154 cout << "before cor, dedx " << pre_dedx << " with cos(theta) " << costheta << endl;
155 dedx = exsvc->StandardTrackCorrec( 0, runFlag, 2, runNO, evtNO, dedx, costheta, tes );
156 if ( m_debug && iCos >= m_debug_min && iCos <= m_debug_max )
157 cout << "after cor, dedx " << dedx << " with gain " << pre_dedx / dedx << endl;
158 m_costheta[iCos]->Fill( dedx );
159 if ( charge > 0 ) m_pos_costheta[iCos]->Fill( dedx );
160 if ( charge < 0 ) m_neg_costheta[iCos]->Fill( dedx );
161
162 usedhit = ndedxhit * truncate;
163 set_dEdx( 1, dedx, recalg, runFlag, usedhit, ptrk, acos( costheta ), 1.5, vFlag, tes );
164 double chi = m_chi_ex[ParticleFlag];
165 m_chi[iCos]->Fill( chi );
166 if ( charge > 0 ) m_pos_chi[iCos]->Fill( chi );
167 if ( charge < 0 ) m_neg_chi[iCos]->Fill( chi );
168
169 Vec_dedx.push_back( dedx );
170 Vec_costheta.push_back( costheta );
171 }
172 }
173}
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
#define Size
const double StepCostheta
double m_chi_ex[5]
Definition DedxCalib.h:48
void set_dEdx(int landau, float dEdx, int trkalg, int runflag, int dedxhit_use, float ptrk, float theta, float pl_rp, int vflag[3], double t0)
int ParticleFlag
Definition DedxCalib.h:52
int vFlag[3]
Definition DedxCalib.h:41
IDedxCorrecSvc * exsvc
Definition DedxCalib.h:29
double cal_dedx_bitrunc(float truncate, std::vector< double > phlist)
vector< string > m_recFileVector
Definition DedxCalib.h:50
float truncate
Definition DedxCalib.h:32

◆ genNtuple()

void DedxCalibCostheta::genNtuple ( )
inlinevirtual

Implements DedxCalib.

Definition at line 18 of file DedxCalibCostheta.h.

18{}

◆ initializing()

void DedxCalibCostheta::initializing ( )
inlinevirtual

Implements DedxCalib.

Definition at line 16 of file DedxCalibCostheta.h.

16{}

◆ WriteHists()

void DedxCalibCostheta::WriteHists ( )
virtual

Implements DedxCalib.

Definition at line 192 of file DedxCalibCostheta.cxx.

192 {
193 MsgStream log( msgSvc(), name() );
194 log << MSG::INFO << "DedxCalibCostheta::WriteHists()" << endmsg;
195
196 double chientryNo[NumBinCostheta] = { 0 }, chimean[NumBinCostheta] = { 0 },
197 chimeanerr[NumBinCostheta] = { 0 }, chisigma[NumBinCostheta] = { 0 },
198 chisq[NumBinCostheta] = { 0 };
199 double pos_chientryNo[NumBinCostheta] = { 0 }, pos_chimean[NumBinCostheta] = { 0 },
200 pos_chimeanerr[NumBinCostheta] = { 0 }, pos_chisigma[NumBinCostheta] = { 0 },
201 pos_chisq[NumBinCostheta] = { 0 };
202 double neg_chientryNo[NumBinCostheta] = { 0 }, neg_chimean[NumBinCostheta] = { 0 },
203 neg_chimeanerr[NumBinCostheta] = { 0 }, neg_chisigma[NumBinCostheta] = { 0 },
204 neg_chisq[NumBinCostheta] = { 0 };
205 double fitentryNo[NumBinCostheta] = { 0 }, fitmean[NumBinCostheta] = { 0 },
206 fitmeanerr[NumBinCostheta] = { 0 }, fitsigma[NumBinCostheta] = { 0 },
207 gain[NumBinCostheta] = { 0 }, fitchisq[NumBinCostheta] = { 0 };
208 double pos_fitentryNo[NumBinCostheta] = { 0 }, pos_fitmean[NumBinCostheta] = { 0 },
209 pos_fitmeanerr[NumBinCostheta] = { 0 }, pos_fitsigma[NumBinCostheta] = { 0 },
210 pos_gain[NumBinCostheta] = { 0 }, pos_fitchisq[NumBinCostheta] = { 0 };
211 double neg_fitentryNo[NumBinCostheta] = { 0 }, neg_fitmean[NumBinCostheta] = { 0 },
212 neg_fitmeanerr[NumBinCostheta] = { 0 }, neg_fitsigma[NumBinCostheta] = { 0 },
213 neg_gain[NumBinCostheta] = { 0 }, neg_fitchisq[NumBinCostheta] = { 0 };
214 double cosBin[NumBinCostheta] = { 0 };
215
216 for ( int i = 0; i < NumBinCostheta; i++ )
217 {
218 cosBin[i] = ( i + 0.5 ) * StepCostheta + CosthetaMin;
219
220 chientryNo[i] = m_chi[i]->GetEntries();
221 if ( m_debug )
222 cout << "get results at " << i << " bin with chi entries " << chientryNo[i] << endl;
223 if ( m_chi[i]->GetFunction( "gaus" ) )
224 {
225 chimean[i] = m_chi[i]->GetFunction( "gaus" )->GetParameter( 1 );
226 chimeanerr[i] = m_chi[i]->GetFunction( "gaus" )->GetParError( 1 );
227 chisigma[i] = m_chi[i]->GetFunction( "gaus" )->GetParameter( 2 );
228 chisq[i] = ( m_chi[i]->GetFunction( "gaus" )->GetChisquare() ) /
229 ( m_chi[i]->GetFunction( "gaus" )->GetNDF() );
230 }
231 pos_chientryNo[i] = m_pos_chi[i]->GetEntries();
232 if ( m_debug )
233 cout << "get results at " << i << " bin with pos_chi entries " << pos_chientryNo[i]
234 << endl;
235 if ( m_pos_chi[i]->GetFunction( "gaus" ) )
236 {
237 pos_chimean[i] = m_pos_chi[i]->GetFunction( "gaus" )->GetParameter( 1 );
238 pos_chimeanerr[i] = m_pos_chi[i]->GetFunction( "gaus" )->GetParError( 1 );
239 pos_chisigma[i] = m_pos_chi[i]->GetFunction( "gaus" )->GetParameter( 2 );
240 pos_chisq[i] = ( m_pos_chi[i]->GetFunction( "gaus" )->GetChisquare() ) /
241 ( m_pos_chi[i]->GetFunction( "gaus" )->GetNDF() );
242 }
243 neg_chientryNo[i] = m_neg_chi[i]->GetEntries();
244 if ( m_debug )
245 cout << "get results at " << i << " bin with neg_chi entries " << neg_chientryNo[i]
246 << endl;
247 if ( m_neg_chi[i]->GetFunction( "gaus" ) )
248 {
249 neg_chimean[i] = m_neg_chi[i]->GetFunction( "gaus" )->GetParameter( 1 );
250 neg_chimeanerr[i] = m_neg_chi[i]->GetFunction( "gaus" )->GetParError( 1 );
251 neg_chisigma[i] = m_neg_chi[i]->GetFunction( "gaus" )->GetParameter( 2 );
252 neg_chisq[i] = ( m_neg_chi[i]->GetFunction( "gaus" )->GetChisquare() ) /
253 ( m_neg_chi[i]->GetFunction( "gaus" )->GetNDF() );
254 }
255
256 fitentryNo[i] = m_costheta[i]->GetEntries();
257 if ( m_debug )
258 cout << "get results at " << i << " bin with fit entries " << fitentryNo[i] << endl;
259 if ( m_costheta[i]->GetFunction( "gaus" ) )
260 {
261 fitmean[i] = m_costheta[i]->GetFunction( "gaus" )->GetParameter( 1 );
262 fitmeanerr[i] = m_costheta[i]->GetFunction( "gaus" )->GetParError( 1 );
263 fitsigma[i] = m_costheta[i]->GetFunction( "gaus" )->GetParameter( 2 );
264 gain[i] = fitmean[i] / NormalMean;
265 fitchisq[i] = ( m_costheta[i]->GetFunction( "gaus" )->GetChisquare() ) /
266 ( m_costheta[i]->GetFunction( "gaus" )->GetNDF() );
267 }
268 pos_fitentryNo[i] = m_pos_costheta[i]->GetEntries();
269 if ( m_debug )
270 cout << "get results at " << i << " bin with pos_fit entries " << pos_fitentryNo[i]
271 << endl;
272 if ( m_pos_costheta[i]->GetFunction( "gaus" ) )
273 {
274 pos_fitmean[i] = m_pos_costheta[i]->GetFunction( "gaus" )->GetParameter( 1 );
275 pos_fitmeanerr[i] = m_pos_costheta[i]->GetFunction( "gaus" )->GetParError( 1 );
276 pos_fitsigma[i] = m_pos_costheta[i]->GetFunction( "gaus" )->GetParameter( 2 );
277 pos_gain[i] = pos_fitmean[i] / NormalMean;
278 pos_fitchisq[i] = ( m_pos_costheta[i]->GetFunction( "gaus" )->GetChisquare() ) /
279 ( m_pos_costheta[i]->GetFunction( "gaus" )->GetNDF() );
280 }
281 neg_fitentryNo[i] = m_neg_costheta[i]->GetEntries();
282 if ( m_debug )
283 cout << "get results at " << i << " bin with neg_fit entries " << neg_fitentryNo[i]
284 << endl;
285 if ( m_neg_costheta[i]->GetFunction( "gaus" ) )
286 {
287 neg_fitmean[i] = m_neg_costheta[i]->GetFunction( "gaus" )->GetParameter( 1 );
288 neg_fitmeanerr[i] = m_neg_costheta[i]->GetFunction( "gaus" )->GetParError( 1 );
289 neg_fitsigma[i] = m_neg_costheta[i]->GetFunction( "gaus" )->GetParameter( 2 );
290 neg_gain[i] = neg_fitmean[i] / NormalMean;
291 neg_fitchisq[i] = ( m_neg_costheta[i]->GetFunction( "gaus" )->GetChisquare() ) /
292 ( m_neg_costheta[i]->GetFunction( "gaus" )->GetNDF() );
293 }
294 // use unweighted average to reduce the asymmetry in Bhabha
295 if ( m_ave )
296 {
297 fitmean[i] = ( pos_fitmean[i] + neg_fitmean[i] ) / 2;
298 fitmeanerr[i] = sqrt( pow( pos_fitmeanerr[i], 2 ) + pow( neg_fitmeanerr[i], 2 ) );
299 fitsigma[i] = sqrt( pow( pos_fitsigma[i], 2 ) + pow( neg_fitsigma[i], 2 ) );
300 gain[i] = fitmean[i] / NormalMean;
301 fitchisq[i] = ( pos_fitchisq[i] + neg_fitchisq[i] ) / 2;
302 }
303
304 if ( fitentryNo[i] > 100 ) m_dEdxVsCostheta->SetBinContent( i + 1, fitmean[i] );
305 if ( pos_fitentryNo[i] > 100 )
306 m_pos_dEdxVsCostheta->SetBinContent( i + 1, pos_fitmean[i] );
307 if ( neg_fitentryNo[i] > 100 )
308 m_neg_dEdxVsCostheta->SetBinContent( i + 1, neg_fitmean[i] );
309 }
310
311 double dedx1[Size] = { 0 };
312 double costheta1[Size] = { 0 };
313 cout << "Vec_dedx.size() = " << Vec_dedx.size() << endl;
314 for ( unsigned int i = 0; i < Size && i < Vec_dedx.size(); i++ )
315 {
316 dedx1[i] = Vec_dedx[i];
317 costheta1[i] = Vec_costheta[i];
318 // cout<<"i= "<<i<<"dedx= "<<dedx1[i]<<" costheta= "<<costheta1[i]<<endl;
319 }
320
321 log << MSG::INFO << "begin generating root file!!! " << endmsg;
322 TFile* f = new TFile( m_rootfile.c_str(), "RECREATE" );
323 for ( int i = 0; i < NumBinCostheta; i++ )
324 {
325 m_chi[i]->Write();
326 m_pos_chi[i]->Write();
327 m_neg_chi[i]->Write();
328 m_costheta[i]->Write();
329 m_pos_costheta[i]->Write();
330 m_neg_costheta[i]->Write();
331 }
332 m_dEdxVsCostheta->Write();
333 m_pos_dEdxVsCostheta->Write();
334 m_neg_dEdxVsCostheta->Write();
335
336 TTree* costhetacalib = new TTree( "costhetacalib", "costhetacalib" );
337 costhetacalib->Branch( "chientryNo", chientryNo, "chientryNo[80]/D" );
338 costhetacalib->Branch( "chimean", chimean, "chimean[80]/D" );
339 costhetacalib->Branch( "chimeanerr", chimeanerr, "chimeanerr[80]/D" );
340 costhetacalib->Branch( "chisigma", chisigma, "chisigma[80]/D" );
341 costhetacalib->Branch( "chisq", chisq, "chisq[80]/D" );
342 costhetacalib->Branch( "pos_chientryNo", pos_chientryNo, "pos_chientryNo[80]/D" );
343 costhetacalib->Branch( "pos_chimean", pos_chimean, "pos_chimean[80]/D" );
344 costhetacalib->Branch( "pos_chimeanerr", pos_chimeanerr, "pos_chimeanerr[80]/D" );
345 costhetacalib->Branch( "pos_chisigma", pos_chisigma, "pos_chisigma[80]/D" );
346 costhetacalib->Branch( "pos_chisq", pos_chisq, "pos_chisq[80]/D" );
347 costhetacalib->Branch( "neg_chientryNo", neg_chientryNo, "neg_chientryNo[80]/D" );
348 costhetacalib->Branch( "neg_chimean", neg_chimean, "neg_chimean[80]/D" );
349 costhetacalib->Branch( "neg_chimeanerr", neg_chimeanerr, "neg_chimeanerr[80]/D" );
350 costhetacalib->Branch( "neg_chisigma", neg_chisigma, "neg_chisigma[80]/D" );
351 costhetacalib->Branch( "neg_chisq", neg_chisq, "neg_chisq[80]/D" );
352 costhetacalib->Branch( "fitentryNo", fitentryNo, "fitentryNo[80]/D" );
353 costhetacalib->Branch( "fitmean", fitmean, "fitmean[80]/D" );
354 costhetacalib->Branch( "fitmeanerr", fitmeanerr, "fitmeanerr[80]/D" );
355 costhetacalib->Branch( "fitsigma", fitsigma, "fitsigma[80]/D" );
356 costhetacalib->Branch( "costheta", gain, "gain[80]/D" );
357 costhetacalib->Branch( "fitchisq", fitchisq, "fitchisq[80]/D" );
358 costhetacalib->Branch( "pos_fitentryNo", pos_fitentryNo, "pos_fitentryNo[80]/D" );
359 costhetacalib->Branch( "pos_fitmean", pos_fitmean, "pos_fitmean[80]/D" );
360 costhetacalib->Branch( "pos_fitmeanerr", pos_fitmeanerr, "pos_fitmeanerr[80]/D" );
361 costhetacalib->Branch( "pos_fitsigma", pos_fitsigma, "pos_fitsigma[80]/D" );
362 costhetacalib->Branch( "pos_gain", pos_gain, "pos_gain[80]/D" );
363 costhetacalib->Branch( "pos_fitchisq", pos_fitchisq, "pos_fitchisq[80]/D" );
364 costhetacalib->Branch( "neg_fitentryNo", neg_fitentryNo, "neg_fitentryNo[80]/D" );
365 costhetacalib->Branch( "neg_fitmean", neg_fitmean, "neg_fitmean[80]/D" );
366 costhetacalib->Branch( "neg_fitmeanerr", neg_fitmeanerr, "neg_fitmeanerr[80]/D" );
367 costhetacalib->Branch( "neg_fitsigma", neg_fitsigma, "neg_fitsigma[80]/D" );
368 costhetacalib->Branch( "neg_gain", neg_gain, "neg_gain[80]/D" );
369 costhetacalib->Branch( "neg_fitchisq", neg_fitchisq, "neg_fitchisq[80]/D" );
370 costhetacalib->Branch( "cosBin", cosBin, "cosBin[80]/D" );
371 costhetacalib->Branch( "costheta1", costheta1, "costheta1[700000]/D" );
372 costhetacalib->Branch( "dedx1", dedx1, "dedx1[700000]/D" );
373 costhetacalib->Fill();
374 costhetacalib->Write();
375
376 TCanvas c1( "c1", "canvas", 500, 400 );
377 c1.Print( ( m_rootfile + ".ps[" ).c_str() );
378 costhetacalib->Draw( "dedx1:costheta1", "dedx1>200 && dedx1<1000" );
379 c1.Print( ( m_rootfile + ".ps" ).c_str() );
380 m_dEdxVsCostheta->Draw();
381 c1.Print( ( m_rootfile + ".ps" ).c_str() );
382 m_pos_dEdxVsCostheta->Draw();
383 c1.Print( ( m_rootfile + ".ps" ).c_str() );
384 m_neg_dEdxVsCostheta->Draw();
385 c1.Print( ( m_rootfile + ".ps" ).c_str() );
386 for ( Int_t i = 0; i < NumBinCostheta; i++ )
387 {
388 m_chi[i]->Draw();
389 c1.Print( ( m_rootfile + ".ps" ).c_str() );
390 }
391 for ( Int_t i = 0; i < NumBinCostheta; i++ )
392 {
393 m_pos_chi[i]->Draw();
394 c1.Print( ( m_rootfile + ".ps" ).c_str() );
395 }
396 for ( Int_t i = 0; i < NumBinCostheta; i++ )
397 {
398 m_neg_chi[i]->Draw();
399 c1.Print( ( m_rootfile + ".ps" ).c_str() );
400 }
401 for ( Int_t i = 0; i < NumBinCostheta; i++ )
402 {
403 m_costheta[i]->Draw();
404 c1.Print( ( m_rootfile + ".ps" ).c_str() );
405 }
406 for ( Int_t i = 0; i < NumBinCostheta; i++ )
407 {
408 m_pos_costheta[i]->Draw();
409 c1.Print( ( m_rootfile + ".ps" ).c_str() );
410 }
411 for ( Int_t i = 0; i < NumBinCostheta; i++ )
412 {
413 m_neg_costheta[i]->Draw();
414 c1.Print( ( m_rootfile + ".ps" ).c_str() );
415 }
416 c1.Print( ( m_rootfile + ".ps]" ).c_str() );
417 f->Close();
418
419 for ( Int_t i = 0; i < NumBinCostheta; i++ )
420 {
421 delete m_chi[i];
422 delete m_pos_chi[i];
423 delete m_neg_chi[i];
424 delete m_costheta[i];
425 delete m_pos_costheta[i];
426 delete m_neg_costheta[i];
427 }
428 delete m_dEdxVsCostheta;
429 delete m_pos_dEdxVsCostheta;
430 delete m_neg_dEdxVsCostheta;
431}
std::string m_rootfile
Definition DedxCalib.h:57
char * c_str(Index i)

The documentation for this class was generated from the following files: