]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFcalib.cxx
Rearrangement of Calibration objects for simulation
[u/mrichter/AliRoot.git] / TOF / AliTOFcalib.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 /*
17 $Log$
18 Revision 1.16  2007/10/08 10:13:26  zampolli
19 First Run and Last Run members added, infinite validity of calib obj implemented.
20
21 Revision 1.15  2007/10/04 13:23:28  zampolli
22 Updates to handle functionalities in TOF online/offline calibration according to the latest schema
23
24 Revision 1.14  2007/06/06 16:26:30  arcelli
25 remove fall-back call to local CDB storage
26
27 Revision 1.13  2007/04/20 13:59:40  arcelli
28 make protections agains failed retrieval of the CDB object in a proper way
29
30 Revision 1.12  2007/03/23 11:31:16  arcelli
31 CDB Entry for TOF Reconstruction Parameters
32
33 Revision 1.11  2007/02/28 18:08:26  arcelli
34 Add protection against failed retrieval of the CDB cal object
35
36 Revision 1.10  2006/08/22 13:30:49  arcelli
37 removal of effective c++ warnings (C.Zampolli)
38
39 Revision 1.9  2006/04/20 22:30:50  hristov
40 Coding conventions (Annalisa)
41
42 Revision 1.8  2006/04/16 22:29:05  hristov
43 Coding conventions (Annalisa)
44
45 Revision 1.7  2006/04/16 20:12:46  hristov
46 Removing memory leak in case of cached CDB entries
47
48 Revision 1.6  2006/04/11 15:28:32  hristov
49 Checks on cache status before deleting calibration objects (A.Colla)
50
51 Revision 1.5  2006/04/05 08:35:38  hristov
52 Coding conventions (S.Arcelli, C.Zampolli)
53
54 Revision 1.4  2006/03/31 11:26:46  arcelli
55  changing CDB Ids according to standard convention
56
57 Revision 1.3  2006/03/28 14:57:02  arcelli
58 updates to handle new V5 geometry & some re-arrangements
59
60 Revision 1.2  2006/02/13 17:22:26  arcelli
61 just Fixing Log info
62
63 Revision 1.1  2006/02/13 16:10:48  arcelli
64 Add classes for TOF Calibration (C.Zampolli)
65
66 author: Chiara Zampolli, zampolli@bo.infn.it
67 */  
68
69 ///////////////////////////////////////////////////////////////////////////////
70 //                                                                           //
71 // class for TOF calibration                                                 //
72 //                                                                           //
73 ///////////////////////////////////////////////////////////////////////////////
74
75 #include "TF1.h"
76 #include "TFile.h"
77 #include "TH1F.h"
78 #include "TH2F.h"
79 #include "TList.h"
80 #include "TROOT.h"
81 #include "TStyle.h"
82 #include "TTree.h"
83 #include "TProfile.h"
84 #include "TGrid.h"
85
86 #include "AliCDBEntry.h"
87 #include "AliCDBRunRange.h"
88 #include "AliCDBId.h"
89 #include "AliCDBManager.h"
90 #include "AliCDBStorage.h"
91 #include "AliCDBMetaData.h"
92 #include "AliESDtrack.h"
93 #include "AliESD.h"
94 #include "AliLog.h"
95
96 #include "AliTOFcalib.h"
97 #include "AliTOFChannelOnline.h"
98 #include "AliTOFChannelOffline.h"
99 #include "AliTOFGeometry.h"
100 #include "AliTOFRecoParam.h"
101
102 extern TROOT *gROOT;
103 extern TStyle *gStyle;
104
105 ClassImp(AliTOFcalib)
106
107 //_______________________________________________________________________
108 AliTOFcalib::AliTOFcalib():
109   TTask("AliTOFcalib",""),
110   fNChannels(-1),
111   fTOFCalOnline(0x0),
112   fTOFCalOffline(0x0),
113   fTOFSimToT(0x0),
114   fkValidity(0x0),
115   fTree(0x0),
116   fNruns(0),
117   fFirstRun(0),
118   fLastRun(AliCDBRunRange::Infinity())
119
120   //TOF Calibration Class ctor
121   fNChannels = AliTOFGeometry::NSectors()*(2*(AliTOFGeometry::NStripC()+AliTOFGeometry::NStripB())+AliTOFGeometry::NStripA())*AliTOFGeometry::NpadZ()*AliTOFGeometry::NpadX();
122 }
123 //____________________________________________________________________________ 
124
125 AliTOFcalib::AliTOFcalib(const AliTOFcalib & calib):
126   TTask("AliTOFcalib",""),
127   fNChannels(calib.fNChannels),
128   fTOFCalOnline(0x0),
129   fTOFCalOffline(0x0),
130   fTOFSimToT(calib.fTOFSimToT),
131   fkValidity(calib.fkValidity),
132   fTree(calib.fTree),
133   fNruns(calib.fNruns),
134   fFirstRun(calib.fFirstRun),
135   fLastRun(calib.fLastRun)
136 {
137   //TOF Calibration Class copy ctor
138   for (Int_t iarray = 0; iarray<fNChannels; iarray++){
139     AliTOFChannelOnline * calChOnline = (AliTOFChannelOnline*)calib.fTOFCalOnline->At(iarray);
140     AliTOFChannelOffline * calChOffline = (AliTOFChannelOffline*)calib.fTOFCalOffline->At(iarray);
141     fTOFCalOnline->AddAt(calChOnline,iarray);
142     fTOFCalOffline->AddAt(calChOffline,iarray);
143
144   }
145 }
146
147 //____________________________________________________________________________ 
148
149 AliTOFcalib& AliTOFcalib::operator=(const AliTOFcalib &calib)
150 {
151   //TOF Calibration Class assignment operator
152   this->fNChannels = calib.fNChannels;
153   this->fTOFSimToT = calib.fTOFSimToT;
154   this->fkValidity = calib.fkValidity;
155   this->fTree = calib.fTree;
156   this->fNruns = calib.fNruns;
157   this->fFirstRun = calib.fFirstRun;
158   this->fLastRun = calib.fLastRun;
159   for (Int_t iarray = 0; iarray<fNChannels; iarray++){
160     AliTOFChannelOnline * calChOnline = (AliTOFChannelOnline*)calib.fTOFCalOnline->At(iarray);
161     AliTOFChannelOffline * calChOffline = (AliTOFChannelOffline*)calib.fTOFCalOffline->At(iarray);
162     this->fTOFCalOnline->AddAt(calChOnline,iarray);
163     this->fTOFCalOffline->AddAt(calChOffline,iarray);
164   }
165   return *this;
166 }
167
168 //____________________________________________________________________________ 
169
170 AliTOFcalib::~AliTOFcalib()
171 {
172   //TOF Calibration Class dtor
173   if(!(AliCDBManager::Instance()->GetCacheFlag())){ // CDB objects must NOT be deleted if cache is active!
174     if (fTOFCalOnline){
175       //      fTOFCalOnline->Clear();
176       delete fTOFCalOnline;
177     }
178     if (fTOFCalOffline){
179       //fTOFCalOffline->Clear();
180       delete fTOFCalOffline;
181     }
182   }
183   if (fTree!=0x0) delete fTree;
184 }
185 //_____________________________________________________________________________
186 void AliTOFcalib::CreateCalArrays(){
187
188   // creating arrays for online/offline calibration objs
189
190   fTOFCalOnline = new TObjArray(fNChannels);
191   fTOFCalOffline = new TObjArray(fNChannels);
192   fTOFCalOnline->SetOwner();
193   fTOFCalOffline->SetOwner();
194   for (Int_t iarray = 0; iarray<fNChannels; iarray++){
195     AliTOFChannelOnline * calChOnline = new AliTOFChannelOnline();
196     AliTOFChannelOffline * calChOffline = new AliTOFChannelOffline();
197     fTOFCalOnline->AddAt(calChOnline,iarray);
198     fTOFCalOffline->AddAt(calChOffline,iarray);
199   }
200 }
201 //_____________________________________________________________________________
202 void AliTOFcalib::WriteParOnlineOnCDB(Char_t *sel, Int_t minrun, Int_t maxrun)
203 {
204   //Write calibration parameters to the CDB
205   SetFirstRun(minrun);
206   SetLastRun(maxrun);
207   AliCDBManager *man = AliCDBManager::Instance();
208   Char_t *sel1 = "ParOnline" ;  // to be consistent with TOFPreprocessor
209   Char_t  out[100];
210   sprintf(out,"%s/%s",sel,sel1); 
211   AliDebug(2,Form("Writing TOF online calib obj on CDB with run range [%i, %i] ",fFirstRun,fLastRun));
212   AliCDBId id(out,fFirstRun,fLastRun);
213   AliCDBMetaData *md = new AliCDBMetaData();
214   md->SetResponsible("Chiara Zampolli");
215   if (!fTOFCalOnline) {
216     // deve uscire!!
217   }
218   man->Put(fTOFCalOnline,id,md);
219   delete md;
220 }
221 //_____________________________________________________________________________
222
223 void AliTOFcalib::WriteParOnlineOnCDB(Char_t *sel)
224 {
225   //Write calibration parameters to the CDB with infinite validity
226   AliCDBManager *man = AliCDBManager::Instance();
227   Char_t *sel1 = "ParOnline" ;  // to be consistent with TOFPreprocessor
228   Char_t  out[100];
229   sprintf(out,"%s/%s",sel,sel1); 
230   AliCDBRunRange runrange(fFirstRun,fLastRun);
231   AliDebug(2,Form("Writing TOF online calib obj on CDB with run range [%i, %i] ",runrange.GetFirstRun(),runrange.GetLastRun()));
232   AliCDBId id(out,runrange);
233   AliCDBMetaData *md = new AliCDBMetaData();
234   md->SetResponsible("Chiara Zampolli");
235   if (!fTOFCalOnline) {
236     // deve uscire!!
237   }
238   man->Put(fTOFCalOnline,id,md);
239   delete md;
240 }
241 //_____________________________________________________________________________
242
243 void AliTOFcalib::WriteParOfflineOnCDB(Char_t *sel, const Char_t *validity, Int_t minrun, Int_t maxrun)
244 {
245   //Write calibration parameters to the CDB
246   SetFirstRun(minrun);
247   SetLastRun(maxrun);
248   AliCDBManager *man = AliCDBManager::Instance();
249   Char_t *sel1 = "ParOffline" ;
250   Char_t  out[100];
251   sprintf(out,"%s/%s",sel,sel1); 
252   AliDebug(2,Form("Writing TOF offline calib obj on CDB with run range [%i, %i] ",fFirstRun,fLastRun));
253   AliCDBId id(out,fFirstRun,fLastRun);
254   AliCDBMetaData *md = new AliCDBMetaData();
255   md->SetResponsible("Chiara Zampolli");
256   md->SetComment(validity);
257   man->Put(fTOFCalOffline,id,md);
258   delete md;
259 }
260 //_____________________________________________________________________________
261
262 void AliTOFcalib::WriteParOfflineOnCDB(Char_t *sel, const Char_t *validity)
263 {
264   //Write calibration parameters to the CDB with infinite validity
265   AliCDBManager *man = AliCDBManager::Instance();
266   Char_t *sel1 = "ParOffline" ;
267   Char_t  out[100];
268   sprintf(out,"%s/%s",sel,sel1); 
269   AliCDBRunRange runrange(fFirstRun,fLastRun);
270   AliDebug(2,Form("Writing TOF offline calib obj on CDB with run range [%i, %i] ",runrange.GetFirstRun(),runrange.GetLastRun()));
271   AliCDBId id(out,runrange);
272   AliCDBMetaData *md = new AliCDBMetaData();
273   md->SetResponsible("Chiara Zampolli");
274   md->SetComment(validity);
275   man->Put(fTOFCalOffline,id,md);
276   delete md;
277 }
278 //_____________________________________________________________________________
279
280 Bool_t AliTOFcalib::ReadParOnlineFromCDB(Char_t *sel, Int_t nrun)
281 {
282   //Read calibration parameters from the CDB
283   AliCDBManager *man = AliCDBManager::Instance();
284   Char_t *sel1 = "ParOnline" ;
285   Char_t  out[100];
286   sprintf(out,"%s/%s",sel,sel1); 
287   if (!man->Get(out,nrun)) { 
288     return kFALSE;
289   }
290   AliCDBEntry *entry = man->Get(out,nrun);
291   if(!entry->GetObject()){
292     return kFALSE;
293   }  
294   
295   fTOFCalOnline =(TObjArray*)entry->GetObject();
296
297   return kTRUE; 
298    
299 }
300 //_____________________________________________________________________________
301
302 Bool_t AliTOFcalib::ReadParOfflineFromCDB(Char_t *sel, Int_t nrun)
303 {
304   //Read calibration parameters from the CDB
305   AliCDBManager *man = AliCDBManager::Instance();
306   Char_t *sel1 = "ParOffline" ;
307   Char_t  out[100];
308   sprintf(out,"%s/%s",sel,sel1); 
309   if (!man->Get(out,nrun)) { 
310     return kFALSE;
311   }
312   AliCDBEntry *entry = man->Get(out,nrun);
313   if(!entry->GetObject()){
314     return kFALSE;
315   }  
316   AliCDBMetaData * md = entry->GetMetaData();
317   fkValidity = md->GetComment();  
318   fTOFCalOffline =(TObjArray*)entry->GetObject();
319
320   return kTRUE; 
321    
322 }
323 //_____________________________________________________________________________
324 void AliTOFcalib::WriteSimHistoOnCDB(Char_t *sel, Int_t minrun, Int_t maxrun, TH1F *histo){
325   //Write Sim miscalibration parameters to the CDB
326
327   fTOFSimToT=histo;
328   AliCDBManager *man = AliCDBManager::Instance();
329   Char_t *sel1 = "SimHisto" ;
330   Char_t  out[100];
331   sprintf(out,"%s/%s",sel,sel1); 
332   AliCDBMetaData *mdhisto = new AliCDBMetaData();
333   mdhisto->SetResponsible("Chiara Zampolli");
334   AliCDBId id(out,minrun,maxrun);
335   man->Put(fTOFSimToT,id,mdhisto);
336   delete mdhisto;
337 }
338 //_____________________________________________________________________________
339 Bool_t AliTOFcalib::ReadSimHistoFromCDB(Char_t *sel, Int_t nrun)
340 {
341   //Read miscalibration parameters from the CDB
342   AliCDBManager *man = AliCDBManager::Instance();
343
344   // The Tot Histo
345
346   Char_t *sel1 = "SimHisto" ;
347   Char_t  out[100];
348   sprintf(out,"%s/%s",sel,sel1); 
349   if (!man->Get(out,nrun)) { 
350     AliFatal("Exiting, no CDB object (SimHisto) found!!!");
351     exit(0);  
352   }
353   AliCDBEntry *entry = man->Get(out,nrun);
354   if(!entry->GetObject()){
355     AliFatal("Exiting, no CDB object (SimHisto) found!!!");
356     exit(0);  
357   }  
358   TH1F *histo =(TH1F*)entry->GetObject();
359   fTOFSimToT=histo;
360 }
361 //_____________________________________________________________________________
362 void AliTOFcalib::WriteRecParOnCDB(Char_t *sel, Int_t minrun, Int_t maxrun, AliTOFRecoParam *param){
363   //Write reconstruction parameters to the CDB
364
365   AliCDBManager *man = AliCDBManager::Instance();
366   AliCDBMetaData *md = new AliCDBMetaData();
367   md->SetResponsible("Silvia Arcelli");
368   Char_t *sel1 = "RecPar" ;
369   Char_t  out[100];
370   sprintf(out,"%s/%s",sel,sel1); 
371   AliCDBId id(out,minrun,maxrun);
372   man->Put(param,id,md);
373   delete md;
374 }
375 //_____________________________________________________________________________
376 AliTOFRecoParam * AliTOFcalib::ReadRecParFromCDB(Char_t *sel, Int_t nrun)
377 {
378   //Read reconstruction parameters from the CDB
379   AliCDBManager *man = AliCDBManager::Instance();
380   Char_t *sel1 = "RecPar" ;
381   Char_t  out[100];
382   sprintf(out,"%s/%s",sel,sel1); 
383   if (!man->Get(out,nrun)) { 
384     AliFatal("Exiting, no CDB object (RecPar) found!!!");
385     exit(0);  
386   }  
387   AliCDBEntry *entry = man->Get(out,nrun);
388   if(!entry->GetObject()){
389     AliFatal("Exiting, no CDB object (RecPar) found!!!");
390     exit(0);  
391   }  
392
393   AliTOFRecoParam *param=(AliTOFRecoParam*)entry->GetObject();
394   return param;
395 }
396 //-----------------------------------------------------------------------------
397 // Calibration methods
398 //-----------------------------------------------------------------------------
399 void AliTOFcalib::CreateTreeFromCDB(Int_t minrun, Int_t maxrun){
400
401   // creating the chain with the trees for calibration
402   // collecting them from reference data 
403   // from minrun to maxrun
404
405   Float_t p[CHENTRIESSMALL];
406   Int_t nentries;
407   fTree = new TTree("TOFCalib","Tree for TOF Calibration");
408   fTree->Branch("nentries",&nentries,"nentries/I");
409   fTree->Branch("TOFentries",p,"TOFentries[nentries]/F");
410   AliCDBManager *man = AliCDBManager::Instance();
411   AliCDBStorage *aStorage = man->GetStorage("local://$ALICE_ROOT");
412   for (Int_t irun = minrun;irun<=maxrun;irun++){
413     AliCDBEntry *entry = aStorage->Get("TOF/RefData/TreeForCalib",irun);
414     if (!entry){
415       AliInfo(Form("No entry found for run %i",irun));
416     }
417     else{
418       TTree *tree = new TTree();
419       tree = (TTree*)entry->GetObject();
420       tree->SetBranchAddress("nentries",&nentries);
421       tree->SetBranchAddress("TOFentries",p);      
422       fTree->CopyEntries(tree);
423       delete tree;
424       fNruns++;
425     }
426   }
427   AliInfo(Form("Number of runs being analyzed %i",fNruns));
428 }
429 //-----------------------------------------------------------------------------
430 void AliTOFcalib::CreateTreeFromGrid(Int_t minrun, Int_t maxrun){
431
432   // creating the chain with the trees for calibration
433   // collecting them from the Grid 
434   // from minrun to maxrun
435
436   Float_t p[CHENTRIESSMALL];
437   Int_t nentries;
438   fTree = new TTree("TOFCalib","Tree for TOF Calibration");
439   fTree->SetDirectory(0);
440   fTree->Branch("nentries",&nentries,"nentries/I");
441   fTree->Branch("TOFentries",p,"TOFentries[nentries]/F");
442   AliInfo("connected to alien");
443   TGrid::Connect("alien://");
444   
445   Char_t filename[100];
446   for (Int_t irun = minrun;irun<=maxrun;irun++){
447     sprintf(filename,"alien:///alice/cern.ch/user/c/czampolli/TOFCalibReference_%i.root",irun);
448     TFile *filegrid = TFile::Open(filename,"READ");
449     TTree *tree = (TTree*)filegrid->Get("T");
450     tree->SetBranchAddress("nentries",&nentries);
451     tree->SetBranchAddress("TOFentries",p);      
452     fTree->CopyEntries(tree);
453     delete tree;
454     fNruns++;    
455   }
456   
457   AliInfo(Form("Number of runs being analyzed %i",fNruns));
458 }
459 //-----------------------------------------------------------------------------
460 void AliTOFcalib::CreateTreeFromFile(Int_t minrun, Int_t maxrun){
461
462   // creating the tree with the trees for calibration
463   // collecting them from reference data (from file)
464   // from minrun to maxrun
465
466   Float_t p[CHENTRIESSMALL];
467   Int_t nentries;
468   fTree = new TTree("TOFCalib","Tree for TOF Calibration");
469   fTree->SetDirectory(0);
470   fTree->Branch("nentries",&nentries,"nentries/I");
471   fTree->Branch("TOFentries",p,"TOFentries[nentries]/F");
472   Char_t filename[100];
473   for (Int_t irun = minrun;irun<=maxrun;irun++){
474     sprintf(filename,"$ALICE_ROOT/TOF/RefData/TreeForCalib/fileout_%i.root",irun);
475     TFile *file = new TFile(filename,"READ");
476     TTree *tree = (TTree*)file->Get("T");
477     tree->SetBranchAddress("nentries",&nentries);
478     tree->SetBranchAddress("TOFentries",p);      
479     fTree->CopyEntries(tree);
480     delete tree;
481     delete file;
482     file = 0x0;
483     fNruns++;
484   }
485
486   AliInfo(Form("Number of runs being analyzed %i",fNruns));
487 }
488 //-----------------------------------------------------------------------------
489 Int_t AliTOFcalib::Calibrate(Int_t ichmin, Int_t ichmax, Option_t *optionSave, Option_t *optionFit){
490
491   // calibrating summing more than one channels
492   // computing calibration parameters
493   // Returning codes:
494   // 0 -> everything was ok
495   // 1 -> no tree for calibration found
496   // 2 -> not enough statistics to perform calibration
497   // 3 -> problems with arrays
498   
499   TH1::AddDirectory(0);
500
501   AliInfo(Form("*** Calibrating Histograms %s, summing more channels, from channel %i, to channel %i, storing Calib Pars in channel %i", GetName(),ichmin,ichmax,ichmin)) ; 
502   AliInfo(Form("Option for Saving histos = %s",optionSave )) ; 
503   AliInfo(Form("Option for Fitting Profile histos = %s",optionFit )) ; 
504
505   Float_t p[CHENTRIESSMALL];
506   Int_t nentries;
507   fTree->SetBranchAddress("nentries",&nentries);
508   fTree->SetBranchAddress("TOFentries",p);
509
510   Float_t ntracksTotalmean =0;
511   for (Int_t i=ichmin; i<ichmax; i++){
512     Int_t ientry = -1;
513     for (Int_t irun=0;irun<fNruns;irun++){
514       ientry = i+irun*fNChannels;
515       fTree->GetEntry(ientry);
516       Int_t ntracksRun=nentries/3;
517       ntracksTotalmean+=ntracksRun;
518     }
519   }
520   
521   if (ntracksTotalmean < MEANENTRIES) {
522     AliInfo(Form(" Too small mean number of entires per channel (mean number = %f) not calibrating and exiting.....",ntracksTotalmean));
523     return 2;
524   }
525
526   //filling ToT and Time arrays
527
528   Int_t nbinToT = 100;  // ToT bin width in Profile = 48.8 ps 
529   Float_t minToT = 0;   // ns
530   Float_t maxToT = 4.88;  // ns
531
532   TH1F *hToT = new TH1F("htot","htot",nbinToT, minToT, maxToT);
533   TH1F *hdeltaTime = new TH1F("hdeltaTime","hdeltaTime",200,2,4);
534   Int_t ntracksTotal = 0;
535   Int_t ntracksRun = 0;
536   Double_t binsProfile[101]; // sized larger than necessary, the correct 
537                              // dim being set in the booking of the profile
538   Int_t nusefulbins=0;
539   Float_t meantime=0;
540   for (Int_t i = ichmin;i<ichmax;i++){
541     Int_t ientry = -1;
542     for (Int_t irun=0;irun<fNruns;irun++){
543       ientry = i+irun*fNChannels;
544       fTree->GetEntry(ientry);
545       ntracksTotal+=nentries/3;
546       ntracksRun=nentries/3;
547       AliDebug(2,Form("run %i, channel %i, nentries = %i, ntracks = %i",irun,i,nentries, ntracksRun));
548       for (Int_t j=0;j<ntracksRun;j++){
549         Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
550         Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
551         Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
552         Float_t tot = p[idxexToT];
553         hdeltaTime->Fill(p[idxexTime]-p[idxexExTime]);
554         meantime+=p[idxexTime]-p[idxexExTime];
555         hToT->Fill(tot);
556       }
557     }
558   }
559   nusefulbins = FindBins(hToT,&binsProfile[0]);
560   meantime/=ntracksTotal;
561   AliDebug(2, Form("meantime = %f",meantime));
562   
563   for (Int_t j=1;j<=nusefulbins;j++) {
564     AliDebug(2,Form(" summing channels from %i to %i, nusefulbins = %i, bin %i = %f",ichmin,ichmax,nusefulbins,j,binsProfile[j])); 
565   }
566
567   TProfile* hSlewingProf = new TProfile("hSlewingProf", "hSlewingProf",nusefulbins, binsProfile, "G");  // CHECK THE BUILD OPTION, PLEASE!!!!!!
568   TH2F * htimetot = new TH2F("htimetot","htimetot",nbinToT, minToT, maxToT,600,-5,10);
569
570   for (Int_t irun=0;irun<fNruns;irun++){
571     Int_t ientry = -1;
572     for (Int_t i=ichmin; i<ichmax; i++){
573       ientry = i+irun*fNChannels;
574       fTree->GetEntry(ientry);
575       ntracksRun=nentries/3;
576       for (Int_t j=0;j<ntracksRun;j++){
577         Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
578         Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
579         Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
580         Float_t tot = p[idxexToT];
581         Float_t time = p[idxexTime]-p[idxexExTime];
582         AliDebug (2, Form("track = %i, time = %f, tot = %f, time-meantime = %f",j,time, tot, time-meantime));
583         hSlewingProf->Fill(tot,time);  // if meantime is not used, the fill may be moved in the loop above
584         htimetot->Fill(tot,time-meantime);  // if meantime is not used, the fill may be moved in the loop above
585       }
586     }
587   }
588
589   hSlewingProf->Fit("pol5",optionFit,"",0,4);
590   TF1 * calibfunc = (TF1*)hSlewingProf->GetFunction("pol5");
591   Float_t par[6];    
592   for(Int_t kk=0;kk<6;kk++){
593     par[kk]=calibfunc->GetParameter(kk);
594     AliDebug(2,Form("parameter %i = %f",kk,par[kk]));
595   }
596
597   if(strstr(optionSave,"save")){
598     TFile * fileProf = new TFile("TOFCalibSave.root","recreate");
599     fileProf->cd(); 
600     TString profName=Form("Profile%06i_%06i",ichmin,ichmax);
601     TString timeTotName=Form("TimeTot%06i_%06i",ichmin,ichmax);
602     TString totName=Form("Tot%06i_%06i",ichmin,ichmax);
603     TString deltaName=Form("Delta%06i_%06i",ichmin,ichmax);
604     hSlewingProf->Write(profName);
605     htimetot->Write(timeTotName);
606     hToT->Write(totName);
607     hdeltaTime->Write(deltaName);
608     fileProf->Close();
609     delete fileProf;
610     fileProf=0x0;
611   }
612
613   delete hToT;
614   hToT=0x0;
615   delete hSlewingProf;
616   hSlewingProf=0x0;
617   delete htimetot;
618   htimetot=0x0;
619   delete hdeltaTime;
620   hdeltaTime=0x0;
621
622   AliTOFChannelOffline * calChannel = (AliTOFChannelOffline*)fTOFCalOffline->At(ichmin);
623   calChannel->SetSlewPar(par);
624   WriteParOfflineOnCDB("TOF/Calib","valid");
625   return 0;
626 }
627 //----------------------------------------------------------------------------
628 Int_t AliTOFcalib::Calibrate(Int_t i, Option_t *optionSave, Option_t *optionFit){
629
630   // computing calibration parameters for channel i
631   // Returning codes:
632   // 0 -> everything was ok
633   // 1 -> no tree for calibration found
634   // 2 -> not enough statistics to perform calibration
635   // 3 -> problems with arrays
636
637   TH1::AddDirectory(0);
638   
639   AliInfo(Form("*** Calibrating Histograms (one channel) %s", GetName())) ; 
640   AliInfo(Form("Option for Saving histos = %s",optionSave )) ; 
641   AliInfo(Form("Option for Fitting Profile histos = %s",optionFit )) ; 
642
643   Float_t p[MAXCHENTRIESSMALL];
644   Int_t nentries;
645   fTree->SetBranchAddress("nentries",&nentries);
646   fTree->SetBranchAddress("TOFentries",p);
647
648   Float_t ntracksTotal =0;
649   for (Int_t irun=0;irun<fNruns;irun++){
650     Int_t ientry = -1;
651     ientry = i+irun*fNChannels;
652     fTree->GetEntry(ientry);
653     ntracksTotal+=nentries/3;    
654   }
655   
656   if (ntracksTotal < MEANENTRIES) {  
657     AliInfo(Form(" Too small mean number of entires per channel (mean number = %f) not calibrating and exiting.....",ntracksTotal));
658     return 2;
659   }
660
661   //filling ToT and Time arrays
662
663   Int_t nbinToT = 100;  // ToT bin width in Profile = 48.8 ps 
664   Float_t minToT = 0;   // ns
665   Float_t maxToT = 4.88;  // ns
666
667   TH1F *hToT = new TH1F("htot","htot",nbinToT, minToT, maxToT);
668   TH1F *hdeltaTime = new TH1F("hdeltaTime","hdeltaTime",200,2,4);
669   Int_t ntracksRun = 0;
670   Double_t binsProfile[101]; // sized larger than necessary, the correct 
671                              // dim being set in the booking of the profile
672   Int_t nusefulbins=0;
673   Float_t meantime=0;
674   for (Int_t irun=0;irun<fNruns;irun++){
675     Int_t ientry = -1;
676     ientry = i+irun*fNChannels;
677     fTree->GetEntry(ientry);
678     ntracksRun=nentries/3;
679     AliDebug(2,Form("run %i, channel %i, nentries = %i, ntracksRun = %i",irun, i ,nentries, ntracksRun));
680     for (Int_t j=0;j<ntracksRun;j++){
681       Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
682       Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
683       Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
684       Float_t tot = p[idxexToT];
685       meantime+=p[idxexTime]-p[idxexExTime];
686       hdeltaTime->Fill(p[idxexTime]-p[idxexExTime]);
687       hToT->Fill(tot);
688     }
689   }
690
691   nusefulbins = FindBins(hToT,&binsProfile[0]);
692   meantime/=ntracksTotal;
693   AliDebug(2,Form("meantime = %f",meantime));
694   
695   for (Int_t j=1;j<=nusefulbins;j++) {
696     AliDebug(2,Form(" channel %i, nusefulbins = %i, bin %i = %f",i,nusefulbins,j,binsProfile[j])); 
697   }
698
699   TProfile* hSlewingProf = new TProfile("hSlewingProf", "hSlewingProf",nusefulbins, binsProfile, "G");  // CHECK THE BUILD OPTION, PLEASE!!!!!!
700   TH2F * htimetot = new TH2F("htimetot","htimetot",nbinToT, minToT, maxToT,600,-5,10);
701   for (Int_t irun=0;irun<fNruns;irun++){
702     Int_t ientry = -1;
703     ientry = i+irun*fNChannels;
704     fTree->GetEntry(ientry);
705     ntracksRun=nentries/3;
706     for (Int_t j=0;j<ntracksRun;j++){
707       Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
708       Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
709       Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
710       Float_t tot = p[idxexToT];
711       Float_t time = p[idxexTime]-p[idxexExTime];
712       AliDebug (2,Form("track = %i, time = %f, tot = %f, time-meantime = %f",j,time, tot, time-meantime));
713       hSlewingProf->Fill(tot,time);  // if meantime is not used, the fill may be moved in the loop above
714       htimetot->Fill(tot,time-meantime);  // if meantime is not used, the fill may be moved in the loop above
715     }
716   }
717
718   hSlewingProf->Fit("pol5",optionFit,"",0,4);
719   TF1 * calibfunc = (TF1*)hSlewingProf->GetFunction("pol5");
720   Float_t par[6];    
721   for(Int_t kk=0;kk<6;kk++){
722     par[kk]=calibfunc->GetParameter(kk);
723     AliDebug(2,Form("parameter %i = %f",kk,par[kk]));
724   }
725
726
727   if(strstr(optionSave,"save")){
728     TFile * fileProf = new TFile("TOFCalibSave.root","recreate");
729     fileProf->cd();   
730     TString profName=Form("Profile%06i",i);
731     TString timeTotName=Form("TimeTot%06i",i);
732     TString totName=Form("Tot%06i",i);
733     TString deltaName=Form("Delta%06i",i);
734     hSlewingProf->Write(profName);
735     htimetot->Write(timeTotName);
736     hToT->Write(totName);
737     hdeltaTime->Write(deltaName);
738     fileProf->Close();
739     delete fileProf;
740     fileProf=0x0;
741   }
742
743   delete hToT;
744   hToT=0x0; 
745   delete hSlewingProf;
746   hSlewingProf=0x0;
747   delete htimetot;
748   htimetot=0x0;
749   delete hdeltaTime;
750   hdeltaTime=0x0;
751
752   AliTOFChannelOffline * calChannel = (AliTOFChannelOffline*)fTOFCalOffline->At(i);
753   calChannel->SetSlewPar(par);
754   WriteParOfflineOnCDB("TOF/Calib","valid");
755   return 0;
756 }
757 //----------------------------------------------------------------------------
758 Int_t AliTOFcalib::Calibrate(Int_t nch, Int_t *ch, Option_t *optionSave, Option_t *optionFit){
759
760   // calibrating an array of channels
761   // computing calibration parameters
762   // Returning codes:
763   // 0 -> everything was ok
764   // 1 -> no tree for calibration found
765   // 2 -> not enough statistics to perform calibration
766   // 3 -> problems with arrays
767   
768   TH1::AddDirectory(0);
769
770   AliInfo(Form("*** Calibrating Histograms %s, number of channels = %i", GetName(),nch)) ; 
771   AliInfo(Form("Option for Saving histos = %s",optionSave )) ; 
772   AliInfo(Form("Option for Fitting Profile histos = %s",optionFit )) ; 
773   for (Int_t ich=0; ich<nch; ich++){
774     Int_t i = ch[ich];
775     AliInfo(Form("Calibrating channel = %i",i )) ; 
776   }
777   Float_t p[MAXCHENTRIESSMALL];
778   Int_t nentries;
779   fTree->SetBranchAddress("nentries",&nentries);
780   fTree->SetBranchAddress("TOFentries",p);
781
782   Float_t ntracksTotalmean =0;
783   for (Int_t ich=0; ich<nch; ich++){
784     Int_t ientry = -1;
785       Int_t i = ch[ich];
786       for (Int_t irun=0;irun<fNruns;irun++){
787       ientry = i+irun*fNChannels;
788       fTree->GetEntry(ientry);
789       ntracksTotalmean+=nentries/3;
790     }
791   }
792
793   ntracksTotalmean/=nch;
794   if (ntracksTotalmean < MEANENTRIES) { 
795     AliInfo(Form(" Too small mean number of entires per channel (mean number = %f) not calibrating and exiting.....",ntracksTotalmean));
796     return 2;
797   }
798
799   //filling ToT and Time arrays
800
801   Int_t nbinToT = 100;  // ToT bin width in Profile = 48.8 ps 
802   Float_t minToT = 0;   // ns
803   Float_t maxToT = 4.88;  // ns
804   TFile * fileProf=0x0;
805   if(strstr(optionSave,"save")){
806     fileProf = new TFile("TOFCalibSave.root","recreate");
807   }
808   for (Int_t ich=0; ich<nch; ich++) {
809     TH1F *hToT = new TH1F("htot","htot",nbinToT, minToT, maxToT);
810     TH1F *hdeltaTime = new TH1F("hdeltaTime","hdeltaTime",200,2,4);
811     Double_t binsProfile[101]; // sized larger than necessary, the correct 
812     // dim being set in the booking of the profile
813     TH2F * htimetot = new TH2F("htimetot","htimetot",nbinToT, minToT, maxToT,600,-5,10);
814     Int_t ntracksTotal = 0;
815     Int_t ntracksRun = 0;
816     Int_t nusefulbins=0;
817     Float_t meantime=0;
818     Int_t i=-1;
819     for (Int_t irun=0;irun<fNruns;irun++){
820       i = ch[ich]+irun*fNChannels;
821       AliDebug(2,Form("Calibrating channel %i",i));
822       fTree->GetEntry(i);
823       ntracksTotal+=nentries/3;
824     }
825     if (ntracksTotal < MEANENTRIES) {
826       AliInfo(Form(" Too small mean number of entires in channel %i (number of tracks = %f), not calibrating channel and continuing.....",i,ntracksTotal));
827       continue;
828     }
829   
830     for (Int_t irun=0;irun<fNruns;irun++){
831       Int_t i = ch[ich]+irun*fNChannels;
832       fTree->GetEntry(i);
833       ntracksRun=nentries/3;
834       for (Int_t j=0;j<ntracksRun;j++){
835         Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
836         Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
837         Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
838         Float_t tot = p[idxexToT];
839         hdeltaTime->Fill(p[idxexTime]-p[idxexExTime]);
840         meantime+=p[idxexTime]-p[idxexExTime];
841         hToT->Fill(tot);
842       }
843     }
844
845     nusefulbins = FindBins(hToT,&binsProfile[0]);
846     meantime/=ntracksTotal;
847     for (Int_t j=1;j<=nusefulbins;j++) {
848       AliDebug(2,Form(" channel %i, nusefulbins = %i, bin %i = %f",i,nusefulbins,j,binsProfile[j])); 
849     }
850
851     TProfile* hSlewingProf = new TProfile("hSlewingProf", "hSlewingProf",nusefulbins, binsProfile, "G");  // CHECK THE BUILD OPTION, PLEASE!!!!!!
852     for (Int_t irun=0;irun<fNruns;irun++){
853       Int_t i = ch[ich]+irun*fNChannels;
854       fTree->GetEntry(i);
855       ntracksRun=nentries/3;
856       for (Int_t j=0;j<ntracksRun;j++){
857         Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
858         Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
859         Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
860         Float_t tot = p[idxexToT];
861         Float_t time = p[idxexTime]-p[idxexExTime];
862         AliDebug(2,Form("track = %i, time = %f, tot = %f, time-meantime = %f",j,time, tot, time-meantime));
863         hSlewingProf->Fill(tot,time);  // if meantime is not used, the fill may be moved in the loop above
864         htimetot->Fill(tot,time-meantime);  // if meantime is not used, the fill may be moved in the loop above
865       }
866     }
867     
868     hSlewingProf->Fit("pol5",optionFit,"",1,4);
869     TF1 * calibfunc = (TF1*)hSlewingProf->GetFunction("pol5");
870     Float_t par[6];    
871     for(Int_t kk=0;kk<6;kk++){
872       par[kk]=calibfunc->GetParameter(kk);
873       AliDebug(2,Form("parameter %i = %f",kk,par[kk]));
874     }
875     
876     if(strstr(optionSave,"save") && fileProf){
877       TString profName=Form("Profile%06i",i);
878       TString timeTotName=Form("TimeTot%06i",i);
879       TString totName=Form("Tot%06i",i);
880       TString deltaName=Form("Delta%06i",i);
881       fileProf->cd();
882       hSlewingProf->Write(profName);
883       htimetot->Write(timeTotName);
884       hToT->Write(totName);
885       hdeltaTime->Write(deltaName);
886     }
887
888     AliTOFChannelOffline * calChannel = (AliTOFChannelOffline*)fTOFCalOffline->At(i);
889     calChannel->SetSlewPar(par);
890     delete hToT;
891     hToT=0x0;
892     delete hSlewingProf;
893     hSlewingProf=0x0;
894     delete htimetot;
895     htimetot=0x0;
896     delete hdeltaTime;
897     hdeltaTime=0x0;
898   }
899
900   if(strstr(optionSave,"save") && fileProf){
901     fileProf->Close();
902     delete fileProf;
903     fileProf=0x0;
904   }
905   WriteParOfflineOnCDB("TOF/Calib","valid");
906
907   return 0;
908 }
909 //----------------------------------------------------------------------------
910 Int_t AliTOFcalib::CalibrateFromProfile(Int_t ich, Option_t *optionSave, Option_t *optionFit){
911
912   // computing calibration parameters using the old profiling algo
913   // Returning codes:
914   // 0 -> everything was ok
915   // 1 -> no tree for calibration found
916   // 2 -> not enough statistics to perform calibration
917   // 3 -> problems with arrays
918
919   TH1::AddDirectory(0);
920
921   AliInfo(Form("*** Calibrating Histograms From Profile %s", GetName())) ; 
922   AliInfo(Form("Option for Saving histos = %s",optionSave )) ; 
923   AliInfo(Form("Option for Fitting Profile histos = %s",optionFit )) ; 
924   Float_t p[MAXCHENTRIESSMALL];
925   Int_t nentries;
926   Int_t ntracksTotal=0;
927   fTree->SetBranchAddress("nentries",&nentries);
928   fTree->SetBranchAddress("TOFentries",p);
929
930   for (Int_t irun=0;irun<fNruns;irun++){
931     Int_t i = ich+irun*fNChannels;
932     fTree->GetEntry(i);
933     ntracksTotal+=nentries/3;
934   }
935
936   if (ntracksTotal < MEANENTRIES) {  
937     AliInfo(Form(" Too small mean number of entires per channel (mean number = %f) not calibrating and exiting.....",ntracksTotal));
938     return 2;
939   }
940
941   TH1F * hProf = new TH1F();
942   hProf = Profile(ich);
943   hProf->Fit("pol5",optionFit,"",0,4);
944   TF1 * calibfunc = (TF1*)hProf->GetFunction("pol5");
945   Float_t par[6];    
946   for(Int_t kk=0;kk<6;kk++){
947     par[kk]=calibfunc->GetParameter(kk);
948     AliDebug(2,Form("parameter %i = %f",kk,par[kk]));
949   }
950
951   if(strstr(optionSave,"save")){
952     TFile * fileProf = new TFile("TOFCalibSave.root","recreate");
953     fileProf->cd(); 
954     TString profName=Form("Profile%06i",ich);
955     hProf->Write(profName);
956     fileProf->Close();
957     delete fileProf;
958     fileProf=0x0;
959   }
960
961   delete hProf;
962   hProf=0x0;
963   AliTOFChannelOffline * calChannel = (AliTOFChannelOffline*)fTOFCalOffline->At(ich);
964   calChannel->SetSlewPar(par);
965   WriteParOfflineOnCDB("TOF/Calib","valid");
966   return 0;
967 }
968 //----------------------------------------------------------------------------
969 Int_t AliTOFcalib::Calibrate(Option_t *optionSave, Option_t *optionFit){
970
971   // calibrating the whole TOF
972   // computing calibration parameters
973   // Returning codes:
974   // 0 -> everything was ok
975   // 1 -> no tree for calibration found
976   // 2 -> not enough statistics to perform calibration
977   // 3 -> problems with arrays
978
979   TH1::AddDirectory(0);
980
981   AliInfo(Form("*** Calibrating Histograms %s, all channels", GetName())) ; 
982   AliInfo(Form("Option for Saving histos = %s",optionSave )) ; 
983   AliInfo(Form("Option for Fitting Profile histos = %s",optionFit )) ; 
984
985   TFile * fileProf=0x0;
986   if(strstr(optionSave,"save")){
987     fileProf = new TFile("TOFCalibSave.root","recreate");
988   }
989
990   Float_t p[MAXCHENTRIESSMALL];
991   Int_t nentries;
992   fTree->SetBranchAddress("nentries",&nentries);
993   fTree->SetBranchAddress("TOFentries",p);
994
995   Float_t ntracksTotalmean =0;
996   for (Int_t ii=0; ii<fNChannels; ii++){
997     for (Int_t irun=0;irun<fNruns;irun++){
998       Int_t i = ii+irun*fNChannels;
999       fTree->GetEntry(i);
1000       ntracksTotalmean+=nentries/3;
1001     }
1002   }
1003
1004   ntracksTotalmean/=fNChannels;
1005   if (ntracksTotalmean < MEANENTRIES) {
1006     AliInfo(Form(" Too small mean number of entires per channel (mean number = %f) not calibrating and exiting.....",ntracksTotalmean));
1007     return 2;
1008   }
1009
1010   //filling ToT and Time arrays
1011
1012   Int_t nbinToT = 100;  // ToT bin width in Profile = 50.0 ps 
1013   Float_t minToT = 0;   // ns
1014   Float_t maxToT = 4.88;// ns
1015   for (Int_t ii=0; ii<fNChannels; ii++) {
1016     TH1F *hToT = new TH1F("htot","htot",nbinToT, minToT, maxToT);
1017     TH1F *hdeltaTime = new TH1F("hdeltaTime","hdeltaTime",200,2,4);
1018     TH2F * htimetot = new TH2F("htimetot","htimetot",nbinToT, minToT, maxToT,600,-5,10);
1019     if (ii%1000 == 0) AliDebug(1,Form("Calibrating channel %i ",ii));
1020     //Int_t i = 3;
1021     Int_t nusefulbins=0;
1022     Double_t binsProfile[101]; // sized larger than necessary, the correct 
1023                               // dim being set in the booking of the profile
1024     Int_t ntracksRun = 0;
1025     Int_t ntracksTotal = 0;
1026     for (Int_t irun=0;irun<fNruns;irun++){
1027       Int_t i = ii+irun*fNChannels;
1028       fTree->GetEntry(i);
1029       ntracksTotal+=nentries/3;
1030     }
1031     if (ntracksTotal < MEANENTRIES) {
1032       AliInfo(Form(" Too small mean number of entires in channel %i (number of tracks = %f), not calibrating channel and continuing.....",ii,ntracksTotal));
1033       continue;
1034     }
1035     Float_t meantime=0;
1036     for (Int_t irun=0;irun<fNruns;irun++){
1037       Int_t i = ii+irun*fNChannels;
1038       fTree->GetEntry(i);
1039       ntracksRun=nentries/3;
1040       for (Int_t j=0;j<ntracksRun;j++){
1041         Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
1042         Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
1043         Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
1044         Float_t tot = p[idxexToT];
1045         hdeltaTime->Fill(p[idxexTime]-p[idxexExTime]);
1046         meantime+=p[idxexTime]-p[idxexExTime];
1047         hToT->Fill(tot);
1048       }
1049     }
1050     nusefulbins = FindBins(hToT,&binsProfile[0]);
1051     meantime/=ntracksTotal;
1052     for (Int_t j=0;j<nusefulbins;j++) {
1053       AliDebug(2,Form(" channel %i, usefulbin = %i, low edge = %f",ii,j,binsProfile[j])); 
1054     }
1055     TProfile* hSlewingProf = new TProfile("hSlewingProf", "hSlewingProf",nusefulbins, binsProfile, "G");  // CHECK THE BUILD OPTION, PLEASE!!!!!!
1056     for (Int_t irun=0;irun<fNruns;irun++){
1057       Int_t i = ii+irun*fNChannels;
1058       fTree->GetEntry(i);
1059       ntracksRun=nentries/3;
1060       for (Int_t j=0;j<ntracksRun;j++){
1061         Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
1062         Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
1063         Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
1064         Float_t tot = p[idxexToT];
1065         Float_t time = p[idxexTime]-p[idxexExTime];
1066         AliDebug (2,Form("track = %i, time = %f, tot = %f, time-meantime = %f",j,time, tot, time-meantime));
1067         hSlewingProf->Fill(tot,time);  // if meantime is not used, the fill may be moved in the loop above
1068         htimetot->Fill(tot,time-meantime);  // if meantime is not used, the fill may be moved in the loop above
1069       }
1070     }
1071     hSlewingProf->Fit("pol5",optionFit,"",1,4);
1072     TF1 * calibfunc = (TF1*)hSlewingProf->GetFunction("pol5");
1073     Float_t par[6];    
1074     for(Int_t kk=0;kk<6;kk++){
1075       par[kk]=calibfunc->GetParameter(kk);
1076       AliDebug(2,Form("parameter %i = %f",kk,par[kk]));
1077     }
1078
1079     if(strstr(optionSave,"save") && fileProf){
1080       TString profName=Form("Profile%06i",ii);
1081       TString timeTotName=Form("TimeTot%06i",ii);
1082       TString totName=Form("Tot%06i",ii);
1083       TString deltaName=Form("Delta%06i",ii);
1084       fileProf->cd();
1085       hSlewingProf->Write(profName);
1086       htimetot->Write(timeTotName);
1087       hToT->Write(totName);
1088       hdeltaTime->Write(deltaName);
1089     }
1090     AliTOFChannelOffline * calChannel = (AliTOFChannelOffline*)fTOFCalOffline->At(ii);
1091     calChannel->SetSlewPar(par);
1092
1093     delete hToT;
1094     hToT=0x0;
1095     delete hSlewingProf;
1096     hSlewingProf=0x0;
1097     delete htimetot;
1098     htimetot=0x0;
1099     delete hdeltaTime;
1100     hdeltaTime=0x0;
1101   }
1102
1103   if(strstr(optionSave,"save")){
1104     fileProf->Close();
1105     delete fileProf;
1106     fileProf=0x0;
1107   }
1108   WriteParOfflineOnCDB("TOF/Calib","valid");
1109   return 0;
1110 }
1111
1112 //-----------------------------------------------------------------------
1113 TH1F* AliTOFcalib::Profile(Int_t ich)
1114 {
1115   // profiling algo
1116
1117   Float_t p[MAXCHENTRIESSMALL];
1118   Int_t nentries;
1119   fTree->SetBranchAddress("nentries",&nentries);
1120   fTree->SetBranchAddress("TOFentries",p);
1121
1122   //Prepare histograms for Slewing Correction
1123   const Int_t knbinToT = 100;
1124   Int_t nbinTime = 200;
1125   Float_t minTime = -5.5; //ns
1126   Float_t maxTime = 5.5; //ns
1127   Float_t minToT = 0; //ns
1128   Float_t maxToT = 5.; //ns
1129   Float_t deltaToT = (maxToT-minToT)/knbinToT;
1130   Double_t mTime[knbinToT+1],mToT[knbinToT+1],meanTime[knbinToT+1], meanTime2[knbinToT+1],vToT[knbinToT+1], vToT2[knbinToT+1],meanToT[knbinToT+1],meanToT2[knbinToT+1],vTime[knbinToT+1],vTime2[knbinToT+1],xlow[knbinToT+1],sigmaTime[knbinToT+1];
1131   Int_t n[knbinToT+1], nentrx[knbinToT+1];
1132   Double_t sigmaToT[knbinToT+1];
1133   for (Int_t i = 0; i < knbinToT+1 ; i++){
1134     mTime[i]=0;
1135     mToT[i]=0;
1136     n[i]=0;
1137     meanTime[i]=0;
1138     meanTime2[i]=0;
1139     vToT[i]=0;
1140     vToT2[i]=0;
1141     meanToT[i]=0;
1142     meanToT2[i]=0;
1143     vTime[i]=0;
1144     vTime2[i]=0;
1145     xlow[i]=0;
1146     sigmaTime[i]=0;
1147     sigmaToT[i]=0;
1148     n[i]=0;
1149     nentrx[i]=0;
1150   }
1151   TH2F* hSlewing = new TH2F("hSlewing", "hSlewing", knbinToT, minToT, maxToT, nbinTime, minTime, maxTime);
1152   Int_t ntracksRun = 0;
1153   TH1F *histo = new TH1F("histo", "1D Time vs ToT", knbinToT, minToT, maxToT);
1154   for (Int_t irun=0;irun<fNruns;irun++){
1155     Int_t i = ich+irun*fNChannels;
1156     fTree->GetEntry(i);
1157     ntracksRun=nentries/3;
1158     for (Int_t j=0;j<ntracksRun;j++){
1159       Int_t idxexToT = (j* NIDXSMALL)+DELTAIDXTOT; 
1160       Int_t idxexTime = (j* NIDXSMALL)+DELTAIDXTIME; 
1161       Int_t idxexExTime = (j* NIDXSMALL)+DELTAIDXPID; 
1162       Float_t tot = p[idxexToT];
1163       Float_t time = p[idxexTime]-p[idxexExTime];
1164       Int_t nx = (Int_t)((tot-minToT)/deltaToT)+1;
1165       if ((tot != 0) && ( time!= 0)){
1166         vTime[nx]+=time;
1167         vTime2[nx]+=time*time;
1168         vToT[nx]+=tot;
1169         vToT2[nx]+=tot*tot;
1170         nentrx[nx]++;
1171         hSlewing->Fill(tot,time);
1172       }
1173     }
1174   }
1175   Int_t nbinsToT=hSlewing->GetNbinsX();
1176   if (nbinsToT != knbinToT) {
1177     AliError("Profile :: incompatible numbers of bins");
1178     return 0x0;
1179   }
1180   
1181   Int_t usefulBins=0;
1182   for (Int_t i=1;i<=nbinsToT;i++){
1183     if (nentrx[i]!=0){
1184       n[usefulBins]+=nentrx[i];
1185       if (n[usefulBins]==0 && i == nbinsToT) {
1186         break;
1187       }
1188       meanTime[usefulBins]+=vTime[i];
1189       meanTime2[usefulBins]+=vTime2[i];
1190       meanToT[usefulBins]+=vToT[i];
1191       meanToT2[usefulBins]+=vToT2[i];
1192       if (n[usefulBins]<10 && i!=nbinsToT) continue; 
1193       mTime[usefulBins]=meanTime[usefulBins]/n[usefulBins];
1194       mToT[usefulBins]=meanToT[usefulBins]/n[usefulBins];
1195       sigmaTime[usefulBins]=TMath::Sqrt(1./n[usefulBins]/n[usefulBins]
1196                             *(meanTime2[usefulBins]-meanTime[usefulBins]
1197                             *meanTime[usefulBins]/n[usefulBins]));
1198       if ((1./n[usefulBins]/n[usefulBins]
1199            *(meanToT2[usefulBins]-meanToT[usefulBins]
1200              *meanToT[usefulBins]/n[usefulBins]))< 0) {
1201         AliError(" too small radical" );
1202         sigmaToT[usefulBins]=0;
1203       }
1204       else{       
1205         sigmaToT[usefulBins]=TMath::Sqrt(1./n[usefulBins]/n[usefulBins]
1206                              *(meanToT2[usefulBins]-meanToT[usefulBins]
1207                              *meanToT[usefulBins]/n[usefulBins]));
1208       }
1209       usefulBins++;
1210     }
1211   }
1212   for (Int_t i=0;i<usefulBins;i++){
1213     Int_t binN = (Int_t)((mToT[i]-minToT)/deltaToT)+1;
1214     histo->Fill(mToT[i],mTime[i]);
1215     histo->SetBinError(binN,sigmaTime[i]);
1216   } 
1217   delete hSlewing;
1218   hSlewing=0x0;
1219
1220   return histo;
1221 }
1222 //----------------------------------------------------------------------------
1223 Int_t AliTOFcalib::FindBins(TH1F* h, Double_t *binsProfile) const{
1224
1225   // to determine the bins for ToT histo
1226
1227   Int_t cont = 0;
1228   Int_t startBin = 1;
1229   Int_t nbin = h->GetNbinsX();
1230   Int_t nentries = (Int_t)h->GetEntries();
1231   Float_t max = h->GetBinLowEdge(nbin);
1232   Int_t nusefulbins=0;
1233   Int_t maxcont=0;
1234   // setting maxvalue of entries per bin
1235   if (nentries <= 60) maxcont = 2;
1236   else  if (nentries <= 100) maxcont = 5;
1237   else  if (nentries <= 500) maxcont = 10;
1238   else  maxcont = 20;
1239   for (Int_t j=1;j<=nbin;j++) {
1240     cont += (Int_t)h->GetBinContent(j);
1241     if (j<nbin){
1242       if (cont>=maxcont){
1243         nusefulbins++;
1244         binsProfile[nusefulbins-1]=h->GetBinLowEdge(startBin);
1245         cont=0;
1246         startBin=j+1;
1247         continue;
1248       }
1249     }
1250     else{
1251       if (cont>=maxcont){
1252         nusefulbins++;
1253         binsProfile[nusefulbins-1]=h->GetBinLowEdge(startBin);
1254         binsProfile[nusefulbins]=max;
1255       }
1256       else {
1257         binsProfile[nusefulbins]=h->GetBinLowEdge(startBin);
1258       }
1259     }
1260   }
1261   return nusefulbins;
1262 }