]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0CalibTimeEq.cxx
Moving the classes that belong to the following libraries: STEERBase, ESD, CDB, AOD...
[u/mrichter/AliRoot.git] / T0 / AliT0CalibTimeEq.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /* $Id$ */
18
19 ///////////////////////////////////////////////////////////////////////////////
20 //                                                                           //
21 // class for T0 calibration                       TM-AC-AM_6-02-2006  
22 // equalize time shift for each time CFD channel
23 //                                                                           //
24 ///////////////////////////////////////////////////////////////////////////////
25
26 #include "AliT0CalibTimeEq.h"
27 #include "AliLog.h"
28 #include <TFile.h>
29 #include <TMath.h>
30 #include <TF1.h>
31 #include <TSpectrum.h>
32 #include <TProfile.h>
33 #include <iostream>
34
35 ClassImp(AliT0CalibTimeEq)
36
37 //________________________________________________________________
38   AliT0CalibTimeEq::AliT0CalibTimeEq():TNamed(),
39                                        fMeanVertex(0),        
40                                        fRmsVertex(0)      
41 {
42   //
43
44 }
45
46 //________________________________________________________________
47 AliT0CalibTimeEq::AliT0CalibTimeEq(const char* name):TNamed(),
48                                        fMeanVertex(0),        
49                                        fRmsVertex(0)      
50 {
51   //constructor
52
53   TString namst = "Calib_";
54   namst += name;
55   SetName(namst.Data());
56   SetTitle(namst.Data());
57 }
58
59 //________________________________________________________________
60 AliT0CalibTimeEq::AliT0CalibTimeEq(const AliT0CalibTimeEq& calibda):TNamed(calibda),            
61                                        fMeanVertex(0),        
62                                        fRmsVertex(0)      
63 {
64 // copy constructor
65   SetName(calibda.GetName());
66   SetTitle(calibda.GetName());
67
68
69 }
70
71 //________________________________________________________________
72 AliT0CalibTimeEq &AliT0CalibTimeEq::operator =(const AliT0CalibTimeEq& calibda)
73 {
74 // assignment operator
75   SetName(calibda.GetName());
76   SetTitle(calibda.GetName());
77  
78   return *this;
79 }
80
81 //________________________________________________________________
82 AliT0CalibTimeEq::~AliT0CalibTimeEq()
83 {
84   //
85   // destrictor
86 }
87 //________________________________________________________________
88 void AliT0CalibTimeEq::Reset()
89 {
90   //reset values
91
92   memset(fCFDvalue,0,120*sizeof(Float_t));
93   memset(fTimeEq,1,24*sizeof(Float_t));
94 }
95
96
97 //________________________________________________________________
98 void  AliT0CalibTimeEq::Print(Option_t*) const
99 {
100   // print time values
101
102   printf("\n    ----    PM Arrays       ----\n\n");
103   printf(" Time delay CFD \n");
104   for (Int_t i=0; i<24; i++) printf(" CFD  %f ",fTimeEq[i]);
105   printf("\n Mean Vertex %f \n", fMeanVertex);
106
107
108
109 //________________________________________________________________
110 Bool_t AliT0CalibTimeEq::ComputeOnlineParams(const char* filePhys)
111 {
112   // compute online equalized time
113   Float_t meandiff, sigmadiff, meanver, meancfdtime, sigmacfdtime;
114   meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
115   Double_t rmsver=0;
116   Int_t nent=0;
117   Bool_t ok=false;
118   gFile = TFile::Open(filePhys);
119   if(!gFile) {
120     AliError("No input PHYS data found ");
121   }
122   else
123     {
124       ok=true;
125       for (Int_t i=0; i<24; i++)
126         {
127
128           meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
129           TH1F *cfd = (TH1F*) gFile->Get(Form("CFD1minCFD%d",i+1));
130           TH1F *cfdtime = (TH1F*) gFile->Get(Form("CFD%d",i+1));
131           if(!cfd) {
132             AliWarning(Form("no Diff histograms collected by PHYS DA for channel %i", i));
133           }
134           if(!cfdtime) {
135             AliWarning(Form("no CFD histograms collected by PHYS DA for channel %i", i));
136           }
137           if(cfd) {
138             nent = Int_t(cfd->GetEntries());
139             if(nent>50)  { 
140               if(cfd->GetRMS()>1.5 )
141                 GetMeanAndSigma(cfd, meandiff, sigmadiff);
142               if(cfd->GetRMS()<=1.5) 
143                 {
144                   meandiff = cfd->GetMean();
145                   sigmadiff=cfd->GetRMS();
146                 }
147               Int_t   maxBin = cfd->GetMaximumBin(); 
148               Double_t  meanEstimate = cfd->GetBinCenter( maxBin); 
149               if(TMath::Abs(meanEstimate - meandiff) > 20 ) meandiff = meanEstimate; 
150             }
151             else 
152               {
153                 //      ok=false;
154                 AliWarning(Form(" Not  enouph data in PMT %i- PMT1:  %i ", i, nent));
155               }
156           }
157           if(cfdtime) {
158             nent = Int_t(cfdtime->GetEntries());
159             if(nent > 50 )  { //!!!!!!!!!!
160               if(cfdtime->GetRMS()>1.5 )
161                 GetMeanAndSigma(cfdtime,meancfdtime, sigmacfdtime);
162               if(cfdtime->GetRMS()<=1.5) 
163                 {
164                   meancfdtime = cfdtime->GetMean();
165                   sigmacfdtime = cfdtime->GetRMS();
166                 }
167               Int_t   maxBin = cfdtime->GetMaximumBin(); 
168               Double_t  meanEstimate = cfdtime->GetBinCenter( maxBin); 
169               if(TMath::Abs(meanEstimate - meancfdtime) > 20 ) meancfdtime = meanEstimate; 
170             }
171             else 
172               {
173                 //      ok=false;
174                 AliWarning(Form(" Not  enouph data in PMT in CFD peak %i - %i ", i, nent));
175               }
176           }
177           SetTimeEq(i,meandiff);
178           SetTimeEqRms(i,sigmadiff);
179           SetCFDvalue(i,0,meancfdtime);
180           if (cfd) delete cfd;
181           if (cfdtime) delete cfdtime;
182           
183         }
184       TH1F *ver = (TH1F*) gFile->Get("hVertex");
185       if(ver) {
186         meanver = ver->GetMean();
187         rmsver = ver->GetRMS();
188       }
189       SetMeanVertex(meanver);
190       SetRmsVertex(rmsver);
191       
192       gFile->Close();
193       delete gFile;
194       
195     }
196     return ok; 
197 }
198
199 //________________________________________________________________
200   Bool_t AliT0CalibTimeEq::ComputeOfflineParams(const char* filePhys, Float_t *timecdb, Float_t *cfdvalue, Int_t badpmt)
201 {
202   // compute offline equalized time
203   Float_t meandiff, sigmadiff, meanver, meancfdtime, sigmacfdtime;
204   meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
205   Int_t nent=0;
206   Bool_t ok=false;
207   TH1F *cfddiff = NULL; 
208   TH1F *cfdtime = NULL;
209   TObjArray * TzeroObj = NULL;
210
211   gFile = TFile::Open(filePhys);
212   if(!gFile) {
213     AliError("No input PHYS data found ");
214     return ok;
215   }
216   else
217     {
218       meandiff = sigmadiff =  meanver = meancfdtime = sigmacfdtime =0;
219       ok=true;
220       TDirectory *dr = (TDirectory*) gFile->Get("T0Calib");
221       if (dr)   TzeroObj = (TObjArray*) dr->Get("T0Calib");
222       
223       for (Int_t i=0; i<24; i++)
224         {
225           if (i != badpmt) {
226             if(TzeroObj) {
227               cfddiff = (TH1F*)TzeroObj->FindObject(Form("CFD1minCFD%d",i+1));
228               cfdtime = (TH1F*)TzeroObj->FindObject(Form("CFD%d",i+1));
229             }
230             else
231               {
232                 cfddiff = (TH1F*)gFile->Get(Form("CFD1minCFD%d",i+1));
233                 cfdtime = (TH1F*)gFile->Get(Form("CFD%d",i+1));
234                 
235               }
236             if(!cfddiff ) {
237               AliWarning(Form("no  histograms collected by pass0 for diff channel %i", i));
238               SetTimeEq(i,timecdb[i]);
239               SetTimeEqRms(i,0);
240             }
241             if(cfddiff) {
242               nent = Int_t(cfddiff->GetEntries());
243               if(nent>100 )  { //!!!!!
244                 if(cfddiff->GetRMS()>1.5 )
245                   GetMeanAndSigma(cfddiff, meandiff, sigmadiff);
246                 if(cfddiff->GetRMS()<=1.5) 
247                   {
248                     meandiff = cfddiff->GetMean();
249                     sigmadiff = cfddiff->GetRMS();
250                   }
251                 Int_t   maxBin = cfddiff->GetMaximumBin(); 
252                 Double_t  meanEstimate = cfddiff->GetBinCenter( maxBin); 
253                 if(TMath::Abs(meanEstimate - meandiff) > 20 ) meandiff = meanEstimate;        
254               }
255               else 
256                 {
257                   AliWarning(Form(" Not  enouph data in PMT %i- PMT1:  %i ", i, nent));
258                   SetTimeEq(i,timecdb[i]);
259                   SetTimeEqRms(i,0);
260                   
261                 }
262             }       
263             
264             if(!cfdtime ) {
265               AliWarning(Form("no  histograms collected by pass0 for time channel %i", i));
266               SetTimeEq(i, cfdvalue[i]);
267               SetTimeEqRms(i, 0);
268             }
269             if(cfdtime) {
270               nent = Int_t(cfdtime->GetEntries());
271               if(nent>100 )  { //!!!!!
272                 if(cfdtime->GetRMS()>1.5 )
273                   GetMeanAndSigma(cfdtime,meancfdtime, sigmacfdtime);
274                 if(cfdtime->GetRMS()<=1.5) 
275                   {
276                     meancfdtime = cfdtime->GetMean();
277                     sigmacfdtime=cfdtime->GetRMS();
278                 }
279                 Int_t   maxBin = cfdtime->GetMaximumBin(); 
280                 Double_t  meanEstimate = cfdtime->GetBinCenter( maxBin); 
281                 if(TMath::Abs(meanEstimate - meancfdtime) > 20 ) meancfdtime = meanEstimate; 
282             }
283             else 
284               {
285                 AliWarning(Form(" Not  enouph data in PMT in CFD peak %i - %i ", i, nent));
286                 SetTimeEq(i, cfdvalue[i]);
287                 SetTimeEqRms(i, 0);
288               }
289           }
290           
291           SetTimeEq(i,timecdb[i] + meandiff);
292           SetTimeEqRms(i,sigmadiff);
293           SetCFDvalue(i,0, meancfdtime );
294           if (cfddiff) cfddiff->Reset();
295           if (cfdtime) cfdtime->Reset();
296           } //bad pmt
297         }
298       
299       gFile->Close();
300       delete gFile;
301
302     }
303     return ok; 
304    }
305
306 //________________________________________________________________________
307 void AliT0CalibTimeEq::GetMeanAndSigma(TH1F* hist,  Float_t &mean, Float_t &sigma) {
308   
309   const double window = 2.;  //fit window 
310   
311   double meanEstimate, sigmaEstimate; 
312   int maxBin;
313   maxBin        =  hist->GetMaximumBin(); //position of maximum
314   meanEstimate  =  hist->GetBinCenter( maxBin); // mean of gaussian sitting in maximum
315   // sigmaEstimate = hist->GetRMS();
316   sigmaEstimate = 10;
317   TF1* fit= new TF1("fit","gaus", meanEstimate - window*sigmaEstimate, meanEstimate + window*sigmaEstimate);
318   fit->SetParameters(hist->GetBinContent(maxBin), meanEstimate, sigmaEstimate);
319   hist->Fit("fit","RQ","Q");
320
321   mean  = (Float_t) fit->GetParameter(1);
322   sigma = (Float_t) fit->GetParameter(2);
323
324   delete fit;
325 }
326
327