BOSS 7.1.2
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.

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

◆ ~DedxCalibCostheta()

DedxCalibCostheta::~DedxCalibCostheta ( )
inline

Definition at line 16 of file DedxCalibCostheta.h.

16{};

Member Function Documentation

◆ AnalyseHists()

void DedxCalibCostheta::AnalyseHists ( )
virtual

Implements DedxCalib.

Definition at line 164 of file DedxCalibCostheta.cxx.

165{
166 MsgStream log(msgSvc(), name());
167 log<<MSG::INFO<<"DedxCalibCostheta::AnalyseHists()"<<endreq;
168
169 gStyle->SetOptFit(1111);
170 for(int i=0; i<NumBinCostheta; i++)
171 {
172 if(m_debug) cout << "num of bin " << i << endl;
173 if(m_costheta[i]->GetEntries()>100) m_costheta[i]->Fit("gaus", "Q" );
174 if(m_pos_costheta[i]->GetEntries()>100) m_pos_costheta[i]->Fit("gaus", "Q" );
175 if(m_neg_costheta[i]->GetEntries()>100) m_neg_costheta[i]->Fit("gaus", "Q" );
176 if(m_chi[i]->GetEntries()>100) m_chi[i]->Fit("gaus", "Q" );
177 if(m_pos_chi[i]->GetEntries()>100) m_pos_chi[i]->Fit("gaus", "Q" );
178 if(m_neg_chi[i]->GetEntries()>100) m_neg_chi[i]->Fit("gaus", "Q" );
179 }
180}
#define NumBinCostheta
IMessageSvc * msgSvc()

◆ BookHists()

void DedxCalibCostheta::BookHists ( )
virtual

Implements DedxCalib.

Definition at line 30 of file DedxCalibCostheta.cxx.

31{
32 MsgStream log(msgSvc(), name());
33 log<<MSG::INFO<<"DedxCalibCostheta::BookHists()"<<endreq;
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, MinHistValue1, MaxHistValue1);
49 hlabel.str("");
50 hlabel<<"pos_dEdx_costheta"<<i;
51 m_pos_costheta[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue1, MaxHistValue1);
52 hlabel.str("");
53 hlabel<<"neg_dEdx_costheta"<<i;
54 m_neg_costheta[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue1, MaxHistValue1);
55 hlabel.str("");
56 hlabel<<"chi_costheta"<<i;
57 m_chi[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinChiValue, MaxChiValue);
58 hlabel.str("");
59 hlabel<<"pos_chi_costheta"<<i;
60 m_pos_chi[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinChiValue, MaxChiValue);
61 hlabel.str("");
62 hlabel<<"neg_chi_costheta"<<i;
63 m_neg_chi[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinChiValue, MaxChiValue);
64 }
65 hlabel.str("");
66 hlabel<<"dEdxVsCostheta";
67 m_dEdxVsCostheta = new TH1F(hlabel.str().c_str(),"dEdx vs costheta",NumBinCostheta, CosthetaMin, CosthetaMax);
68 m_dEdxVsCostheta->GetXaxis()->SetTitle("cos#theta");
69 m_dEdxVsCostheta->GetYaxis()->SetTitle("dE/dx");
70 hlabel.str("");
71 hlabel<<"pos_dEdxVsCostheta";
72 m_pos_dEdxVsCostheta = new TH1F(hlabel.str().c_str(),"positron dEdx vs costheta",NumBinCostheta, CosthetaMin, CosthetaMax);
73 m_pos_dEdxVsCostheta->GetXaxis()->SetTitle("cos#theta");
74 m_pos_dEdxVsCostheta->GetYaxis()->SetTitle("dE/dx");
75 hlabel.str("");
76 hlabel<<"neg_dEdxVsCostheta";
77 m_neg_dEdxVsCostheta = new TH1F(hlabel.str().c_str(),"electron dEdx vs costheta",NumBinCostheta, CosthetaMin, CosthetaMax);
78 m_neg_dEdxVsCostheta->GetXaxis()->SetTitle("cos#theta");
79 m_neg_dEdxVsCostheta->GetYaxis()->SetTitle("dE/dx");
80
81 Vec_dedx.clear();
82 Vec_costheta.clear();
83}
#define MinHistValue1
#define MaxChiValue
#define MaxHistValue1
#define CosthetaMin
#define MinChiValue
#define CosthetaMax
#define NumHistBins
void ReadRecFileList()
Definition DedxCalib.cxx:89

◆ FillHists()

void DedxCalibCostheta::FillHists ( )
virtual

Implements DedxCalib.

Definition at line 85 of file DedxCalibCostheta.cxx.

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

◆ genNtuple()

void DedxCalibCostheta::genNtuple ( )
inlinevirtual

Implements DedxCalib.

Definition at line 19 of file DedxCalibCostheta.h.

19{}

◆ initializing()

void DedxCalibCostheta::initializing ( )
inlinevirtual

Implements DedxCalib.

Definition at line 17 of file DedxCalibCostheta.h.

17{}

◆ WriteHists()

void DedxCalibCostheta::WriteHists ( )
virtual

Implements DedxCalib.

Definition at line 182 of file DedxCalibCostheta.cxx.

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

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