BOSS 7.0.5
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxCalibWireGain.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
12#include "DedxCalibAlg/DedxCalibWireGain.h"
13
14using namespace std;
15
16DedxCalibWireGain::DedxCalibWireGain(const std::string& name, ISvcLocator* pSvcLocator): DedxCalib(name, pSvcLocator){
17 declareProperty("PMin", pMin=0.2);
18 declareProperty("PMax", pMax=2.3);
19}
20
21const int wireNo = 6796;
22
24{
25 MsgStream log(msgSvc(), name());
26 log<<MSG::INFO<<"DedxCalibWireGain::BookHists()"<<endreq;
27
29
30 m_wiregain = new TH1F*[wireNo];
31 stringstream hlabel;
32 for(int i=0; i<wireNo; i++)
33 {
34 hlabel.str("");
35 hlabel<<"dEdx_Wire_"<<i;
36 m_wiregain[i] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue, MaxHistValue);
37 }
38}
39
41{
42 MsgStream log(msgSvc(), name());
43 log<<MSG::INFO<<"DedxCalibWireGain::FillHists()"<<endreq;
44
45 string runlist;
46
47 double dedx=0;
48 float runNO=0,evtNO=0,runFlag=0,pathlength=0,wid=0,layid=0,dd_in=0,driftdist=0,eangle=0,zhit=0,costheta=0,tes=0,ptrk=0;
49 cout<<"m_recFileVector.size()= "<<m_recFileVector.size()<<endl;
50 for(int i=0; i<(int)m_recFileVector.size(); i++)
51 {
52 runlist = m_recFileVector[i];
53 cout<<"runlist: "<<runlist.c_str()<<endl;
54 TFile* f;
55 TTree* n102;
56 f = new TFile(runlist.c_str());
57 n102 = (TTree*)f->Get("n102");
58 n102-> SetBranchAddress("adc_raw",&dedx);
59 n102-> SetBranchAddress("path_rphi",&pathlength);
60 n102-> SetBranchAddress("wire",&wid);
61 n102-> SetBranchAddress("layer",&layid);
62 n102-> SetBranchAddress("doca_in",&dd_in);
63 n102-> SetBranchAddress("driftdist",&driftdist);
64 n102-> SetBranchAddress("eangle",&eangle);
65 n102-> SetBranchAddress("zhit",&zhit);
66 n102-> SetBranchAddress("runNO",&runNO);
67 n102-> SetBranchAddress("evtNO",&evtNO);
68 n102-> SetBranchAddress("runFlag",&runFlag);
69 n102-> SetBranchAddress("costheta1",&costheta);
70 n102-> SetBranchAddress("t01",&tes);
71 n102-> SetBranchAddress("ptrk1",&ptrk);
72
73 cout <<"entries in this file " << n102->GetEntries() << endl;
74
75 for(int j=0;j<n102->GetEntries();j++)
76 {
77 n102->GetEntry(j);
78 if(tes>1400) continue;
79 if(ptrk>pMax || ptrk<pMin) continue;
80 if(layid <8)
81 {if(dedx<MinADCValuecut || dedx>MaxADCValuecut || pathlength>RPhi_PathMaxCut || pathlength<Iner_RPhi_PathMinCut || driftdist>Iner_DriftDistCut) continue;}
82 else
83 {if(dedx<MinADCValuecut || dedx>MaxADCValuecut || pathlength>RPhi_PathMaxCut || pathlength<Out_RPhi_PathMinCut || driftdist>Out_DriftDistCut) continue;}
84 dedx = exsvc->StandardHitCorrec(0,runFlag,2,runNO,evtNO,pathlength,wid,layid,dedx,dd_in,eangle,costheta);
85 dedx = exsvc->StandardTrackCorrec(0, runFlag, 2, runNO, evtNO, dedx, costheta, tes);
86 m_wiregain[(int)wid]->Fill(dedx);
87 }
88 f->Delete();
89 }
90}
91
93{
94 MsgStream log(msgSvc(), name());
95 log<<MSG::INFO<<"DedxCalibWireGain::AnalyseHists()"<<endreq;
96
97 TF1* func = new TF1();
98 Double_t entry=0,mean=0,rms=0;
99 Double_t binmax=0;
100 Double_t lp[3]={0};
101 gStyle->SetOptFit(1111);
102 for(int wire=0; wire<wireNo; wire++)
103 {
104 entry = m_wiregain[wire]->Integral();
105 if(entry<10) continue;
106 mean = m_wiregain[wire]->GetMean();
107 rms = m_wiregain[wire]->GetRMS();
108 binmax = m_wiregain[wire]->GetMaximumBin();
109 lp[1] = m_wiregain[wire]->GetBinCenter(binmax);
110 lp[2] = 200;
111 lp[0] = (m_wiregain[wire]->GetBinContent(binmax)+m_wiregain[wire]->GetBinContent(binmax-1)+m_wiregain[wire]->GetBinContent(binmax+1))/3;
112
113 if(m_phShapeFlag==0)
114 {
115 func = new TF1("mylan",mylan,10,3000,4);
116 func->SetParameters(entry, mean, rms, -0.8);
117 func->SetParLimits(0, 0, entry);
118 func->SetParLimits(1, 0, mean);
119 func->SetParLimits(2, 10, rms);
120 func->SetParLimits(3, -3, 3);
121 }
122 else if(m_phShapeFlag==1)
123 {
124 func = new TF1("landaun",landaun,10,3000,3);
125 func->SetParameters(lp[0], lp[1], lp[2] );
126 }
127 else if(m_phShapeFlag==2)
128 {
129 func = new TF1("Landau",Landau,10,3000,2);
130 func->SetParameters(0.7*mean, rms );
131 }
132 else if(m_phShapeFlag==3)
133 {
134 func = new TF1("Vavilov",Vavilov,10,3000,2);
135 func->SetParameters(0.05, 0.6);
136 }
137 else if(m_phShapeFlag==4)
138 {
139 func = new TF1("AsymGauss",AsymGauss,10,3000,4);
140 func->SetParameters(entry, mean, rms, 1.0 );
141 }
142 func->SetLineWidth(2.1);
143 func->SetLineColor(2);
144
145 m_wiregain[wire]->Fit(func, "rq", "", mean-rms, mean+rms/2 );
146 }
147}
148
150{
151 MsgStream log(msgSvc(), name());
152 log<<MSG::INFO<<"DedxCalibWireGain::WriteHists()"<<endreq;
153
154
155 double entryNo[wireNo]={0},mean[wireNo]={0},sigma[wireNo]={0},proper[wireNo]={0},fitmean[wireNo]={0},fitmeanerr[wireNo]={0},fitsigma[wireNo]={0},wireg[wireNo]={0},wire[wireNo]={0},chisq[wireNo]={0};
156 for(int i=0; i<wireNo; i++)
157 {
158 entryNo[i] = m_wiregain[i]->Integral();
159 mean[i] = m_wiregain[i]->GetMean();
160 sigma[i] = m_wiregain[i]->GetRMS();
161 proper[i] = m_wiregain[i]->GetBinCenter(m_wiregain[i]->GetMaximumBin());
162 if(entryNo[i]<10) continue;
163
164 if(m_phShapeFlag==0)
165 {
166 fitmean[i] = m_wiregain[i]->GetFunction("mylan")->GetParameter(1);
167 fitmeanerr[i] = m_wiregain[i]->GetFunction("mylan")->GetParError(1);
168 fitsigma[i] = m_wiregain[i]->GetFunction("mylan")->GetParameter(2);
169 chisq[i] = (m_wiregain[i]->GetFunction("mylan")-> GetChisquare())/(m_wiregain[i]->GetFunction("mylan")-> GetNDF());
170 }
171 else if(m_phShapeFlag==1)
172 {
173 fitmean[i] = m_wiregain[i]->GetFunction("landaun")->GetParameter(1);
174 fitmeanerr[i] = m_wiregain[i]->GetFunction("landaun")->GetParError(1);
175 fitsigma[i] = m_wiregain[i]->GetFunction("landaun")->GetParameter(2);
176 chisq[i] = (m_wiregain[i]->GetFunction("landaun")-> GetChisquare())/(m_wiregain[i]->GetFunction("landaun")-> GetNDF());
177 }
178 else if(m_phShapeFlag==2)
179 {
180 fitmean[i] = m_wiregain[i]->GetFunction("Landau")->GetParameter(1);
181 fitmeanerr[i] = m_wiregain[i]->GetFunction("Landau")->GetParError(1);
182 fitsigma[i] = m_wiregain[i]->GetFunction("Landau")->GetParameter(2);
183 chisq[i] = (m_wiregain[i]->GetFunction("Landau")-> GetChisquare())/(m_wiregain[i]->GetFunction("Landau")-> GetNDF());
184 }
185 else if(m_phShapeFlag==3)
186 {
187 fitmean[i] = m_wiregain[i]->GetFunction("Vavilov")->GetParameter(1);
188 fitmeanerr[i] = m_wiregain[i]->GetFunction("Vavilov")->GetParError(1);
189 fitsigma[i] = m_wiregain[i]->GetFunction("Vavilov")->GetParameter(2);
190 chisq[i] = (m_wiregain[i]->GetFunction("Vavilov")-> GetChisquare())/(m_wiregain[i]->GetFunction("Vavilov")-> GetNDF());
191 }
192 else if(m_phShapeFlag==4)
193 {
194 fitmean[i] = m_wiregain[i]->GetFunction("AsymGauss")->GetParameter(1);
195 fitmeanerr[i] = m_wiregain[i]->GetFunction("AsymGauss")->GetParError(1);
196 fitsigma[i] = m_wiregain[i]->GetFunction("AsymGauss")->GetParameter(2);
197 chisq[i] = (m_wiregain[i]->GetFunction("AsymGauss")-> GetChisquare())/(m_wiregain[i]->GetFunction("AsymGauss")-> GetNDF());
198 }
199 //cout<<"fitmean[i]= "<<fitmean[i]<<" fitmeanerr[i]= "<<fitmeanerr[i]<<" fitsigma[i]= "<<fitsigma[i]<<endl;
200 }
201
202 double Norm=0, count=0;
203 for(int i=188; i<wireNo; i++)
204 {
205 if(fitmean[i]>0 && fitmeanerr[i]>0 && fitmeanerr[i]<10 && fitsigma[i]>0 && fitsigma[i]<200 && entryNo[i]>1500)
206 {
207 Norm += fitmean[i];
208 count++;
209 }
210 }
211 Norm/=count;
212 cout<<"count= "<<count<<" average value of good wire: "<<Norm<<endl;
213 // double Norm(550); // use 550 instead of average because we will handle wireg later outside of the program
214
215 for(int i=0;i<wireNo;i++)
216 {
217 wireg[i] = fitmean[i]/Norm;
218 wire[i] = i;
219 }
220
221
222 log<<MSG::INFO<<"begin generating root file!!! "<<endreq;
223 TFile* f = new TFile(m_rootfile.c_str(), "recreate");
224 for(int i=0;i<wireNo;i++) m_wiregain[i]->Write();
225
226 TTree* wiregcalib = new TTree("wiregcalib", "wiregcalib");
227 wiregcalib->Branch("wireg", wireg, "wireg[6796]/D");
228 wiregcalib->Branch("wire", wire, "wire[6796]/D");
229 wiregcalib->Branch("entryNo", entryNo, "entryNo[6796]/D");
230 wiregcalib->Branch("proper", proper, "proper[6796]/D");
231 wiregcalib->Branch("mean", mean, "mean[6796]/D");
232 wiregcalib->Branch("sigma", sigma, "sigma[6796]/D");
233 wiregcalib->Branch("fitmean", fitmean, "fitmean[6796]/D");
234 wiregcalib->Branch("fitmeanerr", fitmeanerr, "fitmeanerr[6796]/D");
235 wiregcalib->Branch("fitsigma", fitsigma, "fitsigma[6796]/D");
236 wiregcalib->Branch("chisq", chisq, "chisq[6796]/D");
237 wiregcalib->Fill();
238 wiregcalib->Write();
239
240 f->Close();
241
242 TCanvas c1("c1", "canvas", 500, 400);
243 c1.Print((m_rootfile+".ps[").c_str());
244 for(int i=0;i<wireNo;i++)
245 {
246 m_wiregain[i]->Draw();
247 c1.Print((m_rootfile+".ps").c_str());
248 }
249 c1.Print((m_rootfile+".ps]").c_str());
250
251 for(int i=0;i<wireNo;i++) delete m_wiregain[i];
252
253}
TTree * sigma
curve Fill()
curve Write()
data SetBranchAddress("time",&time)
const int wireNo
DOUBLE_PRECISION count[3]
double Landau(double *x, double *par)
double Vavilov(double *x, double *par)
double landaun(double *x, double *par)
double AsymGauss(double *x, double *par)
double mylan(double *x, double *par)
DedxCalibWireGain(const std::string &name, ISvcLocator *pSvcLocator)
void ReadRecFileList()
Definition: DedxCalib.cxx:89
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
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
float costheta
float ptrk