CGEM BOSS 6.6.5.g
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxCalibEAng Class Reference

#include <DedxCalibEAng.h>

+ Inheritance diagram for DedxCalibEAng:

Public Member Functions

 DedxCalibEAng (const std::string &name, ISvcLocator *pSvcLocator)
 
 ~DedxCalibEAng ()
 
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
 

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 10 of file DedxCalibEAng.h.

Constructor & Destructor Documentation

◆ DedxCalibEAng()

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

Definition at line 18 of file DedxCalibEAng.cxx.

18: DedxCalib(name, pSvcLocator){}

◆ ~DedxCalibEAng()

DedxCalibEAng::~DedxCalibEAng ( )
inline

Definition at line 14 of file DedxCalibEAng.h.

14{};

Member Function Documentation

◆ AnalyseHists()

void DedxCalibEAng::AnalyseHists ( )
virtual

Implements DedxCalib.

Definition at line 114 of file DedxCalibEAng.cxx.

115{
116 MsgStream log(msgSvc(), name());
117 log<<MSG::INFO<<"DedxCalibEAng::AnalyseHists()"<<endreq;
118
119 TF1* func;
120 Double_t entry=0,mean=0,rms=0;
121 Double_t binmax=0;
122 Double_t lp[3]={0};
123 gStyle->SetOptFit(1111);
124 for(Int_t ip=0;ip<NumSlices;ip++)
125 {
126 entry = m_eangle[ip]->Integral();
127 if(entry<10) continue;
128 mean = m_eangle[ip]->GetMean();
129 rms = m_eangle[ip]->GetRMS();
130 binmax = m_eangle[ip]->GetMaximumBin();
131 lp[1] = m_eangle[ip]->GetBinCenter(binmax);
132 lp[2] = 200;
133 lp[0] = (m_eangle[ip]->GetBinContent(binmax)+m_eangle[ip]->GetBinContent(binmax-1)+m_eangle[ip]->GetBinContent(binmax+1))/3;
134
135 if(m_phShapeFlag==0)
136 {
137 func = new TF1("mylan",mylan,10,3000,4);
138 func->SetParLimits(0, 0, entry);
139 func->SetParLimits(1, 0, mean);
140 func->SetParLimits(2, 10, rms);
141 func->SetParLimits(3, -1, 0);
142 func->SetParameters(entry/36., mean-200, rms/2.3, -0.5);
143 }
144 else if(m_phShapeFlag==1)
145 {
146 func = new TF1("landaun",landaun,10,3000,3);
147 func->SetParameters(lp[0], lp[1], lp[2] );
148 }
149 else if(m_phShapeFlag==2)
150 {
151 func = new TF1("Landau",Landau,10,3000,2);
152 func->SetParameters(0.7*mean, rms );
153 }
154 else if(m_phShapeFlag==3)
155 {
156 func = new TF1("Vavilov",Vavilov,10,3000,2);
157 func->SetParameters(0.05, 0.6);
158 }
159 else if(m_phShapeFlag==4)
160 {
161 func = new TF1("AsymGauss",AsymGauss,10,3000,4);
162 func->SetParameters(entry, mean, rms, 1.0 );
163 }
164 func->SetLineWidth(2.1);
165 func->SetLineColor(2);
166
167 m_eangle[ip]->Fit(func, "rq" );// , "", mean-rms, mean+rms/2);
168 m_pos_eangle[ip]->Fit(func, "rq");// , "", mean-rms, mean+rms/2);
169 m_neg_eangle[ip]->Fit(func, "rq");// , "", mean-rms, mean+rms/2);
170 }
171}
double Landau(double *x, double *par)
double Vavilov(double *x, double *par)
#define NumSlices
double landaun(double *x, double *par)
double AsymGauss(double *x, double *par)
double mylan(double *x, double *par)
IMessageSvc * msgSvc()
int m_phShapeFlag
Definition: DedxCalib.h:52

◆ BookHists()

void DedxCalibEAng::BookHists ( )
virtual

Implements DedxCalib.

Definition at line 20 of file DedxCalibEAng.cxx.

