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