]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0CalibTimeEq.cxx
new reconstruction with pass0 calibration
[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       ok=true;
124       for (Int_t i=0; i<24; i++)
125         {
126           meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
127           TH1F *cfd = (TH1F*) gFile->Get(Form("CFD1minCFD%d",i+1));
128           TH1F *cfdtime = (TH1F*) gFile->Get(Form("CFD%d",i+1));
129           if(!cfd) AliWarning(Form("no Diff histograms collected by PHYS DA for channel %i", i));
130           if(!cfdtime) AliWarning(Form("no CFD 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 )
135                 GetMeanAndSigma(cfd, meandiff, sigmadiff);
136               if(cfd->GetRMS()<=1.5) 
137                 {
138                   meandiff = cfd->GetMean();
139                   sigmadiff=cfd->GetRMS();
140                 }
141               if(TMath::Abs(cfd->GetMean() - meandiff) >10 ) meandiff = cfd->GetMean(); 
142             }
143             else 
144               {
145                 ok=false;
146                 AliWarning(Form(" Not  enouph data in PMT %i- PMT1:  %i ", i, nent));
147               }
148             if(!cfd) AliWarning(Form("no CFD histograms collected by PHYS DA for channel %i", i));
149           }
150           //      printf(" i = %d buf1 = %s\n", i, buf1);
151           if(cfdtime) {
152             nent = Int_t(cfdtime->GetEntries());
153             if(nent > 500 )  { //!!!!!!!!!!
154               if(cfdtime->GetRMS()>1.5 )
155                 GetMeanAndSigma(cfdtime,meancfdtime, sigmacfdtime);
156               if(cfdtime->GetRMS()<=1.5) 
157                 {
158                   meancfdtime = cfdtime->GetMean();
159                   sigmacfdtime = cfdtime->GetRMS();
160                 }
161               if(TMath::Abs(cfdtime->GetMean() - meancfdtime) >20 ) meancfdtime = cfdtime->GetMean(); 
162             }
163           }
164           else 
165             {
166               ok=false;
167               AliWarning(Form(" Not  enouph data in PMT in CFD peak %i - %i ", i, nent));
168             }
169           printf(" %i %f %f %f %f \n",i, meandiff, sigmadiff, meancfdtime, sigmacfdtime);
170           SetTimeEq(i,meandiff);
171           SetTimeEqRms(i,sigmadiff);
172           SetCFDvalue(i,0,meancfdtime);
173           if (cfd) delete cfd;
174           if (cfdtime) delete cfdtime;
175
176         }
177       TH1F *ver = (TH1F*) gFile->Get("hVertex");
178       if(!ver) AliWarning("no T0 histogram collected by PHYS DA ");
179       if(ver) {
180         meanver = ver->GetMean();
181         rmsver = ver->GetRMS();
182       }
183       SetMeanVertex(meanver);
184       SetRmsVertex(rmsver);
185       
186       gFile->Close();
187       delete gFile;
188
189     }
190     return ok; 
191 }
192
193 //________________________________________________________________
194 Bool_t AliT0CalibTimeEq::ComputeOfflineParams(const char* filePhys, Float_t *timecdb, Float_t *cfdvalue, Int_t badpmt)
195 {
196   // compute online equalized time
197   Float_t meandiff, sigmadiff, meanver, meancfdtime, sigmacfdtime;
198   meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
199   Double_t rmsver=0;
200   Int_t nent=0;
201   Bool_t ok=false;
202   gFile = TFile::Open(filePhys);
203   if(!gFile) {
204     AliError("No input PHYS data found ");
205     return ok;
206   }
207   else
208     {
209
210       meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
211       ok=true;
212       TDirectory *dr = (TDirectory*) gFile->Get("T0Calib");
213       if (!dr ) { AliError("No input T0calib found "); 
214         return ok;
215       }
216       TObjArray * TzeroObj = (TObjArray*) dr->Get("fTzeroObject");
217       for (Int_t i=0; i<24; i++)
218         {
219           printf("@@@ were in OCDB before %f \n", timecdb[i]);
220           if (i != badpmt) {
221           TH1F *cfddiff = (TH1F*)TzeroObj->At(i);
222           TH1F *cfdtime = (TH1F*)TzeroObj->At(i+24);
223           if(!cfddiff) AliWarning(Form("no Diff histograms collected by PHYS DA for channel %i", i));
224           //      printf(" i = %d buf1 = %s\n", i, buf1);
225           if(!cfdtime) AliWarning(Form("no CFD histograms collected by PHYS DA for channel %i", i));
226           if(cfddiff) {
227             nent = Int_t(cfddiff->GetEntries());
228             if(nent>500 )  { //!!!!!
229               if(cfddiff->GetRMS()>1.5 )
230                 GetMeanAndSigma(cfddiff, meandiff, sigmadiff);
231               if(cfddiff->GetRMS()<=1.5) 
232                 {
233                   meandiff = cfddiff->GetMean();
234                   sigmadiff = cfddiff->GetRMS();
235                 }
236               if(TMath::Abs(cfddiff->GetMean() - meandiff) >10 ) meandiff = cfddiff->GetMean(); 
237               
238             }
239             else 
240               {
241                 ok=false;
242                 AliWarning(Form(" Not  enouph data in PMT %i- PMT1:  %i ", i, nent));
243               }
244           }         
245           if(cfdtime) {
246             nent = Int_t(cfdtime->GetEntries());
247             if(nent>500 )  { //!!!!!
248               if(cfdtime->GetRMS()>1.5 )
249                 GetMeanAndSigma(cfdtime,meancfdtime, sigmacfdtime);
250               if(cfdtime->GetRMS()<=1.5) 
251                 {
252                   meancfdtime = cfdtime->GetMean();
253                   sigmacfdtime=cfdtime->GetRMS();
254                 }
255               if(TMath::Abs(cfdtime->GetMean() - meancfdtime) >20 ) meancfdtime = cfdtime->GetMean(); 
256             }
257             else 
258               {
259                 ok=false;
260                 AliWarning(Form(" Not  enouph data in PMT in CFD peak %i - %i ", i, nent));
261               }
262           }
263           SetTimeEq(i,timecdb[i] + meandiff);
264           SetTimeEqRms(i,sigmadiff);
265           SetCFDvalue(i,0,cfdvalue[i] + meancfdtime );
266           //SetCFDvalue(i,0,cfdvalue[i] );
267           if (cfddiff) delete cfddiff;
268           if (cfdtime) delete cfdtime;
269           } //bad pmt
270         }
271       TH1F *ver = (TH1F*) gFile->Get("hVertex");
272       if(!ver) AliWarning("no T0 histogram collected by PHYS DA ");
273       if(ver) {
274         meanver = ver->GetMean();
275         rmsver = ver->GetRMS();
276       }
277       SetMeanVertex(meanver);
278       SetRmsVertex(rmsver);
279       
280       gFile->Close();
281       delete gFile;
282
283     }
284     return ok; 
285    }
286
287 //________________________________________________________________________
288 void AliT0CalibTimeEq::GetMeanAndSigma(TH1F* hist,  Float_t &mean, Float_t &sigma) {
289   
290   const double window = 2.;  //fit window 
291   
292   double meanEstimate, sigmaEstimate; 
293   int maxBin;
294   maxBin        =  hist->GetMaximumBin(); //position of maximum
295   meanEstimate  =  hist->GetBinCenter( maxBin); // mean of gaussian sitting in maximum
296   sigmaEstimate = hist->GetRMS();
297   TF1* fit= new TF1("fit","gaus", meanEstimate - window*sigmaEstimate, meanEstimate + window*sigmaEstimate);
298   fit->SetParameters(hist->GetBinContent(maxBin), meanEstimate, sigmaEstimate);
299   hist->Fit("fit","RQ","Q");
300
301   mean  = (Float_t) fit->GetParameter(1);
302   sigma = (Float_t) fit->GetParameter(2);
303
304   delete fit;
305 }
306
307