BOSS 6.6.4.p03
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxCalibLayerGain.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2
3#include <sstream>
4#include <string>
5#include "TTree.h"
6#include "TH1F.h"
7#include "TF1.h"
8#include "TFile.h"
9#include "TStyle.h"
10#include "TCanvas.h"
11
13
14using namespace std;
15
16DedxCalibLayerGain::DedxCalibLayerGain(const std::string& name, ISvcLocator* pSvcLocator): DedxCalib(name, pSvcLocator){}
17
18const int layNo = 43;
20void calculate(vector<double> phlist)
21{
22 double mean=0,rms=0;
23 for(int i=0;i<phlist.size();i++)
24 {
25 mean+=phlist[i];
26 }
27 mean/=phlist.size();
28 for(int i=0;i<phlist.size();i++)
29 {
30 rms += pow((phlist[i]-mean),2);
31 }
32 //cout<<"phlist.size()= "<<phlist.size()<<" rms= "<<rms<<endl;
33 rms = sqrt(rms/phlist.size());
34 //cout<<"mean = "<<mean<<" rms= "<<rms<<endl;
35 m_mean = mean;
36 m_rms = rms;
37}
38double getMean() {return m_mean;}
39double getRms() {return m_rms;}
40
42{
43 MsgStream log(msgSvc(), name());
44 log<<MSG::INFO<<"DedxCalibLayerGain::BookHists()"<<endreq;
45
47
48 m_laygain = new TH1F*[layNo];
49 m_laygain_gaus = new TH1F*[layNo];
50 stringstream hlabel;
51 for(int i=0; i<layNo; i++)
52 {
53 hlabel.str("");
54 hlabel<<"dEdx_Lay_"<<i;
55 m_laygain[i] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue, MaxHistValue);
56 hlabel.str("");
57 hlabel<<"dEdx_gaus_Lay_"<<i;
58 m_laygain_gaus[i] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue1, MaxHistValue1);
59 }
60}
61
63{
64 MsgStream log(msgSvc(), name());
65 log<<MSG::INFO<<"DedxCalibLayerGain::FillHists()"<<endreq;
66
67 TFile* f;
68 TTree* n102;
69 string runlist;
70
71 double dedx=0;
72 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;
73 cout<<"m_recFileVector.size()= "<<m_recFileVector.size()<<endl;
74 for(int i=0; i<m_recFileVector.size(); i++)
75 {
76 runlist = m_recFileVector[i];
77 cout<<"runlist: "<<runlist.c_str()<<endl;
78 f = new TFile(runlist.c_str());
79 n102 = (TTree*)f->Get("n102");
80 n102-> SetBranchAddress("adc_raw",&dedx);
81 n102-> SetBranchAddress("path_rphi",&pathlength);
82 n102-> SetBranchAddress("wire",&wid);
83 n102-> SetBranchAddress("layer",&layid);
84 n102-> SetBranchAddress("doca_in",&dd_in);
85 n102-> SetBranchAddress("driftdist",&driftdist);
86 n102-> SetBranchAddress("eangle",&eangle);
87 n102-> SetBranchAddress("zhit",&zhit);
88 n102-> SetBranchAddress("runNO",&runNO);
89 n102-> SetBranchAddress("runFlag",&runFlag);
90 n102-> SetBranchAddress("costheta1",&costheta);
91 n102-> SetBranchAddress("t01",&tes);
92 n102->SetBranchAddress("ptrk1",&ptrk);
93
94 for(int j=0;j<n102->GetEntries();j++)
95 {
96 n102->GetEntry(j);
97 if(tes>1400) continue;
98 if (ptrk<0.1) continue;
99 //if(wid>=6731 && wid<=6739) continue;
100 if(layid <8)
101 {if(dedx<MinADCValuecut || dedx>MaxADCValuecut || pathlength>RPhi_PathMaxCut || pathlength<Iner_RPhi_PathMinCut || driftdist>Iner_DriftDistCut) continue;}
102 else
103 {if(dedx<MinADCValuecut || dedx>MaxADCValuecut || pathlength>RPhi_PathMaxCut || pathlength<Out_RPhi_PathMinCut || driftdist>Out_DriftDistCut) continue;}
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_laygain[(int)layid]->Fill(dedx);
107
108 m_vector[(int)layid].push_back(dedx);
109 }
110 }
111
112}
113
115{
116 MsgStream log(msgSvc(), name());
117 log<<MSG::INFO<<"DedxCalibLayerGain::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
125 stringstream hlabel;
126 double dedxt;
127 vector<double> phlist;
128 vector<double> phlist_gaus;
129 for(int i=0; i<layNo; i++)
130 {
131 entry = m_laygain[i]->GetEntries();
132 if(entry<10) continue;
133 mean = m_laygain[i]->GetMean();
134 rms = m_laygain[i]->GetRMS();
135 binmax = m_laygain[i]->GetMaximumBin();
136 lp[1] = m_laygain[i]->GetBinCenter(binmax);
137 lp[2] = 200;
138 lp[0] = (m_laygain[i]->GetBinContent(binmax)+m_laygain[i]->GetBinContent(binmax-1)+m_laygain[i]->GetBinContent(binmax+1))/3;
139
140 if(m_phShapeFlag==0)
141 {
142 func = new TF1("mylan",mylan,10,3000,4);
143 func->SetParameters(entry, mean, rms, -0.9);
144 }
145 else if(m_phShapeFlag==1)
146 {
147 func = new TF1("landaun",landaun,10,3000,3);
148 func->SetParameters(lp[0], lp[1], lp[2] );
149 }
150 else if(m_phShapeFlag==2)
151 {
152 func = new TF1("Landau",Landau,10,3000,2);
153 func->SetParameters(0.7*mean, rms );
154 }
155 else if(m_phShapeFlag==3)
156 {
157 func = new TF1("Vavilov",Vavilov,10,3000,2);
158 func->SetParameters(0.05, 0.6);
159 }
160 else if(m_phShapeFlag==4)
161 {
162 func = new TF1("AsymGauss",AsymGauss,10,3000,4);
163 func->SetParameters(entry, mean, rms, 1.0 );
164 }
165 func->SetLineWidth(2.1);
166 func->SetLineColor(2);
167
168 m_laygain[i]->Fit(func, "r" );
169
170
171 //******* begin truncated dedx fitting **************************//
172 for(int j=0;j<m_vector[i].size();j+=100)
173 {
174 for(int k=0;k<100;k++) phlist.push_back(m_vector[i][j+k]);
175 dedxt = cal_dedx_bitrunc(truncate, phlist);
176 phlist_gaus.push_back(dedxt);
177 phlist.clear();
178 }
179
180 hlabel.str("");
181 hlabel<<"dEdx_gaus_Lay_"<<i;
182 calculate(phlist_gaus);
183 //cout<<getMean()-10*getRms()<<" "<<getMean()+10*getRms()<<endl;
184 m_laygain_gaus[i] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, getMean()-10*getRms(), getMean()+10*getRms());
185 for(int j=0;j<phlist_gaus.size();j++) m_laygain_gaus[i]->Fill(phlist_gaus[j]);
186 phlist_gaus.clear();
187 m_laygain_gaus[i]->Fit("gaus","Q");
188
189
190 }
191}
192
194{
195 MsgStream log(msgSvc(), name());
196 log<<MSG::INFO<<"DedxCalibLayerGain::WriteHists()"<<endreq;
197
198
199 double entryNo[layNo]={0},mean[layNo]={0},sigma[layNo]={0},proper[layNo]={0},fitmean[layNo]={0},fitmeanerr[layNo]={0},fitsigma[layNo]={0},layerg[layNo]={0},layer[layNo]={0},chisq[layNo]={0};
200 double fitmean_gaus[layNo]={0},fitmeanerr_gaus[layNo]={0},fitsigma_gaus[layNo]={0},layerg_gaus[layNo]={0};
201 for(int i=0; i<layNo; i++)
202 {
203 fitmean_gaus[i] = m_laygain_gaus[i]->GetFunction("gaus")->GetParameter(1);
204 fitmeanerr_gaus[i] = m_laygain_gaus[i]->GetFunction("gaus")->GetParError(1);
205 fitsigma_gaus[i] = m_laygain_gaus[i]->GetFunction("gaus")->GetParameter(2);
206
207 entryNo[i] = m_laygain[i]->GetEntries();
208 mean[i] = m_laygain[i]->GetMean();
209 sigma[i] = m_laygain[i]->GetRMS();
210 proper[i] = m_laygain[i]->GetBinCenter(m_laygain[i]->GetMaximumBin());
211 if(entryNo[i]<10) continue;
212
213 if(m_phShapeFlag==0)
214 {
215 fitmean[i] = m_laygain[i]->GetFunction("mylan")->GetParameter(1);
216 fitmeanerr[i] = m_laygain[i]->GetFunction("mylan")->GetParError(1);
217 fitsigma[i] = m_laygain[i]->GetFunction("mylan")->GetParameter(2);
218 chisq[i] = (m_laygain[i]->GetFunction("mylan")-> GetChisquare())/(m_laygain[i]->GetFunction("mylan")-> GetNDF());
219 }
220 else if(m_phShapeFlag==1)
221 {
222 fitmean[i] = m_laygain[i]->GetFunction("landaun")->GetParameter(1);
223 fitmeanerr[i] = m_laygain[i]->GetFunction("landaun")->GetParError(1);
224 fitsigma[i] = m_laygain[i]->GetFunction("landaun")->GetParameter(2);
225 chisq[i] = (m_laygain[i]->GetFunction("landaun")-> GetChisquare())/(m_laygain[i]->GetFunction("landaun")-> GetNDF());
226 }
227 else if(m_phShapeFlag==2)
228 {
229 fitmean[i] = m_laygain[i]->GetFunction("Landau")->GetParameter(1);
230 fitmeanerr[i] = m_laygain[i]->GetFunction("Landau")->GetParError(1);
231 fitsigma[i] = m_laygain[i]->GetFunction("Landau")->GetParameter(2);
232 chisq[i] = (m_laygain[i]->GetFunction("Landau")-> GetChisquare())/(m_laygain[i]->GetFunction("Landau")-> GetNDF());
233 }
234 else if(m_phShapeFlag==3)
235 {
236 fitmean[i] = m_laygain[i]->GetFunction("Vavilov")->GetParameter(1);
237 fitmeanerr[i] = m_laygain[i]->GetFunction("Vavilov")->GetParError(1);
238 fitsigma[i] = m_laygain[i]->GetFunction("Vavilov")->GetParameter(2);
239 chisq[i] = (m_laygain[i]->GetFunction("Vavilov")-> GetChisquare())/(m_laygain[i]->GetFunction("Vavilov")-> GetNDF());
240 }
241 else if(m_phShapeFlag==4)
242 {
243 fitmean[i] = m_laygain[i]->GetFunction("AsymGauss")->GetParameter(1);
244 fitmeanerr[i] = m_laygain[i]->GetFunction("AsymGauss")->GetParError(1);
245 fitsigma[i] = m_laygain[i]->GetFunction("AsymGauss")->GetParameter(2);
246 chisq[i] = (m_laygain[i]->GetFunction("AsymGauss")-> GetChisquare())/(m_laygain[i]->GetFunction("AsymGauss")-> GetNDF());
247 }
248 //cout<<"fitmean[i]= "<<fitmean[i]<<" fitmeanerr[i]= "<<fitmeanerr[i]<<" fitsigma[i]= "<<fitsigma[i]<<endl;
249 }
250
251 double sum=0, sum_gaus=0, n=0;
252 for(int i=4; i<layNo; i++)
253 {
254 if(fitmean[i]>0 && fitmeanerr[i]>0 && fitmeanerr[i]<10 && fitsigma[i]>0 && fitsigma[i]<200 && entryNo[i]>1500)
255 {
256 sum += fitmean[i];
257 sum_gaus+= fitmean_gaus[i];
258 n++;
259 }
260 }
261 sum/=n;
262 sum_gaus/=n;
263 cout<<"average value of good layer: "<<sum<<endl;
264
265 for(int i=0;i<layNo;i++)
266 {
267 layerg[i] = fitmean[i]/sum;
268 layerg_gaus[i] = fitmean_gaus[i]/sum_gaus;
269 layer[i] = i;
270 }
271
272
273 log<<MSG::INFO<<"begin generating root file!!! "<<endreq;
274 TFile* f = new TFile(m_rootfile.c_str(), "recreate");
275 for(int i=0;i<layNo;i++)
276 {
277 m_laygain[i]->Write();
278 m_laygain_gaus[i]->Write();
279 }
280
281 TTree* layergcalib = new TTree("layergcalib", "layergcalib");
282 layergcalib->Branch("layerg_gaus", layerg_gaus, "layerg_gaus[43]/D");
283 layergcalib->Branch("layerg", layerg, "layerg[43]/D");
284 layergcalib->Branch("layer", layer, "layer[43]/D");
285 layergcalib->Branch("entryNo", entryNo, "entryNo[43]/D");
286 layergcalib->Branch("mean", mean, "mean[43]/D");
287 layergcalib->Branch("sigma", sigma, "sigma[43]/D");
288 layergcalib->Branch("fitmean", fitmean, "fitmean[43]/D");
289 layergcalib->Branch("fitmeanerr", fitmeanerr, "fitmeanerr[43]/D");
290 layergcalib->Branch("fitsigma", fitsigma, "fitsigma[43]/D");
291 layergcalib->Branch("chisq", chisq, "chisq[43]/D");
292 layergcalib->Branch("fitmean_gaus", fitmean_gaus, "fitmean_gaus[43]/D");
293 layergcalib->Branch("fitmeanerr_gaus", fitmeanerr_gaus, "fitmeanerr_gaus[43]/D");
294 layergcalib->Branch("fitsigma_gaus", fitsigma_gaus, "fitsigma_gaus[43]/D");
295 layergcalib->Fill();
296 layergcalib->Write();
297
298 f->Close();
299
300 TCanvas c1("c1", "canvas", 500, 400);
301 c1.Print((m_rootfile+".ps[").c_str());
302 for(int i=0;i<layNo;i++)
303 {
304 m_laygain[i]->Draw();
305 c1.Print((m_rootfile+".ps").c_str());
306 }
307 c1.Print((m_rootfile+".ps]").c_str());
308
309 for(int i=0;i<layNo;i++)
310 {
311 delete m_laygain[i];
312 delete m_laygain_gaus[i];
313 }
314}
const Int_t n
data SetBranchAddress("time",&time)
double m_rms
const int layNo
double getMean()
double getRms()
void calculate(vector< double > phlist)
double m_mean
double Landau(double *x, double *par)
#define MaxADCValuecut
#define MinHistValue1
double Vavilov(double *x, double *par)
#define MaxHistValue1
#define Iner_DriftDistCut
#define RPhi_PathMaxCut
#define MinHistValue
double landaun(double *x, double *par)
#define Out_DriftDistCut
double AsymGauss(double *x, double *par)
#define NumHistBins
#define MaxHistValue
double mylan(double *x, double *par)
IMessageSvc * msgSvc()
DedxCalibLayerGain(const std::string &name, ISvcLocator *pSvcLocator)
void ReadRecFileList()
Definition: DedxCalib.cxx:89
int m_phShapeFlag
Definition: DedxCalib.h:52
std::string m_rootfile
Definition: DedxCalib.h:55
IDedxCorrecSvc * exsvc
Definition: DedxCalib.h:30
double cal_dedx_bitrunc(float truncate, std::vector< double > phlist)
Definition: DedxCalib.cxx:108
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, 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
TCanvas * c1
Definition: tau_mode.c:75