7ff718443b160eedb5649e96f64a0c7602644b5b
[u/mrichter/AliRoot.git] / T0 / AliT0CalibTimeEq.cxx
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  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // class for T0 calibration                       TM-AC-AM_6-02-2006  
21 // equalize time shift for each time CFD channel
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include "AliT0CalibTimeEq.h"
26 #include "AliLog.h"
27 #include <TFile.h>
28 #include <TMath.h>
29 #include <TF1.h>
30 #include <TSpectrum.h>
31 #include <TProfile.h>
32 #include <iostream>
33
34 ClassImp(AliT0CalibTimeEq)
35
36 //________________________________________________________________
37   AliT0CalibTimeEq::AliT0CalibTimeEq():TNamed(),
38                                        fMeanVertex(0),        
39                                        fRmsVertex(0)      
40 {
41   //
42
43 }
44
45 //________________________________________________________________
46 AliT0CalibTimeEq::AliT0CalibTimeEq(const char* name):TNamed(),
47                                        fMeanVertex(0),        
48                                        fRmsVertex(0)      
49 {
50   //constructor
51
52   TString namst = "Calib_";
53   namst += name;
54   SetName(namst.Data());
55   SetTitle(namst.Data());
56 }
57
58 //________________________________________________________________
59 AliT0CalibTimeEq::AliT0CalibTimeEq(const AliT0CalibTimeEq& calibda):TNamed(calibda),            
60                                        fMeanVertex(0),        
61                                        fRmsVertex(0)      
62 {
63 // copy constructor
64   SetName(calibda.GetName());
65   SetTitle(calibda.GetName());
66
67
68 }
69
70 //________________________________________________________________
71 AliT0CalibTimeEq &AliT0CalibTimeEq::operator =(const AliT0CalibTimeEq& calibda)
72 {
73 // assignment operator
74   SetName(calibda.GetName());
75   SetTitle(calibda.GetName());
76  
77   return *this;
78 }
79
80 //________________________________________________________________
81 AliT0CalibTimeEq::~AliT0CalibTimeEq()
82 {
83   //
84   // destrictor
85 }
86 //________________________________________________________________
87 void AliT0CalibTimeEq::Reset()
88 {
89   //reset values
90
91   memset(fCFDvalue,0,120*sizeof(Float_t));
92   memset(fTimeEq,1,24*sizeof(Float_t));
93 }
94
95
96 //________________________________________________________________
97 void  AliT0CalibTimeEq::Print(Option_t*) const
98 {
99   // print time values
100
101   printf("\n    ----    PM Arrays       ----\n\n");
102   printf(" Time delay CFD \n");
103   for (Int_t i=0; i<24; i++) printf(" CFD  %f ",fTimeEq[i]);
104   printf("\n Mean Vertex %f \n", fMeanVertex);
105
106
107
108 //________________________________________________________________
109 Bool_t AliT0CalibTimeEq::ComputeOnlineParams(const char* filePhys)
110 {
111   // compute online equalized time
112   Float_t meandiff, sigmadiff, meanver, meancfdtime, sigmacfdtime;
113   meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
114   Double_t rmsver=0;
115  Int_t nent=0;
116   Bool_t ok=false;
117   gFile = TFile::Open(filePhys);
118   if(!gFile) {
119     AliError("No input PHYS data found ");
120   }
121   else
122     {
123         gFile->ls();
124       ok=true;
125       for (Int_t i=0; i<24; i++)
126         {
127           meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
128           TH1F *cfd = (TH1F*) gFile->Get(Form("CFD1minCFD%d",i+1));
129           TH1F *cfdtime = (TH1F*) gFile->Get(Form("CFD%d",i+1));
130           if(!cfd) AliWarning(Form("no Diff histograms collected by PHYS DA for channel %i", i));
131           if(cfd) {
132             nent = Int_t(cfd->GetEntries());
133             if(nent>500 )  {
134               if(cfd->GetRMS()>1.5 &&  cfd->GetRMS()<20)
135                 GetMeanAndSigma(cfd, meandiff, sigmadiff);
136               if(cfd->GetRMS()<=1.5) 
137                 {
138                   meandiff = cfd->GetMean();
139                   sigmadiff=cfd->GetRMS();
140                 }
141               
142               if(cfd->GetRMS()>20) 
143                 {
144                   ok=false;
145                   AliWarning(Form("Data is not good  in PMT %i - mean %f rsm %f nentries %i", i,meandiff,sigmadiff , nent));
146                 }
147             }
148             else 
149               {
150                 ok=false;
151                 AliWarning(Form(" Not  enouph data in PMT %i- PMT1:  %i ", i, nent));
152               }
153             if(!cfd) AliWarning(Form("no CFD histograms collected by PHYS DA for channel %i", i));
154           }
155             //      printf(" i = %d buf1 = %s\n", i, buf1);
156           if(cfdtime) {
157             nent = Int_t(cfdtime->GetEntries());
158             if(nent>500 )  {
159               if(cfdtime->GetRMS()>1.5 &&  cfdtime->GetRMS()<30)
160                 GetMeanAndSigma(cfdtime,meancfdtime, sigmacfdtime);
161               if(cfdtime->GetRMS()<=1.5) 
162                 {
163                   meancfdtime = cfdtime->GetMean();
164                   sigmacfdtime = cfdtime->GetRMS();
165                 }
166               if(cfdtime->GetRMS()>30) 
167                 {
168                 ok=false;
169                 AliWarning(Form("Data is not good enouph in PMT %i  - meancfdtime %f rsm %f nentries %i", i,meancfdtime, sigmacfdtime, nent));
170                 }
171             }
172           }
173           else 
174             {
175               ok=false;
176               AliWarning(Form(" Not  enouph data in PMT in CFD peak %i - %i ", i, nent));
177             }
178           //      printf(" %i %f %f %f %f \n",i, meandiff, sigmadiff, meancfdtime, sigmacfdtime);
179           SetTimeEq(i,meandiff);
180           SetTimeEqRms(i,sigmadiff);
181           SetCFDvalue(i,0,meancfdtime);
182           SetCFDvalue(i,0,sigmacfdtime);
183           if (cfd) delete cfd;
184           if (cfdtime) delete cfdtime;
185
186         }
187       TH1F *ver = (TH1F*) gFile->Get("hVertex");
188       if(!ver) AliWarning("no T0 histogram collected by PHYS DA ");
189       if(ver) {
190         meanver = ver->GetMean();
191         rmsver = ver->GetRMS();
192       }
193       SetMeanVertex(meanver);
194       SetRmsVertex(rmsver);
195       
196       gFile->Close();
197       delete gFile;
198
199     }
200     return ok; 
201 }
202
203 //________________________________________________________________
204 Bool_t AliT0CalibTimeEq::ComputeOfflineParams(const char* filePhys)
205 {
206   // compute online equalized time
207   Float_t meandiff, sigmadiff, meanver, meancfdtime, sigmacfdtime;
208   meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
209   // Double_t rms=0, rmscfd=0; 
210   Double_t rmsver=0;
211   Int_t nent=0;
212   Bool_t ok=false;
213   gFile = TFile::Open(filePhys);
214   if(!gFile) {
215     AliError("No input PHYS data found ");
216   }
217   else
218     {
219       meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
220       ok=true;
221       TObjArray * TzeroObj = (TObjArray*) gFile->Get("fTzeroObject");
222       for (Int_t i=0; i<24; i++)
223         {
224           TH1F *cfddiff = (TH1F*)TzeroObj->At(i);
225           TH1F *cfdtime = (TH1F*)TzeroObj->At(i+24);
226           if(!cfddiff) AliWarning(Form("no Diff histograms collected by PHYS DA for channel %i", i));
227           //      printf(" i = %d buf1 = %s\n", i, buf1);
228           if(cfddiff) {
229             nent = Int_t(cfddiff->GetEntries());
230             if(nent>500 )  {
231               if(cfddiff->GetRMS()>1.5 &&  cfddiff->GetRMS()<20)
232                 GetMeanAndSigma(cfddiff, meandiff, sigmadiff);
233               if(cfddiff->GetRMS()<=1.5) 
234                 {
235                   meandiff = cfddiff->GetMean();
236                   sigmadiff = cfddiff->GetRMS();
237                 }
238               
239               if(cfddiff->GetRMS()>20) 
240                 {
241                   ok=false;
242                   AliWarning(Form("Data is not good  in PMT %i - mean %f rsm %f nentries %i", i,meandiff,sigmadiff , nent));
243                 }
244             }
245             else 
246               {
247                 ok=false;
248                 AliWarning(Form(" Not  enouph data in PMT %i- PMT1:  %i ", i, nent));
249               }
250           }         
251             //      printf(" i = %d buf1 = %s\n", i, buf1);
252           if(cfdtime) {
253             nent = Int_t(cfdtime->GetEntries());
254             if(nent>500 )  {
255               if(cfdtime->GetRMS()>1.5 &&  cfdtime->GetRMS()<30)
256                 GetMeanAndSigma(cfdtime,meancfdtime, sigmacfdtime);
257               if(cfdtime->GetRMS()<=1.5) 
258                 {
259                   meancfdtime = cfdtime->GetMean();
260                   sigmacfdtime=cfdtime->GetRMS();
261                 }
262               if(cfdtime->GetRMS()>30) 
263                 {
264                   ok=false;
265                   AliWarning(Form("Data is not good enouph in PMT %i  - meancfdtime %f rsm %f nentries %i", i,meancfdtime, sigmacfdtime, nent));
266                 }
267             }
268             else 
269               {
270                 ok=false;
271                 AliWarning(Form(" Not  enouph data in PMT in CFD peak %i - %i ", i, nent));
272               }
273           }
274
275           SetTimeEq(i,meandiff);
276           SetTimeEqRms(i,sigmadiff);
277           SetCFDvalue(i,0,meancfdtime);
278           SetCFDvalue(i,0,sigmacfdtime);
279           if (cfddiff) delete cfddiff;
280           if (cfdtime) delete cfdtime;
281
282         }
283       TH1F *ver = (TH1F*) gFile->Get("hVertex");
284       if(!ver) AliWarning("no T0 histogram collected by PHYS DA ");
285       if(ver) {
286         meanver = ver->GetMean();
287         rmsver = ver->GetRMS();
288       }
289       SetMeanVertex(meanver);
290       SetRmsVertex(rmsver);
291       
292       gFile->Close();
293       delete gFile;
294
295     }
296     return ok; 
297    }
298
299 //________________________________________________________________________
300 void AliT0CalibTimeEq::GetMeanAndSigma(TH1F* hist,  Float_t &mean, Float_t &sigma) {
301
302   const double window = 5.;  //fit window 
303  
304   double meanEstimate, sigmaEstimate; 
305   int maxBin;
306   maxBin        =  hist->GetMaximumBin(); //position of maximum
307   meanEstimate  =  hist->GetBinCenter( maxBin); // mean of gaussian sitting in maximum
308   sigmaEstimate = hist->GetRMS();
309   TF1* fit= new TF1("fit","gaus", meanEstimate - window*sigmaEstimate, meanEstimate + window*sigmaEstimate);
310   fit->SetParameters(hist->GetBinContent(maxBin), meanEstimate, sigmaEstimate);
311   hist->Fit("fit","R");
312
313   mean  = (Float_t) fit->GetParameter(1);
314   sigma = (Float_t) fit->GetParameter(2);
315
316   delete fit;
317 }
318
319