21{
22 MsgStream log(msgSvc(), name());
23 log<<MSG::INFO<<"DedxCalibEAng::BookHists()"<<endreq;
24
26
27 m_eangle = new TH1F*[NumSlices];
28 m_pos_eangle = new TH1F*[NumSlices];
29 m_neg_eangle = new TH1F*[NumSlices];
30 stringstream hlabel;
31 for(int i=0;i<NumSlices;i++)
32 {
33 hlabel.str("");
34 hlabel<<"dEdx"<<"_eang"<<i;
35 m_eangle[i] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue, MaxHistValue);
36 hlabel.str("");
37 hlabel<<"pos_dEdx"<<"_eang"<<i;
38 m_pos_eangle[i] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue, MaxHistValue);
39 hlabel.str("");
40 hlabel<<"neg_dEdx"<<"_eang"<<i;
41 m_neg_eangle[i] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue, MaxHistValue);
42 }
43 hlabel.str("");
44 hlabel<<"dEdxVsEAng";
45 m_EAngAverdE = new TH1F(hlabel.str().c_str(),"dEdx & EAng",NumSlices,PhiMin,PhiMax);
46 hlabel.str("");
47 hlabel<<"pos_dEdxVsEAng";
48 m_pos_EAngAverdE = new TH1F(hlabel.str().c_str(),"pos_dEdx & EAng",NumSlices,PhiMin,PhiMax);
49 hlabel.str("");
50 hlabel<<"neg_dEdxVsEAng";
51 m_neg_EAngAverdE = new TH1F(hlabel.str().c_str(),"neg_dEdx & EAng",NumSlices,PhiMin,PhiMax);
52
53}
#define MinHistValue
#define PhiMax
#define PhiMin
#define NumHistBins
#define MaxHistValue
void ReadRecFileList()
Definition: DedxCalib.cxx:89

◆ FillHists()

void DedxCalibEAng::FillHists ( )
virtual

Implements DedxCalib.

Definition at line 55 of file DedxCalibEAng.cxx.

56{
57 MsgStream log(msgSvc(), name());
58 log<<MSG::INFO<<"DedxCalibEAng::FillHists()"<<endreq;
59
60 TFile* f;
61 TTree* n102;
62 string runlist;
63
64 double dedx=0;
65 float runNO=0,runFlag=0,pathlength=0,wid=0,layid=0,dd_in=0,driftdist=0,eangle=0,zhit=0,costheta=0,tes=0,ptrk=0,charge=0;
66 cout<<"m_recFileVector.size()= "<<m_recFileVector.size()<<endl;
67 for(int i=0; i<m_recFileVector.size(); i++)
68 {
69 runlist = m_recFileVector[i];
70 cout<<"runlist: "<<runlist.c_str()<<endl;
71 f = new TFile(runlist.c_str());
72 n102 = (TTree*)f->Get("n102");
73 n102-> SetBranchAddress("adc_raw",&dedx);
74 n102-> SetBranchAddress("path_rphi",&pathlength);
75 n102-> SetBranchAddress("wire",&wid);
76 n102-> SetBranchAddress("layer",&layid);
77 n102-> SetBranchAddress("doca_in",&dd_in);
78 n102-> SetBranchAddress("driftdist",&driftdist);
79 n102-> SetBranchAddress("eangle",&eangle);
80 n102-> SetBranchAddress("zhit",&zhit);
81 n102-> SetBranchAddress("runNO",&runNO);
82 n102-> SetBranchAddress("runFlag",&runFlag);
83 n102-> SetBranchAddress("costheta1",&costheta);
84 n102-> SetBranchAddress("t01",&tes);
85 n102->SetBranchAddress("ptrk1",&ptrk);
86 n102->SetBranchAddress("charge",&charge);
87
88 for(int j=0;j<n102->GetEntries();j++)
89 {
90 n102->GetEntry(j);
91 if(eangle >0.25) eangle = eangle -0.5;
92 else if(eangle <-0.25) eangle = eangle +0.5;
93 if(eangle>PhiMin && eangle<PhiMax)
94 {
95 if(tes>1400) continue;
96 if (ptrk<0.1) continue;
97 if(layid <4) continue;
98 else if(layid <8)
99 {if(dedx<MinADCValuecut || dedx>MaxADCValuecut || pathlength>RPhi_PathMaxCut || pathlength<Iner_RPhi_PathMinCut || driftdist>Iner_DriftDistCut) continue;}
100 else
101 {if(dedx<MinADCValuecut || dedx>MaxADCValuecut || pathlength>RPhi_PathMaxCut || pathlength<Out_RPhi_PathMinCut || driftdist>Out_DriftDistCut) continue;}
102 int iEAng = 0;
103 iEAng = (int) ((eangle-PhiMin)/StepEAng);
104 dedx = exsvc->StandardHitCorrec(0,runFlag,2,runNO,pathlength,wid,layid,dedx,dd_in,eangle,costheta);
105 dedx = exsvc->StandardTrackCorrec(0, runFlag, 2, runNO, dedx, costheta, tes);
106 m_eangle[iEAng]->Fill(dedx);
107 if(charge>0) m_pos_eangle[iEAng]->Fill(dedx);
108 if(charge<0) m_neg_eangle[iEAng]->Fill(dedx);
109 }
110 }
111 }
112}
data SetBranchAddress("time",&time)
const double StepEAng
#define MaxADCValuecut
#define Iner_DriftDistCut
#define RPhi_PathMaxCut
#define Out_DriftDistCut
IDedxCorrecSvc * exsvc
Definition: DedxCalib.h:30
vector< string > m_recFileVector
Definition: DedxCalib.h:48
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

