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