T0calib lib added
[u/mrichter/AliRoot.git] / T0 / AliT0CalibWalk.cxx
CommitLineData
a2ad8166 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
271b5bd3 15
a2ad8166 16/* $Id$ */
17
18///////////////////////////////////////////////////////////////////////////////
19// //
20// class T0 walk correction Alla Maevskaya alla@inr.ru 20.11.2007 //
21// //
22///////////////////////////////////////////////////////////////////////////////
23
24#include "AliT0CalibWalk.h"
a2ad8166 25#include "AliLog.h"
a2ad8166 26
27#include <TObjArray.h>
28#include <TGraph.h>
29#include <TFile.h>
a2ad8166 30#include <TH2F.h>
31#include <TMath.h>
32#include <TSystem.h>
33#include <Riostream.h>
34#include <TSpectrum.h>
a2ad8166 35#include <TProfile.h>
8f1790c4 36#include <TF1.h>
a2ad8166 37
cdf48647 38using std::cout;
39using std::endl;
a2ad8166 40ClassImp(AliT0CalibWalk)
41
42//________________________________________________________________
43 AliT0CalibWalk::AliT0CalibWalk(): TNamed(),
44 fWalk(0),
c883fdf2 45 fAmpLEDRec(0),
46 fQTC(0),
d0f0bf2a 47 fAmpLED(0),
48 fCalibByData(kFALSE)
a2ad8166 49{
50 //
51}
52
53//________________________________________________________________
54AliT0CalibWalk::AliT0CalibWalk(const char* name):TNamed(),
55 fWalk(0),
c883fdf2 56 fAmpLEDRec(0), fQTC(0),
d0f0bf2a 57 fAmpLED(0),
58 fCalibByData(kFALSE)
c883fdf2 59
a2ad8166 60{
61 TString namst = "Calib_";
62 namst += name;
63 SetName(namst.Data());
64 SetTitle(namst.Data());
65
66}
67
68//________________________________________________________________
69AliT0CalibWalk::AliT0CalibWalk(const AliT0CalibWalk& calibda) :
70 TNamed(calibda),
71 fWalk(0),
c883fdf2 72 fAmpLEDRec(0),
73 fQTC(0),
d0f0bf2a 74 fAmpLED(0) ,
75 fCalibByData(kFALSE)
a2ad8166 76{
77// copy constructor
78 SetName(calibda.GetName());
79 SetTitle(calibda.GetName());
80
81
82}
83
84//________________________________________________________________
85AliT0CalibWalk &AliT0CalibWalk::operator =(const AliT0CalibWalk& calibda)
86{
87// assignment operator
88 SetName(calibda.GetName());
89 SetTitle(calibda.GetName());
90
91 return *this;
92}
93
94//________________________________________________________________
95AliT0CalibWalk::~AliT0CalibWalk()
96{
97 //
98}
d95247c3 99//________________________________________________________________
a2ad8166 100
d95247c3 101Bool_t AliT0CalibWalk::MakeWalkCorrGraph(const char *laserFile)
a2ad8166 102{
7576945f 103 //make walk corerction for preprocessor
f1799eed 104 Float_t sigma,cfdmean, qtmean, ledmean;
d95247c3 105 Bool_t ok=true;
f1799eed 106 Float_t mips[50];
d0f0bf2a 107 cout<<" @@@@ fCalibByData "<<fCalibByData<<endl;
108
c2337900 109 gFile = TFile::Open(laserFile);
6318dd46 110 if(!gFile) {
111 AliError("No input laser data found ");
112 }
113 else
114 {
8f1790c4 115 // gFile->ls();
d95247c3 116 TH1F* hAmp = (TH1F*) gFile->Get("hAmpLaser");
117 Int_t nmips=0;
f1799eed 118 for (Int_t ibin=0; ibin<2000; ibin++) {
d95247c3 119 Float_t bincont = hAmp->GetBinContent(ibin);
120 if(bincont>0){
121 mips[nmips] = hAmp->GetXaxis()->GetBinCenter(ibin);
122 cout<<ibin<<" bincont "<<bincont<<" amp "<< mips[nmips]<<endl;
123 nmips++;
124 }
c69fcbc1 125 }
d0f0bf2a 126 /*
f1799eed 127 if (nmips<17) {
07342b9d 128 ok=false;
129 return ok;
f1799eed 130 }
d0f0bf2a 131 */
c69fcbc1 132 Float_t x1[50], y1[50];
133 Float_t x2[50], xx2[50],y2[50];
134 Float_t xx1[50],yy1[50], xx[50];
95f63b36 135
f1799eed 136 Float_t cfd0 = 0;
95f63b36 137
f1799eed 138 for (Int_t ii=0; ii<nmips; ii++)
139 x1[ii] = y1[ii] = x2[ii] = y2[ii] = 0;
95f63b36 140
141 for (Int_t i=0; i<24; i++)
c883fdf2 142 {
f1799eed 143 cfd0 = 0;
144 for (Int_t im=0; im<nmips; im++)
8f1790c4 145 {
d95247c3 146 TString cfd = Form("hCFD%i_%i",i+1,im+1);
147 TString qtc = Form("hQTC%i_%i",i+1,im+1);
3740a735 148 TString led = Form("hLED%i_%i",i+1,im+1);
95f63b36 149
d95247c3 150 TH1F *hCFD = (TH1F*) gFile->Get(cfd.Data()) ;
151 TH1F *hLED = (TH1F*) gFile->Get(led.Data());
152 TH1F *hQTC = (TH1F*) gFile->Get(qtc.Data()) ;
f1799eed 153 // hCFD->SetDirectory(0);
154 // hQTC->SetDirectory(0);
155 // hLED->SetDirectory(0);
d95247c3 156 if(!hCFD )
157 AliWarning(Form(" no CFD data in LASER DA for channel %i for amplitude %f MIPs",i,mips[im]));
158 if(!hQTC )
159 AliWarning(Form(" no QTC correction data in LASER DA for channel %i for amplitude %f MIPs",i,mips[im]));
160 if(!hLED)
161 AliWarning(Form(" no LED correction data in LASER DA for channel %i for amplitude %f MIPs",i,mips[im]));
f1799eed 162 if( hCFD && hCFD->GetEntries()<500 ) {
163 ok=false;
164 printf("no peak in CFD spectrum for PMT %i amplitude %i\n",i,im);
165 return ok;
166 }
167 if(hCFD && hCFD->GetEntries()>500 ) {
168 if( hCFD->GetRMS() >= 1.5)
169 GetMeanAndSigma(hCFD, cfdmean, sigma);
8704b209 170 else
f1799eed 171 cfdmean = hCFD->GetMean();
8704b209 172
f1799eed 173 Int_t maxBin = hCFD->GetMaximumBin();
174 Double_t meanEstimate = hCFD->GetBinCenter( maxBin);
175 if(TMath::Abs(meanEstimate - cfdmean) > 20 ) cfdmean = meanEstimate;
176 if (im == 0) cfd0 = cfdmean;
177 y1[im] = cfdmean - cfd0;
178 }
179 if(hQTC && hQTC->GetEntries()>500) {
180 GetMeanAndSigma(hQTC, qtmean, sigma);
181
182 x1[im] = qtmean;
8704b209 183 if( x1[im] == 0) {
07342b9d 184 ok=false;
8704b209 185 printf("no peak in QTC signal for PMT %i amplitude %i\n",i,im);
07342b9d 186 return ok;
95f63b36 187 }
c5a24053 188 }
f1799eed 189
190 if( hLED && hLED->GetEntries()>500) {
191 GetMeanAndSigma(hLED, ledmean, sigma);
192 }
193 else
194 {
d0f0bf2a 195 // ok=false;
f1799eed 196 printf("no peak in LED spectrum for PMT %i amplitude %i\n",i,im);
d0f0bf2a 197 ledmean=0;
198 // return ok;
f1799eed 199 }
200 x2[im] = ledmean;
201 xx2[im] = x2[nmips-im-1];
202
203 xx[im]=mips[im];
204 if (hQTC) delete hQTC;
205 if (hCFD) delete hCFD;
206 if (hLED) delete hLED;
207
e47b556f 208 }
d95247c3 209
8f1790c4 210 for (Int_t imi=0; imi<nmips; imi++)
95f63b36 211 {
8f1790c4 212 yy1[imi] = Float_t (mips[nmips-imi-1]);
213 xx1[imi] = x2[nmips-imi-1];
d95247c3 214 // cout<<" LED rev "<<i<<" "<<imi<<" "<<xx1[imi]<<" "<< yy1[imi]<<"nmips-imi-1 = "<< nmips-imi-1<<endl;
95f63b36 215 }
95f63b36 216
d95247c3 217 if(i==0) cout<<"Making graphs..."<<endl;
c5a24053 218
f1799eed 219 TGraph *grwalkqtc = new TGraph (nmips,x1,y1);
220 grwalkqtc->SetTitle(Form("PMT%i",i));
221 TGraph *grwalkled = new TGraph (nmips,x2,y1);
222 grwalkled->SetTitle(Form("PMT%i",i));
d0f0bf2a 223 if(!fCalibByData)
224 fWalk.AddAtAndExpand(grwalkqtc,i);
f1799eed 225 fAmpLEDRec.AddAtAndExpand(grwalkled,i);
226 // cout<<" add walk "<<i<<endl;
c5a24053 227
f1799eed 228 //fit amplitude graphs to make comparison wth new one
229 TGraph *grampled = new TGraph (nmips,xx1,yy1);
230 TGraph *grqtc = new TGraph (nmips,x1,xx);
231 fQTC.AddAtAndExpand(grqtc,i);
232 fAmpLED.AddAtAndExpand(grampled,i);
233 // cout<<" add amp "<<i<<endl;
234
235 if(i==23)
236 cout<<"Graphs created..."<<endl;
237 }
d0f0bf2a 238 }
239 Float_t xpoint, ypoint, xdata[250], ydata[250];
240 Int_t ipmt;
241 if(fCalibByData) {
242 cout<<" read ingraph "<<endl;
243 ifstream ingraph ("calibfit.txt");
244 for (Int_t i=0; i<24; i++)
245 {
246 for (Int_t ip=0; ip<200; ip++)
247 {
248 ingraph>>ipmt>>xpoint>>ypoint;
249 // cout<<i<<" "<<ipmt<<" "<<ip<<" "<< xpoint<<" "<<ypoint<<endl;
250 xdata[ip]=xpoint;
251 ydata[ip]=ypoint;
252 }
253 for (Int_t ip=200; ip<250; ip++) {
254 xdata[ip] =xdata[ip-1]+10;
255 ydata[ip]=ydata[199];
256 }
257
258 TGraph *grwalkqtc = new TGraph (250,xdata,ydata);
259 grwalkqtc->Print();
260 grwalkqtc->SetTitle(Form("PMT%i",i));
261 fWalk.AddAtAndExpand(grwalkqtc,i);
262 for (Int_t ip=0; ip<250; ip++)
263 {
264 xdata[ip]=0;
265 ydata[ip]=0;
266 }
267 }
268 }
07342b9d 269
d0f0bf2a 270
271return ok;
a2ad8166 272}
f1799eed 273
ff49b4c0 274void AliT0CalibWalk::GetMeanAndSigma(TH1F* hist, Float_t &mean, Float_t &sigma)
275{
276
657abde7 277 const double window = 2.; //fit window
ff49b4c0 278
279 double meanEstimate, sigmaEstimate;
280 int maxBin;
281 maxBin = hist->GetMaximumBin(); //position of maximum
282 meanEstimate = hist->GetBinCenter( maxBin); // mean of gaussian sitting in maximum
283 sigmaEstimate = hist->GetRMS();
284 TF1* fit= new TF1("fit","gaus", meanEstimate - window*sigmaEstimate, meanEstimate + window*sigmaEstimate);
285 fit->SetParameters(hist->GetBinContent(maxBin), meanEstimate, sigmaEstimate);
286 hist->Fit("fit","RQ","Q");
287
288 mean = (Float_t) fit->GetParameter(1);
289 sigma = (Float_t) fit->GetParameter(2);
d0f0bf2a 290 // printf(" mean %f max %f sigma %f \n",mean, meanEstimate , sigma);
ff49b4c0 291
292 delete fit;
293}
c883fdf2 294
295
296