◆ genNtuple()

void DedxCalibEAng::genNtuple ( )
inlinevirtual

Implements DedxCalib.

Definition at line 17 of file DedxCalibEAng.h.

17{}

◆ initializing()

void DedxCalibEAng::initializing ( )
inlinevirtual

Implements DedxCalib.

Definition at line 15 of file DedxCalibEAng.h.

15{}

◆ WriteHists()

void DedxCalibEAng::WriteHists ( )
virtual

Implements DedxCalib.

Definition at line 173 of file DedxCalibEAng.cxx.

174{
175 MsgStream log(msgSvc(), name());
176 log<<MSG::INFO<<"DedxCalibEAng::WriteHists()"<<endreq;
177
178 double entryNo[NumSlices]={0},mean[NumSlices]={0},sigma[NumSlices]={0},fitmean[NumSlices]={0},fitmeanerr[NumSlices]={0},fitsigma[NumSlices]={0},gain[NumSlices]={0},chisq[NumSlices]={0};
179 double pos_entryNo[NumSlices]={0},pos_mean[NumSlices]={0},pos_sigma[NumSlices]={0},pos_fitmean[NumSlices]={0},pos_fitmeanerr[NumSlices]={0},pos_fitsigma[NumSlices]={0},pos_gain[NumSlices]={0},pos_chisq[NumSlices]={0};
180 double neg_entryNo[NumSlices]={0},neg_mean[NumSlices]={0},neg_sigma[NumSlices]={0},neg_fitmean[NumSlices]={0},neg_fitmeanerr[NumSlices]={0},neg_fitsigma[NumSlices]={0},neg_gain[NumSlices]={0},neg_chisq[NumSlices]={0};
181 double Ip_eangle[NumSlices]={0},eangle[NumSlices]={0};
182
183 double Norm=0;
184 int count=0;
185 for(Int_t ip=0;ip<NumSlices;ip++)
186 {
187 Ip_eangle[ip] = ip;
188 eangle[ip] = (ip+0.5)*StepEAng + PhiMin;
189
190 entryNo[ip]= m_eangle[ip]->Integral();
191 mean[ip] = m_eangle[ip]->GetMean();
192 sigma[ip] = m_eangle[ip]->GetRMS();
193 pos_entryNo[ip]= m_pos_eangle[ip]->Integral();
194 pos_mean[ip] = m_pos_eangle[ip]->GetMean();
195 pos_sigma[ip] = m_pos_eangle[ip]->GetRMS();
196 neg_entryNo[ip]= m_neg_eangle[ip]->Integral();
197 neg_mean[ip] = m_neg_eangle[ip]->GetMean();
198 neg_sigma[ip] = m_neg_eangle[ip]->GetRMS();
199 if(entryNo[ip]<10 || pos_entryNo[ip]<10 || neg_entryNo[ip]<10) continue;
200
201 if(m_phShapeFlag==0)
202 {
203 fitmean[ip] = m_eangle[ip]->GetFunction("mylan")->GetParameter(1);
204 fitmeanerr[ip] = m_eangle[ip]->GetFunction("mylan")->GetParError(1);
205 fitsigma[ip] = m_eangle[ip]->GetFunction("mylan")->GetParameter(2);
206 chisq[ip] = (m_eangle[ip]->GetFunction("mylan")-> GetChisquare())/(m_eangle[ip]->GetFunction("mylan")-> GetNDF());
207
208 pos_fitmean[ip] = m_pos_eangle[ip]->GetFunction("mylan")->GetParameter(1);
209 pos_fitmeanerr[ip] = m_pos_eangle[ip]->GetFunction("mylan")->GetParError(1);
210 pos_fitsigma[ip] = m_pos_eangle[ip]->GetFunction("mylan")->GetParameter(2);
211 pos_chisq[ip] = (m_pos_eangle[ip]->GetFunction("mylan")-> GetChisquare())/(m_pos_eangle[ip]->GetFunction("mylan")-> GetNDF());
212
213 neg_fitmean[ip] = m_neg_eangle[ip]->GetFunction("mylan")->GetParameter(1);
214 neg_fitmeanerr[ip] = m_neg_eangle[ip]->GetFunction("mylan")->GetParError(1);
215 neg_fitsigma[ip] = m_neg_eangle[ip]->GetFunction("mylan")->GetParameter(2);
216 neg_chisq[ip] = (m_neg_eangle[ip]->GetFunction("mylan")-> GetChisquare())/(m_neg_eangle[ip]->GetFunction("mylan")-> GetNDF());
217 }
218 if(m_phShapeFlag==1)
219 {
220 fitmean[ip] = m_eangle[ip]->GetFunction("landaun")->GetParameter(1);
221 fitmeanerr[ip] = m_eangle[ip]->GetFunction("landaun")->GetParError(1);
222 fitsigma[ip] = m_eangle[ip]->GetFunction("landaun")->GetParameter(2);
223 chisq[ip] = (m_eangle[ip]->GetFunction("landaun")-> GetChisquare())/(m_eangle[ip]->GetFunction("mylan")-> GetNDF());
224
225 pos_fitmean[ip] = m_pos_eangle[ip]->GetFunction("landaun")->GetParameter(1);
226 pos_fitmeanerr[ip] = m_pos_eangle[ip]->GetFunction("landaun")->GetParError(1);
227 pos_fitsigma[ip] = m_pos_eangle[ip]->GetFunction("landaun")->GetParameter(2);
228 pos_chisq[ip] = (m_pos_eangle[ip]->GetFunction("landaun")-> GetChisquare())/(m_pos_eangle[ip]->GetFunction("landaun")-> GetNDF());
229
230 neg_fitmean[ip] = m_neg_eangle[ip]->GetFunction("landaun")->GetParameter(1);
231 neg_fitmeanerr[ip] = m_neg_eangle[ip]->GetFunction("landaun")->GetParError(1);
232 neg_fitsigma[ip] = m_neg_eangle[ip]->GetFunction("landaun")->GetParameter(2);
233 neg_chisq[ip] = (m_neg_eangle[ip]->GetFunction("landaun")-> GetChisquare())/(m_neg_eangle[ip]->GetFunction("landaun")-> GetNDF());
234 }
235 if(m_phShapeFlag==2)
236 {
237 fitmean[ip] = m_eangle[ip]->GetFunction("Landau")->GetParameter(1);
238 fitmeanerr[ip] = m_eangle[ip]->GetFunction("Landau")->GetParError(1);
239 fitsigma[ip] = m_eangle[ip]->GetFunction("Landau")->GetParameter(2);
240 chisq[ip] = (m_eangle[ip]->GetFunction("Landau")-> GetChisquare())/(m_eangle[ip]->GetFunction("Landau")-> GetNDF());
241
242 pos_fitmean[ip] = m_pos_eangle[ip]->GetFunction("Landau")->GetParameter(1);
243 pos_fitmeanerr[ip] = m_pos_eangle[ip]->GetFunction("Landau")->GetParError(1);
244 pos_fitsigma[ip] = m_pos_eangle[ip]->GetFunction("Landau")->GetParameter(2);
245 pos_chisq[ip] = (m_pos_eangle[ip]->GetFunction("Landau")-> GetChisquare())/(m_pos_eangle[ip]->GetFunction("Landau")-> GetNDF());
246
247 neg_fitmean[ip] = m_neg_eangle[ip]->GetFunction("Landau")->GetParameter(1);
248 neg_fitmeanerr[ip] = m_neg_eangle[ip]->GetFunction("Landau")->GetParError(1);
249 neg_fitsigma[ip] = m_neg_eangle[ip]->GetFunction("Landau")->GetParameter(2);
250 neg_chisq[ip] = (m_neg_eangle[ip]->GetFunction("Landau")-> GetChisquare())/(m_neg_eangle[ip]->GetFunction("Landau")-> GetNDF());
251 }
252 if(m_phShapeFlag==3)
253 {
254 fitmean[ip] = m_eangle[ip]->GetFunction("Vavilov")->GetParameter(1);
255 fitmeanerr[ip] = m_eangle[ip]->GetFunction("Vavilov")->GetParError(1);
256 fitsigma[ip] = m_eangle[ip]->GetFunction("Vavilov")->GetParameter(2);
257 chisq[ip] = (m_eangle[ip]->GetFunction("Vavilov")-> GetChisquare())/(m_eangle[ip]->GetFunction("Vavilov")-> GetNDF());
258
259 pos_fitmean[ip] = m_pos_eangle[ip]->GetFunction("Vavilov")->GetParameter(1);
260 pos_fitmeanerr[ip] = m_pos_eangle[ip]->GetFunction("Vavilov")->GetParError(1);
261 pos_fitsigma[ip] = m_pos_eangle[ip]->GetFunction("Vavilov")->GetParameter(2);
262 pos_chisq[ip] = (m_pos_eangle[ip]->GetFunction("Vavilov")-> GetChisquare())/(m_pos_eangle[ip]->GetFunction("Vavilov")-> GetNDF());
263
264 neg_fitmean[ip] = m_neg_eangle[ip]->GetFunction("Vavilov")->GetParameter(1);
265 neg_fitmeanerr[ip] = m_neg_eangle[ip]->GetFunction("Vavilov")->GetParError(1);
266 neg_fitsigma[ip] = m_neg_eangle[ip]->GetFunction("Vavilov")->GetParameter(2);
267 neg_chisq[ip] = (m_neg_eangle[ip]->GetFunction("Vavilov")-> GetChisquare())/(m_neg_eangle[ip]->GetFunction("Vavilov")-> GetNDF());
268 }
269 if(m_phShapeFlag==4)
270 {
271 fitmean[ip] = m_eangle[ip]->GetFunction("AsymGauss")->GetParameter(1);
272 fitmeanerr[ip] = m_eangle[ip]->GetFunction("AsymGauss")->GetParError(1);
273 fitsigma[ip] = m_eangle[ip]->GetFunction("AsymGauss")->GetParameter(2);
274 chisq[ip] = (m_eangle[ip]->GetFunction("AsymGauss")-> GetChisquare())/(m_eangle[ip]->GetFunction("AsymGauss")-> GetNDF());
275
276 pos_fitmean[ip] = m_pos_eangle[ip]->GetFunction("AsymGauss")->GetParameter(1);
277 pos_fitmeanerr[ip] = m_pos_eangle[ip]->GetFunction("AsymGauss")->GetParError(1);
278 pos_fitsigma[ip] = m_pos_eangle[ip]->GetFunction("AsymGauss")->GetParameter(2);
279 pos_chisq[ip] = (m_pos_eangle[ip]->GetFunction("AsymGauss")-> GetChisquare())/(m_pos_eangle[ip]->GetFunction("AsymGauss")-> GetNDF());
280
281 neg_fitmean[ip] = m_neg_eangle[ip]->GetFunction("AsymGauss")->GetParameter(1);
282 neg_fitmeanerr[ip] = m_neg_eangle[ip]->GetFunction("AsymGauss")->GetParError(1);
283 neg_fitsigma[ip] = m_neg_eangle[ip]->GetFunction("AsymGauss")->GetParameter(2);
284 neg_chisq[ip] = (m_neg_eangle[ip]->GetFunction("AsymGauss")-> GetChisquare())/(m_neg_eangle[ip]->GetFunction("AsymGauss")-> GetNDF());
285 }
286 if(entryNo[ip]>1500 && fitmean[ip]>0 && fitmeanerr[ip]>0 && fitmeanerr[ip]<10 && fitsigma[ip]>0 && fitsigma[ip]<200)
287 {
288 Norm+=fitmean[ip];
289 count++;
290 }
291 }
292 Norm = Norm/count;
293 // Norm = 550; // use an absolute normalization
294 cout<<"count= "<<count<<" average value of eangle bins: "<<Norm<<endl;
295 for(Int_t i=0;i<NumSlices;i++)
296 {
297 if(entryNo[i]>10 && fitmean[i]>0 && fitsigma[i]>0) gain[i]=fitmean[i]/Norm;
298 else{
299 gain[i] = 1;
300 fitmeanerr[i]=10000;
301 }
302 if(pos_entryNo[i]>10 && pos_fitmean[i]>0 && pos_fitsigma[i]>0) pos_gain[i]=pos_fitmean[i]/Norm;
303 else pos_gain[i] = 0;
304 if(neg_entryNo[i]>10 && neg_fitmean[i]>0 && neg_fitsigma[i]>0) neg_gain[i]=neg_fitmean[i]/Norm;
305 else neg_gain[i] = 0;
306 if(gain[i] != 0){
307 m_EAngAverdE -> SetBinContent(i+1, gain[i]);
308 m_EAngAverdE->SetBinError(i+1, fitmeanerr[i]/Norm);
309 }
310 if(pos_gain[i] != 0){
311 m_pos_EAngAverdE -> SetBinContent(i+1, pos_gain[i]);
312 m_pos_EAngAverdE->SetBinError(i+1, pos_fitmeanerr[i]/Norm);
313 }
314 if(neg_gain[i] != 0){
315 m_neg_EAngAverdE -> SetBinContent(i+1, neg_gain[i]);
316 m_neg_EAngAverdE->SetBinError(i+1, neg_fitmeanerr[i]/Norm);
317 }
318 }
319
320 log<<MSG::INFO<<"begin getting calibration constant!!! "<<endreq;
321 double denangle[NumSlices]={0};
322 int denangle_entry[1] = {100};
323 m_EAngAverdE->Fit("pol1","r","",-0.25,0.25);
324 double par0 = m_EAngAverdE->GetFunction("pol1")->GetParameter(0);
325 double par1 = m_EAngAverdE->GetFunction("pol1")->GetParameter(1);
326 // double par2 = m_EAngAverdE->GetFunction("pol2")->GetParameter(2);
327 // par1 = -par1; // a trick is used here, because we don't really understand
328 // for(int i=0;i<NumSlices;i++) denangle[i] = (par0+par1*eangle[i])/par0;//+par2*pow(eangle[i],2))/par0;
329 for(int i=0;i<NumSlices;i++) denangle[i] = gain[i]; // use raw values instead of fit
330
331 log<<MSG::INFO<<"begin generating root file!!! "<<endreq;
332 TFile* f = new TFile(m_rootfile.c_str(), "recreate");
333 for(Int_t ip=0;ip<NumSlices;ip++)
334 {
335 m_eangle[ip]->Write();
336 m_pos_eangle[ip]->Write();
337 m_neg_eangle[ip]->Write();
338 }
339 m_EAngAverdE->Write();
340 m_pos_EAngAverdE->Write();
341 m_neg_EAngAverdE->Write();
342
343 TTree* entracalib = new TTree("entracalib","entracalib");
344 entracalib->Branch("1denangle", denangle, "denangle[100]/D");
345 entracalib->Branch("1denangle_entry", denangle_entry, "denangle_entry[1]/I");
346 entracalib->Fill();
347 entracalib-> Write();
348
349 TTree* ddgcalib = new TTree("ddgcalib", "ddgcalib");
350 ddgcalib -> Branch("hits", entryNo, "entryNo[100]/D");
351 ddgcalib -> Branch("mean", mean, "mean[100]/D");
352 ddgcalib -> Branch("sigma", sigma, "sigma[100]/D");
353 ddgcalib -> Branch("fitmean", fitmean, "fitmean[100]/D");
354 ddgcalib -> Branch("fitmeanerr", fitmeanerr, "fitmeanerr[100]/D");
355 ddgcalib -> Branch("chi", fitsigma, "fitsigma[100]/D");
356 ddgcalib -> Branch("gain", gain, "gain[100]/D");
357 ddgcalib -> Branch("chisq", chisq, "chisq[100]/D");
358 ddgcalib -> Branch("pos_hits", pos_entryNo, "pos_entryNo[100]/D");
359 ddgcalib -> Branch("pos_mean", pos_mean, "pos_mean[100]/D");
360 ddgcalib -> Branch("pos_sigma", pos_sigma, "pos_sigma[100]/D");
361 ddgcalib -> Branch("pos_fitmean", pos_fitmean, "pos_fitmean[100]/D");
362 ddgcalib -> Branch("pos_fitmeanerr", pos_fitmeanerr, "pos_fitmeanerr[100]/D");
363 ddgcalib -> Branch("pos_chi", pos_fitsigma, "pos_fitsigma[100]/D");
364 ddgcalib -> Branch("pos_gain", pos_gain, "pos_gain[100]/D");
365 ddgcalib -> Branch("pos_chisq", pos_chisq, "pos_chisq[100]/D");
366 ddgcalib -> Branch("neg_hits", neg_entryNo, "neg_entryNo[100]/D");
367 ddgcalib -> Branch("neg_mean", neg_mean, "neg_mean[100]/D");
368 ddgcalib -> Branch("neg_sigma", neg_sigma, "neg_sigma[100]/D");
369 ddgcalib -> Branch("neg_fitmean", neg_fitmean, "neg_fitmean[100]/D");
370 ddgcalib -> Branch("neg_fitmeanerr", neg_fitmeanerr, "neg_fitmeanerr[100]/D");
371 ddgcalib -> Branch("neg_chi", neg_fitsigma, "neg_fitsigma[100]/D");
372 ddgcalib -> Branch("neg_gain", neg_gain, "neg_gain[100]/D");
373 ddgcalib -> Branch("neg_chisq", neg_chisq, "neg_chisq[100]/D");
374 ddgcalib -> Branch("Ip_eangle", Ip_eangle, "Ip_eangle[100]/D");
375 ddgcalib -> Branch("eangle", eangle, "eangle[100]/D");
376 ddgcalib -> Fill();
377 ddgcalib -> Write();
378 f->Close();
379
380 TCanvas c1("c1", "canvas", 500, 400);
381 c1.Print((m_rootfile+".ps[").c_str());
382 m_EAngAverdE->Draw();
383 c1.Print((m_rootfile+".ps").c_str());
384 m_pos_EAngAverdE->Draw();
385 c1.Print((m_rootfile+".ps").c_str());
386 m_neg_EAngAverdE->Draw();
387 c1.Print((m_rootfile+".ps").c_str());
388 for(Int_t ip=0;ip<NumSlices;ip++)
389 {
390 m_eangle[ip]->Draw();
391 c1.Print((m_rootfile+".ps").c_str());
392 }
393 for(Int_t ip=0;ip<NumSlices;ip++)
394 {
395 m_pos_eangle[ip]->Draw();
396 c1.Print((m_rootfile+".ps").c_str());
397 }
398 for(Int_t ip=0;ip<NumSlices;ip++)
399 {
400 m_neg_eangle[ip]->Draw();
401 c1.Print((m_rootfile+".ps").c_str());
402 }
403 c1.Print((m_rootfile+".ps]").c_str());
404
405 for(Int_t ip=0;ip<NumSlices;ip++)
406 {
407 delete m_eangle[ip];
408 delete m_pos_eangle[ip];
409 delete m_neg_eangle[ip];
410 }
411 delete m_EAngAverdE;
412 delete m_pos_EAngAverdE;
413 delete m_neg_EAngAverdE;
414
415}
std::string m_rootfile
Definition: DedxCalib.h:55
char * c_str(Index i)
Definition: EvtCyclic3.cc:252
uint32_t count(const node_t &list)
Definition: node.cxx:42
TCanvas * c1
Definition: tau_mode.c:75

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