1cd3cd7d114ff302b65362a1a800540bded11bd4
[u/mrichter/AliRoot.git] / TPC / Calib / AliTPCPreprocessorOffline.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
18 /*
19   Responsible: marian.ivanov@cern.ch 
20   Code to analyze the TPC calibration and to produce OCDB entries  
21
22
23    .x ~/rootlogon.C
24    gSystem->Load("libANALYSIS");
25    gSystem->Load("libTPCcalib");
26
27    AliTPCPreprocessorOffline proces;
28    TString ocdbPath="local:////"
29    ocdbPath+=gSystem->GetFromPipe("pwd");
30
31    proces.CalibTimeGain("CalibObjects.root",run0,run1,ocdbPath);
32    proces.CalibTimeVdrift("CalibObjects.root",run0,run1,ocdbPath);
33   // take the raw calibration data from the file CalibObjects.root 
34   // and make a OCDB entry with run  validity run0-run1
35   // results are stored at the ocdbPath - local or alien ...
36   // default storage ""- data stored at current working directory 
37  
38   e.g.
39   gSystem->Load("libANALYSIS");
40   gSystem->Load("libTPCcalib");
41   AliTPCPreprocessorOffline proces;
42   proces.CalibTimeGain("TPCMultObjects.root",114000,140040,0);
43   TFile oo("OCDB/TPC/Calib/TimeGain/Run114000_121040_v0_s0.root")
44  TObjArray * arr = AliCDBEntry->GetObject()
45   arr->At(4)->Draw("alp")
46
47 */
48 #include "Riostream.h"
49 #include <fstream>
50 #include "TMap.h"
51 #include "TGraphErrors.h"
52 #include "AliExternalTrackParam.h"
53 #include "TFile.h"
54 #include "TDirectory.h"
55 #include "TGraph.h"
56 #include "TMultiGraph.h"
57 #include "TCanvas.h"
58 #include "THnSparse.h"
59 #include "TLegend.h"
60 #include "TPad.h"
61 #include "TH2D.h"
62 #include "TH3D.h"
63 #include "AliTPCROC.h"
64 #include "AliTPCCalROC.h"
65 #include "AliESDfriend.h"
66 #include "AliTPCcalibTime.h"
67 #include "AliSplineFit.h"
68 #include "AliCDBMetaData.h"
69 #include "AliCDBId.h"
70 #include "AliCDBManager.h"
71 #include "AliCDBStorage.h"
72 #include "AliTPCcalibBase.h"
73 #include "AliTPCcalibDB.h"
74 #include "AliTPCcalibDButil.h"
75 #include "AliRelAlignerKalman.h"
76 #include "AliTPCParamSR.h"
77 #include "AliTPCcalibTimeGain.h"
78 #include "AliTPCcalibGainMult.h"
79 #include "AliTPCcalibAlign.h"
80 #include "AliSplineFit.h"
81 #include "AliTPCComposedCorrection.h"
82 #include "AliTPCExBTwist.h"
83 #include "AliTPCCalibGlobalMisalignment.h"
84 #include "TStatToolkit.h"
85 #include "TChain.h"
86 #include "TCut.h"
87 #include "AliTrackerBase.h"
88 #include "AliTracker.h"
89 #include "AliTPCPreprocessorOffline.h"
90 #include "AliTPCCorrectionFit.h"
91
92 #include "AliTPCClusterParam.h"
93 #include "AliTPCRecoParam.h"
94
95 using std::endl;
96 using std::cout;
97
98 ClassImp(AliTPCPreprocessorOffline)
99
100 AliTPCPreprocessorOffline::AliTPCPreprocessorOffline():
101   TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
102   fMinEntries(500),                      // minimal number of entries for fit
103   fStartRun(0),                         // start Run - used to make fast selection in THnSparse
104   fEndRun(0),                           // end   Run - used to make fast selection in THnSparse
105   fStartTime(0),                        // fStartTime - used to make fast selection in THnSparse
106   fEndTime(0),                          // fEndTime   - used to make fast selection in THnSparse
107   fOCDBstorage(0),                       // OCDB storage
108   fVdriftArray(new TObjArray),
109   fTimeDrift(0),
110   fGraphMIP(0),                // graph time dependence of MIP
111   fGraphCosmic(0),             // graph time dependence at Plateu
112   fGraphAttachmentMIP(0),
113   fFitMIP(0),                  // fit of dependence - MIP
114   fFitCosmic(0),               // fit of dependence - Plateu
115   fGainArray(new TObjArray),               // array to be stored in the OCDB
116   fGainMIP(0),          // calibration component for MIP
117   fGainCosmic(0),       // calibration component for cosmic
118   fGainMult(0),
119   fAlignTree(0),        // alignment tree
120   fSwitchOnValidation(kFALSE), // flag to switch on validation of OCDB parameters
121   fMinGain(2.0),
122   fMaxGain(3.0),
123   fMaxVdriftCorr(0.03),
124   fNtracksVdrift(0),
125   fMinTracksVdrift(0),
126   fNeventsVdrift(0),
127   fMinEventsVdrift(0),
128   fCalibrationStatus(0)
129 {
130   //
131   // default constructor
132   //
133 }
134
135 AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
136   //
137   // Destructor
138   //
139 }
140
141
142
143
144 void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const  timeDrift){
145   //
146   // find the fist and last run
147   //
148   TObjArray *hisArray =timeDrift->GetHistoDrift();
149   {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
150     THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
151     if (!addHist) continue;
152     if (addHist->GetEntries()<fMinEntries) continue;
153     TH1D* histo    =addHist->Projection(3);
154     TH1D* histoTime=addHist->Projection(0);
155     printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
156
157     if (fStartRun<=0){ 
158       fStartRun=histo->FindFirstBinAbove(0);
159       fEndRun  =histo->FindLastBinAbove(0);
160     }else{
161       fStartRun=TMath::Min(histo->FindFirstBinAbove(0),fStartRun);
162       fEndRun  =TMath::Max(histo->FindLastBinAbove(0),fEndRun);
163     }
164     if (fStartTime==0){ 
165       fStartTime=histoTime->FindFirstBinAbove(0);
166       fEndTime  =histoTime->FindLastBinAbove(0);
167     }else{
168       fStartTime=TMath::Min(histoTime->FindFirstBinAbove(0),fStartTime);
169       fEndTime  =TMath::Max(histoTime->FindLastBinAbove(0),fEndTime);
170     }
171     delete histo;
172     delete histoTime;
173   }}
174   if (fStartRun<0) fStartRun=0;
175   if (fEndRun<0) fEndRun=100000000;
176   printf("Run range  :\t%d-%d\n", fStartRun, fEndRun);
177   printf("Time range :\t%d-%d\n", fStartTime, fEndTime);
178
179 }
180
181
182
183 void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, AliCDBStorage* pocdbStorage){
184   //
185   // make calibration of the drift velocity
186   // Input parameters:
187   //      file                   - the location of input file
188   //      ustartRun, uendrun     - run validity period 
189   //      pocdbStorage           - path to hte OCDB storage
190   //                             - if empty - local storage 'pwd' uesed
191   if (pocdbStorage) fOCDBstorage=pocdbStorage;
192   else {
193     TString localStorage = "local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; 
194     fOCDBstorage=AliCDBManager::Instance()->GetStorage(localStorage.Data());
195   }
196
197   //
198   // 1. Initialization and run range setting
199   TFile fcalib(file);
200   TObject* obj = dynamic_cast<TObject*>(fcalib.Get("TPCCalib"));
201   TObjArray* array = dynamic_cast<TObjArray*>(obj);
202   TDirectory* dir = dynamic_cast<TDirectory*>(obj);
203   if (dir) {
204     fTimeDrift = dynamic_cast<AliTPCcalibTime*>(dir->Get("calibTime"));
205   }
206   else if (array){
207     fTimeDrift = (AliTPCcalibTime *)array->FindObject("calibTime");
208   } else {
209     fTimeDrift = (AliTPCcalibTime*)fcalib.Get("calibTime");
210   }
211   if(!fTimeDrift) return;
212
213   //extract statistics
214   fNtracksVdrift = TMath::Nint(fTimeDrift->GetResHistoTPCITS(0)->GetEntries());
215   //if we have 0 ITS TPC matches it means we have no ITS tracks and we try to use TPC-TOF matching for calibration
216   if (fNtracksVdrift==0) fNtracksVdrift=TMath::Nint(fTimeDrift->GetResHistoTPCTOF(0)->GetEntries());
217   fNeventsVdrift = TMath::Nint(fTimeDrift->GetTPCVertexHisto(0)->GetEntries());
218
219   fStartRun=ustartRun;
220   fEndRun=ustartRun; 
221   TObjArray *hisArray =fTimeDrift->GetHistoDrift();  
222   GetRunRange(fTimeDrift);
223   for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
224     THnSparse* addHist=(THnSparse*)hisArray->At(i);
225     if (!addHist) continue;
226     if (fStartTime<fEndTime) addHist->GetAxis(0)->SetRange(fStartTime-1,fEndTime+1);
227     if (fStartRun<fEndRun) addHist->GetAxis(3)->SetRange(fStartRun-1,fEndRun+1);
228   }
229   //
230   //
231   // 2. extraction of the information
232   //
233   fVdriftArray = new TObjArray();
234   AddAlignmentGraphs(fVdriftArray,fTimeDrift);
235   AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
236   AddLaserGraphs(fVdriftArray,fTimeDrift);
237   
238   //
239   // 3. Append QA plots
240   //
241   MakeDefaultPlots(fVdriftArray,fVdriftArray);
242
243   //
244   // 4. validate OCDB entries
245   //
246   if(fSwitchOnValidation==kTRUE && ValidateTimeDrift()==kFALSE) { 
247     Printf("TPC time drift OCDB parameters out of range!");
248     return;
249   }
250   //
251   //4.b make alignment
252   //
253   MakeFitTime();
254   TFile * ftime= TFile::Open("fitITSVertex.root");
255   if (ftime){
256     TObject * alignmentTime=ftime->Get("FitCorrectionTime");
257     if (alignmentTime) fVdriftArray->AddLast(alignmentTime);
258   }
259   //
260   // 5.) Add the RecoParam and ClusterParam - for compatibility checks -different sets of parameters can invalidate calibration 
261   //
262   AliTPCClusterParam *clParam =   AliTPCcalibDB::Instance()->GetClusterParam();
263   TObjArray *recoParams = new TObjArray(4) ;
264   for (Int_t i=0;i<4;i++) recoParams->AddAt(AliTPCcalibDB::Instance()->GetRecoParam(i),i);
265   fVdriftArray->AddLast(clParam);
266   fVdriftArray->AddLast(recoParams);
267   //
268   //
269   // 6. update of OCDB
270   //
271   //
272   UpdateOCDBDrift(ustartRun,uendRun,fOCDBstorage);
273 }
274
275 void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun,  AliCDBStorage* storage ){
276   //
277   // Update OCDB 
278   //
279   AliCDBMetaData *metaData= new AliCDBMetaData();
280   metaData->SetObjectClassName("TObjArray");
281   metaData->SetResponsible("Marian Ivanov");
282   metaData->SetBeamPeriod(1);
283   metaData->SetAliRootVersion("05-25-01"); //root version
284   metaData->SetComment("Calibration of the time dependence of the drift velocity");
285   AliCDBId* id1=NULL;
286   id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
287   storage->Put(fVdriftArray, (*id1), metaData);
288 }
289
290 Bool_t AliTPCPreprocessorOffline::ValidateTimeGain()
291 {
292   //
293   // Validate time gain corrections 
294   //
295   Printf("ValidateTimeGain..." );
296   Float_t minGain = fMinGain;
297   Float_t maxGain = fMaxGain;
298
299   TGraphErrors *gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
300   if (!gr) {
301     gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
302     if (!gr) 
303     { 
304       fCalibrationStatus |= kCalibFailedTimeGain;
305       return kFALSE;
306     }
307     Printf("Assuming given run is a cosmic run. Using gain calibration from Fermi-plateau muons.");
308   }
309   if(gr->GetN()<1) 
310   { 
311     fCalibrationStatus |= kCalibFailedTimeGain;
312     return kFALSE;
313   }
314
315   // check whether gain in the range
316   for(Int_t iPoint=0; iPoint<gr->GetN(); iPoint++) 
317   {
318     if(gr->GetY()[iPoint] < minGain || gr->GetY()[iPoint] > maxGain)  
319     { 
320       fCalibrationStatus |= kCalibFailedTimeGain;
321       return kFALSE;
322     }
323   }
324
325 return kTRUE;
326 }
327
328
329 Bool_t AliTPCPreprocessorOffline::ValidateTimeDrift()
330 {
331   //
332   // Validate time drift velocity corrections 
333   //
334   Printf("ValidateTimeDrift..." );
335
336   Float_t maxVDriftCorr = fMaxVdriftCorr;
337
338   TGraphErrors* gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
339   Printf("ALIGN_ITSB_TPC_DRIFTVD graph = %p",gr);
340   if (!gr)
341   {
342     gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_TOFB_TPC_DRIFTVD");
343     Printf("ALIGN_TOFB_TPC_DRIFTVD graph = %p",gr);
344   }
345
346   if(!gr) 
347   {
348     fCalibrationStatus|=kCalibFailedTimeDrift;
349     return kFALSE;
350   }
351   
352   // for now we validate even with low statistics
353   ////check if we have enough statistics
354   //if (fNtracksVdrift<fMinTracksVdrift) 
355   //{
356   //  fCalibrationStatus|=kCalibFailedTimeDrift;
357   //  return kFALSE;
358   //}
359
360   if(gr->GetN()<1)  { 
361     Printf("ALIGN_ITSB_TPC_DRIFTVD number of points = %d",gr->GetN());
362     {
363       fCalibrationStatus|=kCalibFailedTimeDrift;
364       return kFALSE;
365     }
366   }
367
368   // check whether drift velocity corrections in the range
369   for(Int_t iPoint = 0; iPoint<gr->GetN(); iPoint++) 
370   {
371     //Printf("Y value from the graph: %f",TMath::Abs(gr->GetY()[iPoint]));
372     if(TMath::Abs(gr->GetY()[iPoint]) > maxVDriftCorr)  
373     {
374       fCalibrationStatus|=kCalibFailedTimeDrift;
375       return kFALSE;
376     }
377   }
378
379 return kTRUE;
380 }
381
382 void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
383   //
384   //  update the OCDB entry for the nominal time0
385   //
386   //
387   //  AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
388   AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
389   TGraphErrors *grT =  (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
390   Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
391   Double_t deltaT   = deltaTcm/param->GetDriftV();
392   paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
393   paramNew->Update();
394
395   AliCDBMetaData *metaData= new AliCDBMetaData();
396   metaData->SetObjectClassName("TObjArray");
397   metaData->SetResponsible("Marian Ivanov");
398   metaData->SetBeamPeriod(1);
399   metaData->SetAliRootVersion("05-25-02"); //root version
400   metaData->SetComment("Updated calibration of nominal time 0");
401   AliCDBId* id1=NULL;
402   id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
403   fOCDBstorage->Put(param, (*id1), metaData);
404
405 }
406
407
408 void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
409   //
410   // Print the names of the entries in array
411   //
412   Int_t entries = array->GetEntries();
413   for (Int_t i=0; i<entries; i++){
414     if (!array->At(i)) continue;
415     printf("%d\t %s\n", i,  array->At(i)->GetName());
416   }
417 }
418
419
420
421 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
422   // 2 filters:
423   //    1. filter graph - error cut errSigmaCut
424   //    2. filter graph - medianCutAbs around median
425   //
426   // errSigmaCut   - cut on error
427   // medianCutAbs  - cut on value around median
428   Double_t dummy=0;               //   
429   //
430   // 1. filter graph - error cut errSigmaCut
431   //              
432   TGraphErrors *graphF; 
433   graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
434   delete graph;
435   if (!graphF) return 0;
436   graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
437   delete graphF;
438   if (!graph) return 0;
439   //
440   // filter graph - kMedianCutAbs around median
441   // 
442   graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
443   delete graph;
444   if (!graphF) return 0;
445   graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
446   delete graphF;
447   if (!graph) return 0;
448   return graph;
449 }
450
451
452
453 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
454   //
455   // filter outlyer measurement
456   // Only points around median +- cut filtered 
457   //
458   if (!graph) return  0;
459   Int_t kMinPoints=2;
460   Int_t npoints0 = graph->GetN();
461   Int_t npoints=0;
462   Float_t  rmsY=0;
463   Double_t *outx=new Double_t[npoints0];
464   Double_t *outy=new Double_t[npoints0];
465   Double_t *errx=new Double_t[npoints0];
466   Double_t *erry=new Double_t[npoints0];
467   //
468   //
469   if (npoints0<kMinPoints) {
470     delete []outx;
471     delete []outy;
472     delete []errx;
473     delete []erry;
474     return 0;
475   }
476   for (Int_t iter=0; iter<3; iter++){
477     npoints=0;
478     for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
479       if (graph->GetY()[ipoint]==0) continue;
480       if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;  
481       outx[npoints]  = graph->GetX()[ipoint];
482       outy[npoints]  = graph->GetY()[ipoint];
483       errx[npoints]  = graph->GetErrorX(ipoint);
484       erry[npoints]  = graph->GetErrorY(ipoint);
485       npoints++;
486     }
487     if (npoints<=1) break;
488     medianY  =TMath::Median(npoints,outy);
489     rmsY   =TMath::RMS(npoints,outy);
490   }
491   TGraphErrors *graphOut=0;
492   if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry); 
493   delete []outx;
494   delete []outy;
495   delete []errx;
496   delete []erry;
497   return graphOut;
498 }
499
500
501 void AliTPCPreprocessorOffline::AddHistoGraphs(  TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
502   //
503   // Add graphs corresponding to the alignment
504   //
505   const Double_t kErrSigmaCut=5;      // error sigma cut - for filtering
506   const Double_t kMedianCutAbs=0.03;  // error sigma cut - for filtering
507   //
508   TObjArray * array=timeDrift->GetHistoDrift();
509   if (array){
510     THnSparse* hist=NULL;
511     // 2.a) cosmics with different triggers
512     for (Int_t i=0; i<array->GetEntriesFast();i++){
513       hist=(THnSparseF*)array->UncheckedAt(i);
514       if(!hist) continue;
515       if (hist->GetEntries()<minEntries) continue;
516       //hist->Print();
517       TString name=hist->GetName();
518       Int_t dim[4]={0,1,2,3};
519       THnSparse* newHist=hist->Projection(4,dim);
520       newHist->SetName(name);
521       TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
522       delete newHist;
523       if (!graph) {
524         printf("Graph =%s filtered out\n", name.Data());
525         continue;
526       }
527       printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
528       Int_t pos=name.Index("_");
529       name=name(pos,name.Capacity()-pos);
530       TString graphName=graph->ClassName();
531       graphName+=name;
532       graphName.ToUpper();
533       //
534       graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
535       //
536       if (graph){
537         graph->SetMarkerStyle(i%8+20);
538         graph->SetMarkerColor(i%7);
539         graph->GetXaxis()->SetTitle("Time");
540         graph->GetYaxis()->SetTitle("v_{dcor}");
541         graph->SetName(graphName);
542         graph->SetTitle(graphName);
543         printf("Graph %d\t=\t%s\n", i, graphName.Data());
544         vdriftArray->Add(graph);
545       }
546     }
547   }
548 }
549
550
551
552
553 void AliTPCPreprocessorOffline::AddAlignmentGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
554   //
555   // Add graphs corresponding to alignment to the object array
556   //
557   TObjArray *arrayITS=0;
558   TObjArray *arrayTOF=0;
559   TObjArray *arrayTRD=0;
560   TMatrixD *mstatITS=0;
561   TMatrixD *mstatTOF=0;
562   TMatrixD *mstatTRD=0;
563   //
564   arrayITS=timeDrift->GetAlignITSTPC();
565   arrayTRD=timeDrift->GetAlignTRDTPC();
566   arrayTOF=timeDrift->GetAlignTOFTPC();
567
568   if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.7,50,fMaxVdriftCorr);
569   if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.7,1000,fMaxVdriftCorr);
570   if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.7,50,fMaxVdriftCorr);
571   //
572   TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
573   TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
574   TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
575   TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
576   TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
577   TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
578
579   TObjArray * arrayTRDP= 0x0;
580   TObjArray * arrayTRDM= 0x0;
581   TObjArray * arrayTRDB= 0x0;
582   arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
583   arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
584   arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
585   //
586   //
587   Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
588   TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
589                          arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
590                          arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
591   TString   grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
592                          "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
593                          "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
594   TString   grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
595                       "DELTAX", "DELTAY", "DELTAZ",
596                       "DRIFTVD", "T0", "VDGY"};
597
598   
599   TVectorD vX(entries);
600   TVectorD vY(entries);
601   TVectorD vEx(entries);
602   TVectorD vEy(entries);
603   TObjArray *arr=0;
604   for (Int_t iarray=0; iarray<12; iarray++){
605     arr = arrays[iarray];
606     if (arr==0) continue;
607     for (Int_t ipar=0; ipar<9; ipar++){      
608       Int_t counter=0;
609       for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
610         AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
611         if (!kalman) continue;
612         vX[counter]=kalman->GetTimeStamp();
613         vY[counter]=(*(kalman->GetState()))[ipar];
614         if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
615         vEx[counter]=0;
616         vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
617         counter++;
618       }
619     
620       TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
621                                           vY.GetMatrixArray(),
622                                           vEx.GetMatrixArray(),
623                                           vEy.GetMatrixArray());
624       TString grName=grnames[iarray];
625       grName+="_TPC_";
626       grName+=grpar[ipar];
627       graph->SetName(grName.Data());
628       vdriftArray->AddLast(graph);
629     }
630     delete arrays[iarray];
631   }  
632 }
633
634
635
636
637 void AliTPCPreprocessorOffline::AddLaserGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
638   //
639   // add graphs for laser
640   //
641   const Double_t delayL0L1 = 0.071;  //this is hack for 1/2 weeks
642   //THnSparse *hisN=0;
643   TGraphErrors *grLaser[6]={0,0,0,0,0,0};
644   //hisN = timeDrift->GetHistVdriftLaserA(0);
645   if (timeDrift->GetHistVdriftLaserA(0)){
646     grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
647     grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
648     vdriftArray->AddLast(grLaser[0]);
649   }    
650   if (timeDrift->GetHistVdriftLaserA(1)){
651     grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
652     grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
653     vdriftArray->AddLast(grLaser[1]);
654   }    
655   if (timeDrift->GetHistVdriftLaserA(2)){
656     grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
657     grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
658     vdriftArray->AddLast(grLaser[2]);
659   }    
660   if (timeDrift->GetHistVdriftLaserC(0)){
661     grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
662     grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
663     vdriftArray->AddLast(grLaser[3]);
664   }    
665   if (timeDrift->GetHistVdriftLaserC(1)){
666     grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
667     grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
668     vdriftArray->AddLast(grLaser[4]);
669   }    
670   if (timeDrift->GetHistVdriftLaserC(2)){
671     grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
672     grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");    
673     vdriftArray->AddLast(grLaser[5]);
674   }    
675   for (Int_t i=0; i<6;i++){
676     if (grLaser[i]) {
677       SetDefaultGraphDrift(grLaser[i], 1,(i+20));
678       grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
679     }
680   }
681 }
682  
683  
684 TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
685   //
686   // Make graph with mean values and rms
687   //
688   hisN->GetAxis(itime)->SetRange(0,100000000);
689   hisN->GetAxis(ival)->SetRange(0,100000000);
690   TH1 * hisT      = hisN->Projection(itime);
691   TH1 * hisV      = hisN->Projection(ival);
692   //
693   Int_t firstBinA = hisT->FindFirstBinAbove(2);
694   Int_t lastBinA  = hisT->FindLastBinAbove(2);    
695   Int_t firstBinV = hisV->FindFirstBinAbove(0);
696   Int_t lastBinV  = hisV->FindLastBinAbove(0);    
697   hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
698   hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
699   Int_t entries=0;
700   for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
701     Double_t cont = hisT->GetBinContent(ibin);
702     if (cont<minEntries) continue;
703     entries++;
704   }
705   TVectorD vecTime(entries);
706   TVectorD vecMean0(entries);
707   TVectorD vecRMS0(entries);
708   TVectorD vecMean1(entries);
709   TVectorD vecRMS1(entries);
710   entries=0;
711   for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
712       Double_t cont = hisT->GetBinContent(ibin);
713       if (cont<minEntries) continue;
714       //hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
715       Int_t minBin = ibin-1;
716       Int_t maxBin = ibin+1;
717       if(minBin <= 0) minBin = 1;
718       if(maxBin >= hisN->GetAxis(itime)->GetNbins()) maxBin = hisN->GetAxis(itime)->GetNbins()-1;
719       hisN->GetAxis(itime)->SetRange(minBin,maxBin);
720
721       Double_t time = hisT->GetBinCenter(ibin);
722       TH1 * his = hisN->Projection(ival);
723       Double_t nentries0= his->GetBinContent(his->FindBin(0));
724       if (cont-nentries0<minEntries) continue;
725       //
726       his->SetBinContent(his->FindBin(0),0);
727       vecTime[entries]=time;
728       vecMean0[entries]=his->GetMean()+offset;
729       vecMean1[entries]=his->GetMeanError();
730       vecRMS0[entries] =his->GetRMS();
731       vecRMS1[entries] =his->GetRMSError();
732       delete his;  
733       entries++;
734   }
735   delete hisT;
736   delete hisV;
737   TGraphErrors * graph =  new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(),                                    0, vecMean1.GetMatrixArray());
738
739   return graph;
740 }
741
742
743
744
745
746
747
748
749 void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
750   //
751   // Set default style for QA views
752   //
753   graph->GetXaxis()->SetTimeDisplay(kTRUE);
754   graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
755   graph->SetMaximum( 0.025);
756   graph->SetMinimum(-0.025);
757   graph->GetXaxis()->SetTitle("Time");
758   graph->GetYaxis()->SetTitle("v_{dcorr}");
759   //
760   graph->GetYaxis()->SetLabelSize(0.03);
761   graph->GetXaxis()->SetLabelSize(0.03);
762   //
763   graph->GetXaxis()->SetNdivisions(10,5,0);
764   graph->GetYaxis()->SetNdivisions(10,5,0);
765   //
766   graph->GetXaxis()->SetLabelOffset(0.02);
767   graph->GetYaxis()->SetLabelOffset(0.005);
768   //
769   graph->GetXaxis()->SetTitleOffset(1.3);
770   graph->GetYaxis()->SetTitleOffset(1.2);
771   //
772   graph->SetMarkerColor(color);
773   graph->SetLineColor(color);
774   graph->SetMarkerStyle(style);
775 }
776
777 void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
778   // 
779   // Set default pad style for QA
780   // 
781   pad->SetTicks(1,1);
782   pad->SetMargin(mx0,mx1,my0,my1);
783 }
784
785
786 void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray * /*picArray*/){
787   //
788   // 0. make a default QA plots
789   // 1. Store them in the array
790   //
791   //
792   Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
793   //
794   TGraphErrors* laserA       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
795   TGraphErrors* laserC       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
796   TGraphErrors* cosmic       =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
797   TGraphErrors* cross        =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
798   TGraphErrors* itstpcP       =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
799   TGraphErrors* itstpcM       =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
800   TGraphErrors* itstpcB       =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
801   //
802   if (laserA)  SetDefaultGraphDrift(laserA,2,25);
803   if (laserC)  SetDefaultGraphDrift(laserC,4,26);
804   if (cosmic)  SetDefaultGraphDrift(cosmic,3,27);
805   if (cross)   SetDefaultGraphDrift(cross,4,28);
806   if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
807   if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
808   if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
809   //
810   //
811   TPad *pad=0;
812   //
813   // Laser-Laser
814   //
815   if (laserA&&laserC){
816     pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
817     laserA->Draw("alp");
818     SetPadStyle(pad,mx0,mx1,my0,my1);
819     laserA->Draw("apl");
820     laserC->Draw("p");
821     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
822     legend->AddEntry(laserA,"Laser A side");
823     legend->AddEntry(laserC,"Laser C side");
824     legend->Draw();    
825     //picArray->AddLast(pad);
826   }
827
828   if (itstpcP&&itstpcM&&itstpcB){
829     pad = new TCanvas("ITSTPC","ITSTPC");
830     itstpcP->Draw("alp");
831     SetPadStyle(pad,mx0,mx1,my0,my1);    
832     itstpcP->Draw("alp");
833     gPad->Clear();
834     itstpcM->Draw("apl");
835     itstpcP->Draw("p");
836     itstpcB->Draw("p");
837     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
838     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
839     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
840     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
841     legend->Draw();    
842     //picArray->AddLast(pad);
843   }
844
845   if (itstpcB&&laserA&&itstpcP&&itstpcM){
846     pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
847     SetPadStyle(pad,mx0,mx1,my0,my1);    
848     laserA->Draw("alp");
849     itstpcP->Draw("p");
850     itstpcM->Draw("p");
851     itstpcB->Draw("p");
852     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
853     legend->AddEntry(laserA,"TPC laser");
854     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");   
855     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");   
856     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
857     legend->Draw();
858     //picArray->AddLast(pad);
859   }
860
861   if (itstpcP&&cross){ 
862     pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
863     SetPadStyle(pad,mx0,mx1,my0,my1);    
864     itstpcP->Draw("alp");
865     pad->Clear();
866     cross->Draw("ap");
867     itstpcP->Draw("p");
868     //
869     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
870
871     legend->AddEntry(cross,"TPC cross tracks");
872     legend->AddEntry(itstpcB,"ITS-TPC smooth");
873     legend->Draw();        
874     //picArray->AddLast(pad);
875   }
876   if (itstpcP&&cosmic){ 
877     pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
878     SetPadStyle(pad,mx0,mx1,my0,my1);    
879     itstpcP->Draw("alp");
880     pad->Clear();
881     cosmic->Draw("ap");
882     itstpcP->Draw("p");
883     //
884     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
885
886     legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
887     legend->AddEntry(itstpcB,"ITS-TPC smooth");
888     legend->Draw();        
889     //picArray->AddLast(pad);
890   }
891 }
892
893
894
895
896 void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber,  AliCDBStorage* pocdbStorage){
897   //
898   // Update OCDB gain
899   //
900   if (pocdbStorage==0) {
901     TString localStorage = "local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
902     pocdbStorage = AliCDBManager::Instance()->GetStorage(localStorage.Data());
903   }
904
905   //
906   // 1. Read gain values
907   //
908   ReadGainGlobal(fileName);
909
910   //
911   // 2. Extract calibration values
912   //
913   AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
914   AnalyzeAttachment(startRunNumber,endRunNumber);
915   AnalyzePadRegionGain();
916   AnalyzeGainMultiplicity();
917   AnalyzeGainChamberByChamber();
918   //
919   AnalyzeGainDipAngle(0); // short pads
920   AnalyzeGainDipAngle(1); // medium pads
921   AnalyzeGainDipAngle(2); // long pads
922   AnalyzeGainDipAngle(3); // absolute calibration on full track
923   //
924   // 3. Make control plots
925   //
926   MakeQAPlot(1.43);  
927
928   //
929   // 4. validate OCDB entries
930   //
931   if(fSwitchOnValidation==kTRUE && ValidateTimeGain()==kFALSE) { 
932     Printf("TPC time gain OCDB parameters out of range!");
933     return;
934   }
935
936   //
937   // 5. Update OCDB
938   //
939   UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage);
940 }
941
942 void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
943   //
944   // read calibration entries from file
945   // 
946   TFile fcalib(fileName);
947   TObject* obj = dynamic_cast<TObject*>(fcalib.Get("TPCCalib"));
948   TObjArray * array = dynamic_cast<TObjArray*>(obj);
949   TDirectory * dir = dynamic_cast<TDirectory*>(obj);
950   if (dir) {
951     fGainMIP    = dynamic_cast<AliTPCcalibTimeGain *>(dir->Get("calibTimeGain"));
952     fGainCosmic = dynamic_cast<AliTPCcalibTimeGain *>(dir->Get("calibTimeGainCosmic"));
953     fGainMult   = dynamic_cast<AliTPCcalibGainMult *>(dir->Get("calibGainMult"));
954   }
955   else if (array){
956     fGainMIP    = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
957     fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
958     fGainMult   = ( AliTPCcalibGainMult *)array->FindObject("calibGainMult");
959   }else{
960     fGainMIP    = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
961     fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
962     fGainMult   = ( AliTPCcalibGainMult *)fcalib.Get("calibGainMult");
963   }
964   if (!fGainMult){
965     TFile calibMultFile("TPCMultObjects.root");
966     fGainMult   = ( AliTPCcalibGainMult *)calibMultFile.Get("calibGainMult");
967   }
968   TH1 * hisT=0;
969   Int_t firstBinA =0, lastBinA=0;
970
971   if (fGainCosmic){ 
972     hisT= fGainCosmic->GetHistGainTime()->Projection(1);
973     firstBinA = hisT->FindFirstBinAbove(2);
974     lastBinA  = hisT->FindLastBinAbove(2);    
975     fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
976     delete hisT;
977   }
978
979   if (fGainMIP){ 
980     hisT= fGainMIP->GetHistGainTime()->Projection(1);
981     firstBinA = hisT->FindFirstBinAbove(2);
982     lastBinA  = hisT->FindLastBinAbove(2);    
983     fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
984     delete hisT;
985   }
986
987 }
988
989 Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit,  Float_t FPtoMIPratio){
990   //
991   // Analyze gain - produce the calibration graphs
992  //
993
994   // 1.) try to create MIP spline
995   if (fGainMIP) 
996   {
997     fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
998     fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
999     fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
1000     //
1001     fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
1002     if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
1003     if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
1004     if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
1005     fGainArray->AddAt(fFitMIP,0);
1006   } 
1007
1008   // 2.) try to create Cosmic spline
1009   if (fGainCosmic)
1010   {
1011     fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
1012     fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100);    // only Fermi-Plateau muons
1013     //
1014     fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
1015     if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
1016     //
1017     if (fGraphCosmic) {
1018       for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
1019         fGraphCosmic->GetY()[i] /= FPtoMIPratio;
1020         fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
1021       }
1022     }
1023     //
1024     if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
1025     if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
1026     fGainArray->AddAt(fFitCosmic,1);
1027   }
1028   // with naming convention and backward compatibility
1029   fGainArray->AddAt(fGraphMIP,2);
1030   fGainArray->AddAt(fGraphCosmic,3);
1031   //
1032   // 3.) Add HV and PT correction parameterization which was used
1033   //
1034   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1035   if (param->GetGainSlopesHV())  fGainArray->AddLast(param->GetGainSlopesHV());
1036   if (param->GetGainSlopesPT())  fGainArray->AddLast(param->GetGainSlopesPT());
1037   //
1038   // 4.) Add the RecoParam and ClusterParam - for compatibility checks -deffrent sets of paramters can invalidate calibration 
1039   //
1040   AliTPCClusterParam *clParam =   AliTPCcalibDB::Instance()->GetClusterParam();
1041   TObjArray *recoParams = new TObjArray(4) ;
1042   for (Int_t i=0;i<4;i++) recoParams->AddAt(AliTPCcalibDB::Instance()->GetRecoParam(i),i);
1043   fGainArray->AddLast(clParam);
1044   fGainArray->AddLast(recoParams);
1045   //
1046   cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
1047   return kTRUE;
1048
1049 }
1050
1051 Bool_t AliTPCPreprocessorOffline::AnalyzeAttachment(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesFit) {
1052   //
1053   // determine slope as a function of mean driftlength
1054   //
1055   if(!fGainMIP) return kFALSE;
1056
1057   fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
1058   //
1059   fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
1060   fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
1061   //
1062   fGainMIP->GetHistGainTime()->GetAxis(3)->SetRangeUser(125,250);// only full tracking region (driftlength)
1063   fGainMIP->GetHistGainTime()->GetAxis(0)->SetRangeUser(1.5,3.5);// only full tracking region (driftlength)
1064   //
1065   TH3D * hist = fGainMIP->GetHistGainTime()->Projection(1, 0, 3);
1066   //
1067   Double_t *xvec = new Double_t[hist->GetNbinsX()];
1068   Double_t *yvec = new Double_t[hist->GetNbinsX()];
1069   Double_t *xerr = new Double_t[hist->GetNbinsX()];
1070   Double_t *yerr = new Double_t[hist->GetNbinsX()];
1071   Int_t counter  = 0;
1072   //
1073   for(Int_t i=1; i < hist->GetNbinsX(); i++) {
1074     Int_t nsum=0;
1075     Int_t imin   =  i;
1076     Int_t imax   =  i;    
1077     for (Int_t idelta=0; idelta<5; idelta++){
1078       //
1079       imin   =  TMath::Max(i-idelta,1);
1080       imax   =  TMath::Min(i+idelta,hist->GetNbinsX());
1081       nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
1082       //if (nsum==0) break;
1083       if (nsum>minEntriesFit) break;
1084     }
1085     if (nsum<minEntriesFit) continue;
1086     //
1087     fGainMIP->GetHistGainTime()->GetAxis(1)->SetRangeUser(hist->GetXaxis()->GetBinCenter(imin-1),hist->GetXaxis()->GetBinCenter(imax+1)); // define time range
1088     TH2D * histZdep = fGainMIP->GetHistGainTime()->Projection(0,3);
1089     TObjArray arr;
1090     histZdep->FitSlicesY(0,0,-1,0,"QNR",&arr);
1091     TH1D * driftDep = (TH1D*)arr.At(1);
1092     delete histZdep;
1093     //TGraphErrors * driftDep = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,3,100,1,0.,1);
1094     /*if (driftDep->GetN() < 4) {
1095       delete driftDep;
1096       continue;
1097       }*/
1098     //
1099     //TObjArray arr;
1100     //
1101     TF1 pol1("polynom1","pol1",125,240);
1102     //driftDep->Fit(&pol1,"QNRROB=0.8");
1103     driftDep->Fit(&pol1,"QNR");
1104     xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin-1)+hist->GetXaxis()->GetBinCenter(imax+1));
1105     yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
1106     xerr[counter] = hist->GetXaxis()->GetBinCenter(imax+1)-hist->GetXaxis()->GetBinCenter(imin-1);
1107     yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
1108     counter++;
1109     //
1110     //delete driftDep;
1111   }
1112   //
1113   fGraphAttachmentMIP = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
1114   if (fGraphAttachmentMIP) fGraphAttachmentMIP->SetName("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");// set proper names according to naming convention
1115   fGainArray->AddLast(fGraphAttachmentMIP);
1116   //
1117   delete [] xvec;
1118   delete [] yvec;
1119   delete [] xerr;
1120   delete [] yerr;
1121   delete hist;
1122   //
1123   if (counter < 1) return kFALSE;
1124   return kTRUE;
1125
1126 }
1127
1128
1129 Bool_t AliTPCPreprocessorOffline::AnalyzePadRegionGain(){
1130   //
1131   // Analyze gain for different pad regions - produce the calibration graphs 0,1,2
1132   //
1133   if (fGainMult) 
1134   {
1135     TH2D * histQmax = (TH2D*) fGainMult->GetHistPadEqual()->Projection(0,2);
1136     TH2D * histQtot = (TH2D*) fGainMult->GetHistPadEqual()->Projection(1,2);
1137     //
1138     TObjArray arr;
1139     histQmax->FitSlicesY(0,0,-1,0,"QNR",&arr);
1140     Double_t xMax[3] = {0,1,2};
1141     Double_t yMax[3]    = {((TH1D*)arr.At(1))->GetBinContent(1),
1142                            ((TH1D*)arr.At(1))->GetBinContent(2),
1143                            ((TH1D*)arr.At(1))->GetBinContent(3)};
1144     Double_t yMaxErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1145                            ((TH1D*)arr.At(1))->GetBinError(2),
1146                            ((TH1D*)arr.At(1))->GetBinError(3)};
1147     TGraphErrors * fitPadRegionQmax = new TGraphErrors(3, xMax, yMax, 0, yMaxErr);
1148     //
1149     histQtot->FitSlicesY(0,0,-1,0,"QNR",&arr);
1150     Double_t xTot[3] = {0,1,2};
1151     Double_t yTot[3]    = {((TH1D*)arr.At(1))->GetBinContent(1),
1152                            ((TH1D*)arr.At(1))->GetBinContent(2),
1153                            ((TH1D*)arr.At(1))->GetBinContent(3)};
1154     Double_t yTotErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1155                            ((TH1D*)arr.At(1))->GetBinError(2),
1156                            ((TH1D*)arr.At(1))->GetBinError(3)};
1157     TGraphErrors * fitPadRegionQtot = new TGraphErrors(3, xTot, yTot, 0, yTotErr);
1158     //
1159     fitPadRegionQtot->SetName("TGRAPHERRORS_MEANQTOT_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1160     fitPadRegionQmax->SetName("TGRAPHERRORS_MEANQMAX_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1161     //
1162     fGainArray->AddLast(fitPadRegionQtot);
1163     fGainArray->AddLast(fitPadRegionQmax);
1164     return kTRUE;
1165   } 
1166   return kFALSE;
1167
1168 }
1169
1170
1171 Bool_t AliTPCPreprocessorOffline::AnalyzeGainDipAngle(Int_t padRegion)  {
1172   //
1173   // Analyze gain as a function of multiplicity and produce calibration graphs
1174   // padRegion -- 0: short, 1: medium, 2: long, 3: absolute calibration of full track
1175   //
1176   Int_t kMarkers[10]={25,24,20,21,22};
1177   Int_t kColors[10]={1,2,4,3,6};
1178   if (!fGainMult) return kFALSE;
1179   if (!(fGainMult->GetHistTopology())) return kFALSE;
1180   const Double_t kMinStat=100;
1181   //
1182   // "dEdxRatioMax","dEdxRatioTot","padType","mult","driftlength"
1183   TObjArray arrMax;
1184   TObjArray arrTot;
1185   //
1186   TH2D * histQmax = 0x0;
1187   TH2D * histQtot = 0x0;
1188   fGainMult->GetHistPadEqual()->GetAxis(4)->SetRangeUser(-0.85,0.85);
1189   fGainMult->GetHistTopology()->GetAxis(2)->SetRangeUser(-0.85,0.85);
1190   if (padRegion < 3) {
1191     fGainMult->GetHistPadEqual()->GetAxis(2)->SetRangeUser(padRegion,padRegion); // short,medium,long
1192     histQmax = (TH2D*) fGainMult->GetHistPadEqual()->Projection(0,4);
1193     histQtot = (TH2D*) fGainMult->GetHistPadEqual()->Projection(1,4);
1194   } else {
1195     fGainMult->GetHistTopology()->GetAxis(1)->SetRangeUser(1,1); //Qmax
1196     histQmax = (TH2D*) fGainMult->GetHistTopology()->Projection(0,2);
1197     histQmax->SetName("fGainMult_GetHistPadEqual_11");
1198     fGainMult->GetHistTopology()->GetAxis(1)->SetRangeUser(0,0); //Qtot
1199     histQtot = (TH2D*) fGainMult->GetHistTopology()->Projection(0,2);
1200     histQtot->SetName("fGainMult_GetHistPadEqual_00");
1201   }
1202   //  
1203   
1204   if (histQmax->GetEntries()<=kMinStat || histQtot->GetEntries()<=kMinStat) {
1205     AliError(Form("hisQtot.GetEntries()=%f",histQtot->GetEntries()));
1206     AliError(Form("hisQmax.GetEntries()=%f",histQmax->GetEntries()));
1207     return kFALSE;
1208   }
1209
1210   TGraphErrors * graphMax = TStatToolkit::MakeStat1D( histQmax,0,0.8,4,kMarkers[padRegion],kColors[padRegion]);
1211   TGraphErrors * graphTot = TStatToolkit::MakeStat1D( histQtot,0,0.8,4,kMarkers[padRegion],kColors[padRegion]);
1212   delete histQmax;
1213   delete histQtot;
1214
1215   //
1216   const char* names[4]={"SHORT","MEDIUM","LONG","ABSOLUTE"};
1217   //
1218   Double_t meanMax = TMath::Mean(graphMax->GetN(), graphMax->GetY());
1219   Double_t meanTot = TMath::Mean(graphTot->GetN(), graphTot->GetY());
1220   if (meanMax<=0 || meanTot<=0){
1221     AliError(Form("meanMax=%f",meanMax));
1222     AliError(Form("meanTot=%f",meanTot));
1223     return kFALSE;
1224   }
1225   //
1226   for (Int_t ipoint=0; ipoint<graphMax->GetN(); ipoint++) {graphMax->GetY()[ipoint]/=meanMax;}
1227   for (Int_t ipoint=0; ipoint<graphTot->GetN(); ipoint++) {graphTot->GetY()[ipoint]/=meanTot;}
1228   //
1229   graphMax->SetNameTitle(Form("TGRAPHERRORS_QMAX_DIPANGLE_%s_BEAM_ALL",names[padRegion]),
1230                         Form("TGRAPHERRORS_QMAX_DIPANGLE_%s_BEAM_ALL",names[padRegion]));
1231   graphTot->SetNameTitle(Form("TGRAPHERRORS_QTOT_DIPANGLE_%s_BEAM_ALL",names[padRegion]),
1232                         Form("TGRAPHERRORS_QTOT_DIPANGLE_%s_BEAM_ALL",names[padRegion]));
1233   //
1234   fGainArray->AddLast(graphMax);
1235   fGainArray->AddLast(graphTot);
1236   //
1237   // Normalization to 1 (mean of the graph.fY --> 1)
1238   //
1239   TF1 * funMax= new TF1("","1++abs(x)++abs(x*x)");
1240   TF1 * funTot= new TF1("","1++abs(x)++abs(x*x)");
1241   graphMax->Fit(funMax,"w","rob=0.9",-0.8,0.8);
1242   graphTot->Fit(funTot,"w","rob=0.9",-0.8,0.8);
1243   funMax->SetNameTitle(Form("TF1_QMAX_DIPANGLE_%s_BEAM_ALL",names[padRegion]),
1244                         Form("TF1_QMAX_DIPANGLE_%s_BEAM_ALL",names[padRegion]));
1245   funTot->SetNameTitle(Form("TF1_QTOT_DIPANGLE_%s_BEAM_ALL",names[padRegion]),
1246                         Form("TF1_QTOT_DIPANGLE_%s_BEAM_ALL",names[padRegion]));
1247
1248   //
1249   fGainArray->AddLast(funMax);
1250   fGainArray->AddLast(funTot);
1251   //
1252   return kTRUE;
1253 }
1254
1255
1256 Bool_t AliTPCPreprocessorOffline::AnalyzeGainMultiplicity() {
1257   //
1258   // Analyze gain as a function of multiplicity and produce calibration graphs
1259   //
1260   if (!fGainMult) return kFALSE;
1261   fGainMult->GetHistGainMult()->GetAxis(3)->SetRangeUser(3,3);
1262   TH2D * histMultMax = fGainMult->GetHistGainMult()->Projection(0,4);
1263   TH2D * histMultTot = fGainMult->GetHistGainMult()->Projection(1,4);
1264   histMultMax->RebinX(4);
1265   histMultTot->RebinX(4);
1266   //
1267   TObjArray arrMax;
1268   TObjArray arrTot;
1269   TF1 fitGaus("fitGaus","gaus(0)",histMultMax->GetYaxis()->GetXmin(),histMultMax->GetYaxis()->GetXmax());
1270   fitGaus.SetParameters(histMultMax->GetEntries()/10., histMultMax->GetMean(2), TMath::Sqrt(TMath::Abs(histMultMax->GetMean(2))));
1271   histMultMax->FitSlicesY(&fitGaus,0,-1,1,"QNRB",&arrMax);
1272   fitGaus.SetParameters(histMultTot->GetEntries()/10., histMultTot->GetMean(2), TMath::Sqrt(TMath::Abs(histMultTot->GetMean(2))));
1273   histMultTot->FitSlicesY(&fitGaus,0,-1,1,"QNRB",&arrTot);
1274   //
1275   TH1D * meanMax = (TH1D*)arrMax.At(1);
1276   TH1D * meanTot = (TH1D*)arrTot.At(1);
1277   Float_t meanMult = histMultMax->GetMean();
1278   delete histMultMax;
1279   delete histMultTot;
1280   if(meanMax->GetBinContent(meanMax->FindBin(meanMult))) {
1281     meanMax->Scale(1./meanMax->GetBinContent(meanMax->FindBin(meanMult)));
1282   }
1283   else {
1284    return kFALSE;
1285   }
1286   if(meanTot->GetBinContent(meanTot->FindBin(meanMult))) {
1287     meanTot->Scale(1./meanTot->GetBinContent(meanTot->FindBin(meanMult)));
1288   }
1289   else {
1290    return kFALSE;
1291   }
1292   Float_t xMultMax[50];
1293   Float_t yMultMax[50];
1294   Float_t yMultErrMax[50];
1295   Float_t xMultTot[50];
1296   Float_t yMultTot[50];
1297   Float_t yMultErrTot[50];
1298   //
1299   Int_t nCountMax = 0;
1300   for(Int_t iBin = 1; iBin < meanMax->GetXaxis()->GetNbins(); iBin++) {
1301     Float_t yValMax = meanMax->GetBinContent(iBin);
1302     if (yValMax < 0.7) continue;
1303     if (yValMax > 1.3) continue;
1304     if (meanMax->GetBinError(iBin)/yValMax > 0.01) continue;
1305     xMultMax[nCountMax] = meanMax->GetXaxis()->GetBinCenter(iBin);
1306     yMultMax[nCountMax] = yValMax;
1307     yMultErrMax[nCountMax] = meanMax->GetBinError(iBin);
1308     nCountMax++;
1309   }
1310   //
1311   if (nCountMax < 10) return kFALSE;
1312   TGraphErrors * fitMultMax = new TGraphErrors(nCountMax, xMultMax, yMultMax, 0, yMultErrMax);
1313   fitMultMax->SetName("TGRAPHERRORS_MEANQMAX_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1314   //
1315   Int_t nCountTot = 0;
1316   for(Int_t iBin = 1; iBin < meanTot->GetXaxis()->GetNbins(); iBin++) {
1317     Float_t yValTot = meanTot->GetBinContent(iBin);
1318     if (yValTot < 0.7) continue;
1319     if (yValTot > 1.3) continue;
1320     if (meanTot->GetBinError(iBin)/yValTot > 0.1) continue;
1321     xMultTot[nCountTot] = meanTot->GetXaxis()->GetBinCenter(iBin);
1322     yMultTot[nCountTot] = yValTot;
1323     yMultErrTot[nCountTot] = meanTot->GetBinError(iBin);
1324     nCountTot++;
1325   }
1326   //
1327   if (nCountTot < 10) return kFALSE;
1328   TGraphErrors *  fitMultTot = new TGraphErrors(nCountTot, xMultTot, yMultTot, 0, yMultErrTot);
1329   fitMultTot->SetName("TGRAPHERRORS_MEANQTOT_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1330   //
1331   fGainArray->AddLast(fitMultMax);
1332   fGainArray->AddLast(fitMultTot);
1333   //
1334   return kTRUE;
1335
1336 }
1337
1338 Bool_t AliTPCPreprocessorOffline::AnalyzeGainChamberByChamber(){
1339   //
1340   // get chamber by chamber gain
1341   //
1342   if (!fGainMult) return kFALSE;
1343   TGraphErrors *grShort  = fGainMult->GetGainPerChamberRobust(0);
1344   TGraphErrors *grMedium = fGainMult->GetGainPerChamberRobust(1);
1345   TGraphErrors *grLong   = fGainMult->GetGainPerChamberRobust(2);
1346   if (grShort==0x0 || grMedium==0x0 || grLong==0x0) {
1347     delete grShort;
1348     delete grMedium;
1349     delete grLong;
1350     return kFALSE;
1351   }
1352
1353   fGainArray->AddLast(grShort);
1354   fGainArray->AddLast(grMedium);
1355   fGainArray->AddLast(grLong);
1356
1357   return kTRUE;
1358 }
1359
1360 void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage *storage){
1361   //
1362   // Update OCDB entry
1363   //
1364   AliCDBMetaData *metaData= new AliCDBMetaData();
1365   metaData->SetObjectClassName("TObjArray");
1366   metaData->SetResponsible("Alexander Kalweit");
1367   metaData->SetBeamPeriod(1);
1368   metaData->SetAliRootVersion("05-24-00"); //root version
1369   metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
1370   AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
1371   storage->Put(fGainArray, id1, metaData);    
1372 }
1373
1374 void AliTPCPreprocessorOffline::MakeQAPlot(Float_t  FPtoMIPratio) {
1375   //
1376   // Make QA plot to visualize results
1377   //
1378   //
1379   //
1380   if (fGraphCosmic) {
1381     TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
1382     canvasCosmic->cd();
1383     TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
1384     gainHistoCosmic->SetDirectory(0);
1385     gainHistoCosmic->SetName("GainHistoCosmic");
1386     gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
1387     gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1388     gainHistoCosmic->Draw("colz");
1389     fGraphCosmic->SetMarkerStyle(25);
1390     fGraphCosmic->Draw("lp");
1391     fGraphCosmic->SetMarkerStyle(25);
1392     TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
1393     if (grfFitCosmic) {
1394       for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
1395         grfFitCosmic->GetY()[i] *= FPtoMIPratio;        
1396       }
1397       for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
1398         fGraphCosmic->GetY()[i] *= FPtoMIPratio;        
1399       }
1400     }
1401     fGraphCosmic->Draw("lp"); 
1402     if (grfFitCosmic) {
1403       grfFitCosmic->SetLineColor(2);
1404       grfFitCosmic->Draw("lu");
1405     }
1406     fGainArray->AddLast(gainHistoCosmic);
1407     //fGainArray->AddLast(canvasCosmic->Clone());
1408     delete canvasCosmic;    
1409   }
1410   if (fFitMIP) {
1411     TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
1412     canvasMIP->cd();
1413     TH2D * gainHistoMIP    = fGainMIP->GetHistGainTime()->Projection(0,1);
1414     gainHistoMIP->SetName("GainHistoCosmic");
1415     gainHistoMIP->SetDirectory(0);
1416     gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
1417     gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1418     gainHistoMIP->Draw("colz");
1419     fGraphMIP->SetMarkerStyle(25);
1420     fGraphMIP->Draw("lp");
1421     TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
1422     grfFitMIP->SetLineColor(2);
1423     grfFitMIP->Draw("lu");    
1424     fGainArray->AddLast(gainHistoMIP);
1425     //fGainArray->AddLast(canvasMIP->Clone());
1426     delete canvasMIP;  
1427     delete grfFitMIP;
1428   }  
1429 }
1430
1431 void AliTPCPreprocessorOffline::MakeFitTime(){
1432   //
1433   // make aligment fit - store results in the file
1434   //
1435   const Int_t kMinEntries=1000;
1436   MakeChainTime();
1437   MakePrimitivesTime();
1438   if (!fAlignTree) return;
1439   if (fAlignTree->GetEntries()<kMinEntries) return;
1440   fAlignTree->SetAlias("ptype","type");
1441   fAlignTree->SetAlias("hasITS","(1+0)");
1442   fAlignTree->SetAlias("dITS","1-2*(refX<40)");
1443   fAlignTree->SetAlias("isITS","refX>10");
1444   fAlignTree->SetAlias("isVertex","refX<10");
1445   // 
1446   Int_t  npointsMax=30000000;
1447   TStatToolkit toolkit;
1448   Double_t chi2=0;
1449   Int_t    npoints=0;
1450   TVectorD param;
1451   TMatrixD covar;
1452
1453   TString fstringFast="";
1454   fstringFast+="FExBTwistX++";
1455   fstringFast+="FExBTwistY++";
1456   fstringFast+="FAlignRot0D++";
1457   fstringFast+="FAlignTrans0D++";
1458   fstringFast+="FAlignTrans1D++";
1459   //
1460   fstringFast+="hasITS*FAlignTrans0++";
1461   fstringFast+="hasITS*FAlignTrans1++";
1462   fstringFast+="hasITS*FAlignRot0++";
1463   fstringFast+="hasITS*FAlignRot1++";
1464   fstringFast+="hasITS*FAlignRot2++";
1465   //
1466   fstringFast+="dITS*FAlignTrans0++";
1467   fstringFast+="dITS*FAlignTrans1++";
1468   fstringFast+="dITS*FAlignRot0++";
1469   fstringFast+="dITS*FAlignRot1++";
1470   fstringFast+="dITS*FAlignRot2++";
1471   
1472   TCut cutFit="entries>10&&abs(mean)>0.00001&&rms>0";
1473   fAlignTree->SetAlias("err","rms");
1474
1475   TString *strDeltaITS = TStatToolkit::FitPlaneConstrain(fAlignTree,"mean:err", fstringFast.Data(),cutFit, chi2,npoints,param,covar,-1,0, npointsMax, 1);
1476   TObjArray* tokArr = strDeltaITS->Tokenize("++");
1477   tokArr->Print();
1478   delete tokArr;
1479   fAlignTree->SetAlias("fitYFast",strDeltaITS->Data());
1480   delete strDeltaITS;
1481
1482   // 
1483   TVectorD paramC= param;
1484   TMatrixD covarC= covar;
1485   TStatToolkit::Constrain1D(fstringFast,"Trans0D",paramC,covarC,0, 0.1);
1486   TStatToolkit::Constrain1D(fstringFast,"Trans1D",paramC,covarC,0, 0.1);
1487   TStatToolkit::Constrain1D(fstringFast,"TwistX",paramC,covarC,0, 0.1);
1488   TStatToolkit::Constrain1D(fstringFast,"TwistY",paramC,covarC,0, 0.1);
1489   TString strFitConst=TStatToolkit::MakeFitString(fstringFast, paramC,covar);
1490   fAlignTree->SetAlias("fitYFastC",strFitConst.Data());
1491   CreateAlignTime(fstringFast,paramC);
1492
1493 }
1494
1495
1496 void AliTPCPreprocessorOffline::MakeChainTime(){
1497   //
1498   //
1499   //
1500   TFile f("CalibObjects.root");
1501   
1502   //  const char *cdtype[7]={"ITS","TRD","Vertex","TOF","TPC","TPC0","TPC1"};
1503   //const char *cptype[5]={"dy","dz","dsnp","dtheta","d1pt"}; 
1504   const char * hname[5]={"dy","dz","dsnp","dtheta","d1pt"};
1505   Int_t run=0;
1506   AliTPCcalibTime  *calibTime = 0;
1507   TObject* obj = dynamic_cast<TObject*>(f.Get("TPCCalib"));
1508   TObjArray * array = dynamic_cast<TObjArray*>(obj);
1509   TDirectory * dir = dynamic_cast<TDirectory*>(obj);
1510   if (dir) {
1511     calibTime = dynamic_cast<AliTPCcalibTime*>(dir->Get("calibTime"));
1512   }
1513   else if (array){
1514     calibTime = (AliTPCcalibTime *)array->FindObject("calibTime");
1515   } else {
1516     calibTime = (AliTPCcalibTime*)f.Get("calibTime");
1517   }
1518   if (!calibTime) return;
1519   AliTPCCorrectionFit::CreateAlignMaps(AliTracker::GetBz(), run);
1520   TTreeSRedirector *pcstream = new TTreeSRedirector("meanITSVertex.root");
1521   //
1522   Int_t ihis=0;
1523   THnSparse *his = calibTime->GetResHistoTPCITS(ihis);
1524   if (his){
1525     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1526     his->GetAxis(2)->SetRange(0,1000000);
1527     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1528     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
1529   }
1530   ihis=1;
1531   his = calibTime->GetResHistoTPCITS(ihis);
1532   if (his){
1533     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1534     his->GetAxis(2)->SetRange(0,1000000);
1535     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1536     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
1537   }
1538   ihis=2;
1539   his = calibTime->GetResHistoTPCITS(ihis);
1540   if (his){
1541     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1542     his->GetAxis(2)->SetRange(0,1000000);
1543     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1544     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
1545   }
1546   ihis=0;
1547   his = calibTime->GetResHistoTPCvertex(ihis);
1548   if (his){
1549     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1550     his->GetAxis(2)->SetRange(0,1000000);
1551     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1552     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
1553   }
1554   ihis=2;
1555   his = calibTime->GetResHistoTPCvertex(ihis);
1556   if (his){
1557     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1558     his->GetAxis(2)->SetRange(0,1000000);
1559     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1560     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
1561
1562   }
1563   ihis=1;
1564   his = calibTime->GetResHistoTPCvertex(ihis);
1565   if (his){
1566     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1567     his->GetAxis(2)->SetRange(0,1000000);
1568     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1569     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
1570
1571   }
1572   ihis=0;
1573   his = calibTime->GetResHistoTPCTOF(ihis);
1574   if (his){
1575     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1576     his->GetAxis(2)->SetRange(0,1000000);
1577     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1578     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("TOF%s",hname[ihis]),run,0.,ihis,3);
1579
1580   }
1581   ihis=0;
1582   his = calibTime->GetResHistoTPCTRD(ihis);
1583   if (his){
1584     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1585     his->GetAxis(2)->SetRange(0,1000000);
1586     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1587     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("TRD%s",hname[ihis]),run,0.,ihis,3);
1588
1589   }
1590   if (dir || (!dir && !array)) { // the object is taken from a directory or a file
1591     delete calibTime;
1592   }
1593   delete pcstream;
1594 }
1595
1596
1597 Double_t AliTPCPreprocessorOffline::EvalAt(Double_t phi, Double_t refX, Double_t theta, Int_t corr, Int_t ptype){
1598   //
1599   //
1600   //
1601   Double_t sector = 9*phi/TMath::Pi();
1602   if (sector<0) sector+=18;
1603   Double_t y85=AliTPCCorrection::GetCorrSector(sector,85,theta,1,corr);
1604   Double_t y245=AliTPCCorrection::GetCorrSector(sector,245,theta,1,corr);
1605   if (ptype==0) return y85+(y245-y85)*(refX-85.)/(245.-85.);
1606   if (ptype==2) return (y245-y85)/(245.-85.);
1607   return 0;
1608 }
1609
1610
1611 Double_t AliTPCPreprocessorOffline::EvalAtPar(Double_t phi0, Double_t snp, Double_t refX, Double_t theta, Int_t corr, Int_t ptype, Int_t nsteps){
1612   //
1613   // Fit the distortion along the line with the parabolic model
1614   // Parameters:
1615   //  phi0 - phi at the entrance of the TPC
1616   //  snp  - local inclination angle at the entrance of the TPC
1617   //  refX - ref X where the distortion is evanluated
1618   //  theta
1619   //  
1620   static TLinearFitter fitter(3,"pol2"); 
1621   fitter.ClearPoints();
1622   if (nsteps<3) nsteps=3;
1623   Double_t deltaX=(245-85)/(nsteps);
1624   for (Int_t istep=0; istep<(nsteps+1); istep++){
1625     //
1626     Double_t localX =85.+deltaX*istep;
1627     Double_t localPhi=phi0+deltaX*snp*istep;
1628     Double_t sector = 9*localPhi/TMath::Pi();
1629     if (sector<0) sector+=18;
1630     Double_t y=AliTPCCorrection::GetCorrSector(sector,localX,theta,1,corr);
1631     Double_t dlocalX=AliTPCCorrection::GetCorrSector(sector,localX,theta,0,corr);
1632     Double_t x[1]={localX-dlocalX};
1633     fitter.AddPoint(x,y);
1634   }
1635   fitter.Eval();
1636   Double_t par[3];
1637   par[0]=fitter.GetParameter(0);
1638   par[1]=fitter.GetParameter(1);
1639   par[2]=fitter.GetParameter(2);
1640
1641   if (ptype==0) return par[0]+par[1]*refX+par[2]*refX*refX;
1642   if (ptype==2) return par[1]+2*par[2]*refX;
1643   if (ptype==4) return par[2];
1644   return 0;
1645 }
1646
1647
1648
1649
1650
1651
1652 void AliTPCPreprocessorOffline::MakePrimitivesTime(){
1653   //
1654   // Create primitive transformation to fit
1655   //
1656   fAlignTree=new TChain("fit","fit");
1657   fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdy");
1658   fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdsnp");
1659   fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdy");
1660   fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdsnp");
1661   // 
1662   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1663   Double_t bzField=AliTrackerBase::GetBz(); 
1664   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
1665   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
1666   Double_t wtP = -10.0 * (bzField) * vdrift /  ezField ; 
1667   AliTPCExBTwist *fitExBTwistX= new  AliTPCExBTwist;
1668   AliTPCExBTwist *fitExBTwistY= new  AliTPCExBTwist;
1669   AliTPCCalibGlobalMisalignment *trans0   =new  AliTPCCalibGlobalMisalignment;
1670   AliTPCCalibGlobalMisalignment *trans1   =new  AliTPCCalibGlobalMisalignment;
1671   AliTPCCalibGlobalMisalignment *trans0D  =new  AliTPCCalibGlobalMisalignment;
1672   AliTPCCalibGlobalMisalignment *trans1D  =new  AliTPCCalibGlobalMisalignment;
1673   AliTPCCalibGlobalMisalignment *rot0     =new  AliTPCCalibGlobalMisalignment;
1674   AliTPCCalibGlobalMisalignment *rot1     =new  AliTPCCalibGlobalMisalignment;
1675   AliTPCCalibGlobalMisalignment *rot2     =new  AliTPCCalibGlobalMisalignment;
1676   AliTPCCalibGlobalMisalignment *rot3     =new  AliTPCCalibGlobalMisalignment;
1677   //
1678   //
1679   fitExBTwistX->SetXTwist(0.001);
1680   fitExBTwistX->SetOmegaTauT1T2(wtP,1,1);  
1681   //
1682   fitExBTwistY->SetYTwist(0.001);
1683   fitExBTwistY->SetOmegaTauT1T2(wtP,1,1);  
1684   //
1685   TGeoHMatrix *matrixRot = new TGeoHMatrix; 
1686   TGeoHMatrix *matrixX = new TGeoHMatrix; 
1687   TGeoHMatrix *matrixY = new TGeoHMatrix; 
1688   matrixX->SetDx(0.1);
1689   matrixY->SetDy(0.1);
1690   Double_t rotAngles0[9]={0};
1691   Double_t rotAngles1[9]={0};
1692   Double_t rotAngles2[9]={0};
1693   //
1694   Double_t rotAngles3[9]={0};
1695
1696   rotAngles0[0]=1; rotAngles0[4]=1; rotAngles0[8]=1;
1697   rotAngles1[0]=1; rotAngles1[4]=1; rotAngles1[8]=1;
1698   rotAngles2[0]=1; rotAngles2[4]=1; rotAngles2[8]=1;
1699   rotAngles3[0]=1; rotAngles3[4]=1; rotAngles3[8]=1;
1700
1701   rotAngles0[1]=-0.001;rotAngles0[3]=0.001;
1702   rotAngles1[5]=-0.001;rotAngles1[7]=0.001;
1703   rotAngles2[2]=0.001;rotAngles2[6]=-0.001;
1704   rotAngles3[1]=0.001;rotAngles3[3]=-0.001;
1705   matrixRot->SetRotation(rotAngles0);
1706   rot0->SetAlignGlobal(matrixRot);
1707   matrixRot->SetRotation(rotAngles1);
1708   rot1->SetAlignGlobal(matrixRot);
1709   matrixRot->SetRotation(rotAngles2);
1710   rot2->SetAlignGlobal(matrixRot); 
1711   matrixRot->SetRotation(rotAngles3);
1712   rot3->SetAlignGlobalDelta(matrixRot); 
1713   //
1714   trans0->SetAlignGlobal(matrixX);
1715   trans1->SetAlignGlobal(matrixY);
1716   trans0D->SetAlignGlobalDelta(matrixX);
1717   trans1D->SetAlignGlobalDelta(matrixY);
1718   fitExBTwistX->Init();
1719   fitExBTwistY->Init();
1720   //
1721   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistX->Clone()),100);
1722   fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistY->Clone()),101);
1723   //
1724   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot0->Clone()),102);
1725   fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(rot1->Clone()),103);
1726   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot2->Clone()),104);
1727   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot3->Clone()),105);
1728
1729   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0->Clone()),106);
1730   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1->Clone()),107);
1731   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0D->Clone()),108);
1732   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1D->Clone()),109);
1733   //
1734   fAlignTree->SetAlias("FExBTwistX", "AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,100,ptype)+0");
1735   fAlignTree->SetAlias("FExBTwistY","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,101,ptype)+0");
1736   fAlignTree->SetAlias("FAlignRot0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,102,ptype)+0");
1737   fAlignTree->SetAlias("FAlignRot0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,105,ptype)+0");
1738   fAlignTree->SetAlias("FAlignRot1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,103,ptype)+0");
1739   fAlignTree->SetAlias("FAlignRot2","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,104,ptype)+0");
1740   fAlignTree->SetAlias("FAlignTrans0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,106,ptype)+0");
1741   fAlignTree->SetAlias("FAlignTrans1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,107,ptype)+0");
1742   fAlignTree->SetAlias("FAlignTrans0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,108,ptype)+0");
1743   fAlignTree->SetAlias("FAlignTrans1D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,109,ptype)+0");
1744   //
1745   // test fast function
1746   //
1747 //   fAlignTree->Draw("FExBTwistX:ExBTwistX","isITS&&ptype==0&&abs(snp)<0.05","");
1748 //   fAlignTree->Draw("FExBTwistY:ExBTwistY","isITS&&ptype==0&&abs(snp)<0.05","");
1749 //   fAlignTree->Draw("FAlignRot0:alignRot0","isITS&&ptype==0&&abs(snp)<0.05","");
1750 //   fAlignTree->Draw("FAlignRot1:alignRot1","isITS&&ptype==0&&abs(snp)<0.05","");
1751 //   fAlignTree->Draw("FAlignRot2:alignRot2","isITS&&ptype==0&&abs(snp)<0.05","");
1752 //   //
1753 //   fAlignTree->Draw("FAlignTrans0:alignTrans0","isITS&&ptype==0&&abs(snp)<0.05","");
1754 //   fAlignTree->Draw("FAlignTrans1:alignTrans1","isITS&&ptype==0&&abs(snp)<0.05","");
1755
1756
1757
1758
1759 void AliTPCPreprocessorOffline::CreateAlignTime(TString fstring, TVectorD paramC){
1760   //
1761   //
1762   //
1763   //
1764   TGeoHMatrix *matrixDelta     = new TGeoHMatrix; 
1765   TGeoHMatrix *matrixGlobal    = new TGeoHMatrix; 
1766   Double_t rAngles[9];
1767   Int_t index=0;
1768   //
1769   index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans0D");
1770   if (index>=0) matrixDelta->SetDx(paramC[index+1]*0.1);
1771   index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans1D");
1772   if (index>=0) matrixDelta->SetDy(paramC[index+1]*0.1);
1773   rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1774   index=TStatToolkit::GetFitIndex(fstring,"FAlignRot0D");
1775   rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1776   rAngles[5]=0; rAngles[7] =0;
1777   rAngles[2]=0; rAngles[6] =0;
1778   matrixDelta->SetRotation(rAngles);
1779   //
1780   //
1781   //
1782   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans0");
1783   if (index>=0) matrixGlobal->SetDx(paramC[index+1]*0.1);
1784   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans1");
1785   if (index>=0) matrixGlobal->SetDy(paramC[index+1]*0.1);
1786   rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1787   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot0");
1788   rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1789   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot1");  
1790   rAngles[5]=-paramC[index+1]*0.001; rAngles[7]=paramC[index+1]*0.001;
1791   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot2");  
1792   rAngles[2]=paramC[index+1]*0.001; rAngles[6] =-paramC[index+1]*0.001;
1793   matrixGlobal->SetRotation(rAngles);
1794   //
1795   AliTPCCalibGlobalMisalignment *fitAlignTime  =0;
1796   fitAlignTime  =new  AliTPCCalibGlobalMisalignment;
1797   fitAlignTime->SetName("FitAlignTime");
1798   fitAlignTime->SetTitle("FitAlignTime");
1799   fitAlignTime->SetAlignGlobalDelta(matrixDelta);
1800   fitAlignTime->SetAlignGlobal(matrixGlobal);
1801   //
1802   AliTPCExBTwist * fitExBTwist= new  AliTPCExBTwist;
1803   Int_t indexX=TStatToolkit::GetFitIndex(fstring,"ExBTwistX");
1804   Int_t indexY=TStatToolkit::GetFitIndex(fstring,"ExBTwistY");  
1805   fitExBTwist->SetXTwist(0.001*paramC[indexX+1]);  // 1 mrad twist in x
1806   fitExBTwist->SetYTwist(0.001*paramC[indexY+1]);  // 1 mrad twist in x
1807   fitExBTwist->SetName("FitExBTwistTime");
1808   fitExBTwist->SetTitle("FitExBTwistTime"); 
1809   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1810   Double_t bzField=AliTrackerBase::GetBz();
1811   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
1812
1813   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
1814   Double_t wt = -10.0 * (bzField) * vdrift /  ezField ; 
1815   //
1816   fitExBTwist->SetOmegaTauT1T2(wt,1,1);  
1817   fitExBTwist->Init();  
1818
1819   AliTPCComposedCorrection *corrTime =  new AliTPCComposedCorrection;
1820   TObjArray *arr = new TObjArray;
1821   corrTime->SetCorrections(arr);
1822   
1823   corrTime->GetCorrections()->Add(fitExBTwist);
1824   corrTime->GetCorrections()->Add(fitAlignTime);
1825   corrTime->SetName("FitCorrectionTime");
1826   corrTime->SetTitle("FitCorrectionTime");
1827
1828   fitExBTwist->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwist->Clone()),1001);
1829   fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(fitAlignTime->Clone()),1002);
1830   fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(corrTime->Clone()),1003);
1831   
1832   
1833   fAlignTree->SetAlias("ExBTwistTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1001,ptype)+0");
1834   fAlignTree->SetAlias("AlignTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1002,ptype)+0");
1835   fAlignTree->SetAlias("FitCorrectionTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1003,ptype)+0");
1836
1837
1838   TFile *f = new TFile("fitITSVertex.root","update");
1839   corrTime->Write("FitCorrectionTime");
1840   f->Close();
1841 }
1842
1843 //_____________________________________________________________________________
1844 Int_t AliTPCPreprocessorOffline::GetStatus()
1845 {
1846   //get the calibration status
1847   // 0 means OK
1848   // positive numbers invalidate for unknown reasons.
1849   // negative numbers invalidate with a known reason (e.g. low statistics).
1850   // the returned integer has one bit set for every component that failed.
1851
1852   Bool_t enoughStatistics = (fNtracksVdrift>fMinTracksVdrift && fNeventsVdrift>fMinEventsVdrift);
1853   
1854   if (!enoughStatistics) 
1855   {
1856     fCalibrationStatus=-TMath::Abs(fCalibrationStatus);
1857   }
1858
1859   return fCalibrationStatus;
1860 }
1861
1862 /*
1863   Short sequence to acces the calbration entry:
1864   TFile *f = TFile::Open("CalibObjects.root");
1865   AliTPCcalibGainMult      * fGainMult = (AliTPCcalibGainMult      *)f->Get("TPCCalib/calibGainMult");
1866   
1867  
1868 */