BOSS 7.0.1
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 ()
 
 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 ()
 
virtual void genNtuple ()=0
 
virtual void initializing ()=0
 
virtual void BookHists ()=0
 
virtual void FillHists ()=0
 
virtual void AnalyseHists ()=0
 
virtual void WriteHists ()=0
 
 DedxCalib (const std::string &name, ISvcLocator *pSvcLocator)
 
 ~DedxCalib ()
 
StatusCode initialize ()
 
StatusCode execute ()
 
StatusCode finalize ()
 
virtual void genNtuple ()=0
 
virtual void initializing ()=0
 
virtual void BookHists ()=0
 
virtual void FillHists ()=0
 
virtual void AnalyseHists ()=0
 
virtual void WriteHists ()=0
 

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 ()
 
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

Constructor & Destructor Documentation

◆ DedxCalibCostheta() [1/2]

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

Definition at line 20 of file DedxCalibCostheta.cxx.

20 : DedxCalib(name, pSvcLocator){
21 declareProperty("Debug", m_debug=false);
22 declareProperty("DebugMin", m_debug_min=2);
23 declareProperty("DebugMax", m_debug_max=2);
24}

◆ ~DedxCalibCostheta() [1/2]

DedxCalibCostheta::~DedxCalibCostheta ( )
inline

◆ DedxCalibCostheta() [2/2]

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

◆ ~DedxCalibCostheta() [2/2]

DedxCalibCostheta::~DedxCalibCostheta ( )
inline

Member Function Documentation

◆ AnalyseHists() [1/2]

void DedxCalibCostheta::AnalyseHists ( )
virtual

Implements DedxCalib.

Definition at line 156 of file DedxCalibCostheta.cxx.

157{
158 MsgStream log(msgSvc(), name());
159 log<<MSG::INFO<<"DedxCalibCostheta::AnalyseHists()"<<endreq;
160
161 gStyle->SetOptFit(1111);
162 for(int i=0; i<NumBinCostheta; i++)
163 {
164 if(m_debug) cout << "num of bin " << i << endl;
165 if(m_costheta[i]->GetEntries()>100) m_costheta[i]->Fit("gaus", "Q" );
166 if(m_pos_costheta[i]->GetEntries()>100) m_pos_costheta[i]->Fit("gaus", "Q" );
167 if(m_neg_costheta[i]->GetEntries()>100) m_neg_costheta[i]->Fit("gaus", "Q" );
168 if(m_chi[i]->GetEntries()>100) m_chi[i]->Fit("gaus", "Q" );
169 if(m_pos_chi[i]->GetEntries()>100) m_pos_chi[i]->Fit("gaus", "Q" );
170 if(m_neg_chi[i]->GetEntries()>100) m_neg_chi[i]->Fit("gaus", "Q" );
171 }
172}

◆ AnalyseHists() [2/2]

void DedxCalibCostheta::AnalyseHists ( )
virtual

Implements DedxCalib.

◆ BookHists() [1/2]

void DedxCalibCostheta::BookHists ( )
virtual

Implements DedxCalib.

Definition at line 26 of file DedxCalibCostheta.cxx.

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

◆ BookHists() [2/2]

void DedxCalibCostheta::BookHists ( )
virtual

Implements DedxCalib.

◆ FillHists() [1/2]

void DedxCalibCostheta::FillHists ( )
virtual

Implements DedxCalib.

Definition at line 81 of file DedxCalibCostheta.cxx.

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

◆ FillHists() [2/2]

void DedxCalibCostheta::FillHists ( )
virtual

Implements DedxCalib.

◆ genNtuple() [1/2]

void DedxCalibCostheta::genNtuple ( )
inlinevirtual

◆ genNtuple() [2/2]

void DedxCalibCostheta::genNtuple ( )
inlinevirtual

◆ initializing() [1/2]

void DedxCalibCostheta::initializing ( )
inlinevirtual

◆ initializing() [2/2]

void DedxCalibCostheta::initializing ( )
inlinevirtual

◆ WriteHists() [1/2]

void DedxCalibCostheta::WriteHists ( )
virtual

Implements DedxCalib.

Definition at line 174 of file DedxCalibCostheta.cxx.

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

◆ WriteHists() [2/2]

void DedxCalibCostheta::WriteHists ( )
virtual

Implements DedxCalib.


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