]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessorOffline.cxx
Empirical cut to increase robustness
[u/mrichter/AliRoot.git] / TPC / 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 "TGraph.h"
55 #include "TMultiGraph.h"
56 #include "TCanvas.h"
57 #include "THnSparse.h"
58 #include "TLegend.h"
59 #include "TPad.h"
60 #include "TH2D.h"
61 #include "TH3D.h"
62 #include "AliTPCROC.h"
63 #include "AliTPCCalROC.h"
64 #include "AliESDfriend.h"
65 #include "AliTPCcalibTime.h"
66 #include "AliSplineFit.h"
67 #include "AliCDBMetaData.h"
68 #include "AliCDBId.h"
69 #include "AliCDBManager.h"
70 #include "AliCDBStorage.h"
71 #include "AliTPCcalibBase.h"
72 #include "AliTPCcalibDB.h"
73 #include "AliTPCcalibDButil.h"
74 #include "AliRelAlignerKalman.h"
75 #include "AliTPCParamSR.h"
76 #include "AliTPCcalibTimeGain.h"
77 #include "AliTPCcalibGainMult.h"
78 #include "AliSplineFit.h"
79 #include "AliTPCComposedCorrection.h"
80 #include "AliTPCExBTwist.h"
81 #include "AliTPCCalibGlobalMisalignment.h"
82 #include "TStatToolkit.h"
83 #include "TChain.h"
84 #include "TCut.h"
85 #include "AliTrackerBase.h"
86 #include "AliTPCPreprocessorOffline.h"
87
88
89 ClassImp(AliTPCPreprocessorOffline)
90
91 AliTPCPreprocessorOffline::AliTPCPreprocessorOffline():
92   TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
93   fMinEntries(500),                      // minimal number of entries for fit
94   startRun(0),                         // start Run - used to make fast selection in THnSparse
95   endRun(0),                           // end   Run - used to make fast selection in THnSparse
96   startTime(0),                        // startTime - used to make fast selection in THnSparse
97   endTime(0),                          // endTime   - used to make fast selection in THnSparse
98   ocdbStorage(""),                   // path to the OCDB storage
99   fVdriftArray(new TObjArray),
100   fTimeDrift(0),
101   fGraphMIP(0),                // graph time dependence of MIP
102   fGraphCosmic(0),             // graph time dependence at Plateu
103   fGraphAttachmentMIP(0),
104   fFitMIP(0),                  // fit of dependence - MIP
105   fFitCosmic(0),               // fit of dependence - Plateu
106   fGainArray(new TObjArray),               // array to be stored in the OCDB
107   fGainMIP(0),          // calibration component for MIP
108   fGainCosmic(0),       // calibration component for cosmic
109   fGainMult(0),
110   fAlignTree(0),        // alignment tree
111   fSwitchOnValidation(kFALSE) // flag to switch on validation of OCDB parameters
112 {
113   //
114   // default constructor
115   //
116 }
117
118 AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
119   //
120   // Destructor
121   //
122 }
123
124
125
126
127 void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const  timeDrift){
128   //
129   // find the fist and last run
130   //
131   TObjArray *hisArray =timeDrift->GetHistoDrift();
132   {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
133     THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
134     if (addHist->GetEntries()<fMinEntries) continue;
135     if (!addHist) continue;
136     TH1D* histo    =addHist->Projection(3);
137     TH1D* histoTime=addHist->Projection(0);
138     printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
139
140     if (startRun<=0){ 
141       startRun=histo->FindFirstBinAbove(0);
142       endRun  =histo->FindLastBinAbove(0);
143     }else{
144       startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun);
145       endRun  =TMath::Max(histo->FindLastBinAbove(0),endRun);
146     }
147     if (startTime==0){ 
148       startTime=histoTime->FindFirstBinAbove(0);
149       endTime  =histoTime->FindLastBinAbove(0);
150     }else{
151       startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime);
152       endTime  =TMath::Max(histoTime->FindLastBinAbove(0),endTime);
153     }
154     delete histo;
155     delete histoTime;
156   }}
157   if (startRun<0) startRun=0;
158   if (endRun<0) endRun=100000000;
159   printf("Run range  :\t%d-%d\n", startRun, endRun);
160   printf("Time range :\t%d-%d\n", startTime, endTime);
161
162 }
163
164
165
166 void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, TString pocdbStorage){
167   //
168   // make calibration of the drift velocity
169   // Input parameters:
170   //      file                   - the location of input file
171   //      ustartRun, uendrun     - run validity period 
172   //      pocdbStorage           - path to hte OCDB storage
173   //                             - if empty - local storage 'pwd' uesed
174   if (pocdbStorage.Length()>0) ocdbStorage=pocdbStorage;
175   else
176   ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
177   //
178   // 1. Initialization and run range setting
179   TFile fcalib(file);
180   TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
181   if (array){
182     fTimeDrift = (AliTPCcalibTime *)array->FindObject("calibTime");
183   } else {
184     fTimeDrift = (AliTPCcalibTime*)fcalib.Get("calibTime");
185   }
186   if(!fTimeDrift) return;
187
188   startRun=ustartRun;
189   endRun=ustartRun; 
190   TObjArray *hisArray =fTimeDrift->GetHistoDrift();  
191   GetRunRange(fTimeDrift);
192   for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
193     THnSparse* addHist=(THnSparse*)hisArray->At(i);
194     if (!addHist) continue;
195     if (startTime<endTime) addHist->GetAxis(0)->SetRange(startTime-1,endTime+1);
196     if (startRun<endRun) addHist->GetAxis(3)->SetRange(startRun-1,endRun+1);
197   }
198   //
199   //
200   // 2. extraction of the information
201   //
202   fVdriftArray = new TObjArray();
203   AddAlignmentGraphs(fVdriftArray,fTimeDrift);
204   AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
205   AddLaserGraphs(fVdriftArray,fTimeDrift);
206   //
207   // 3. Append QA plots
208   //
209   MakeDefaultPlots(fVdriftArray,fVdriftArray);
210
211   //
212   // 4. validate OCDB entries
213   //
214   if(fSwitchOnValidation==kTRUE && ValidateTimeDrift()==kFALSE) { 
215     Printf("TPC time drift OCDB parameters out of range!");
216     return;
217   }
218   //
219   //4.b make alignment
220   //
221   MakeFitTime();
222   TFile * ftime= TFile::Open("fitITSVertex.root");
223   if (ftime){
224     TObject * alignmentTime=ftime->Get("FitCorrectionTime");
225     if (alignmentTime) fVdriftArray->AddLast(alignmentTime);
226   }
227   //
228   //
229   // 5. update of OCDB
230   //
231   //
232   UpdateOCDBDrift(ustartRun,uendRun,ocdbStorage);
233 }
234
235 void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun,  const char* storagePath ){
236   //
237   // Update OCDB 
238   //
239   AliCDBMetaData *metaData= new AliCDBMetaData();
240   metaData->SetObjectClassName("TObjArray");
241   metaData->SetResponsible("Marian Ivanov");
242   metaData->SetBeamPeriod(1);
243   metaData->SetAliRootVersion("05-25-01"); //root version
244   metaData->SetComment("Calibration of the time dependence of the drift velocity");
245   AliCDBId* id1=NULL;
246   id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
247   AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
248   gStorage->Put(fVdriftArray, (*id1), metaData);
249 }
250
251 Bool_t AliTPCPreprocessorOffline::ValidateTimeGain(Double_t minGain, Double_t maxGain)
252 {
253   //
254   // Validate time gain corrections 
255   //
256   Printf("ValidateTimeGain..." );
257
258   TGraphErrors *gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
259   if(!gr) return kFALSE;
260   if(gr->GetN()<1) return kFALSE;
261
262   // check whether gain in the range
263   for(Int_t iPoint=0; iPoint<gr->GetN(); iPoint++) 
264   {
265     if(gr->GetY()[iPoint] < minGain || gr->GetY()[iPoint] > maxGain)  
266       return kFALSE;
267   }
268
269 return kTRUE;
270 }
271
272
273 Bool_t AliTPCPreprocessorOffline::ValidateTimeDrift(Double_t maxVDriftCorr)
274 {
275   //
276   // Validate time drift velocity corrections 
277   //
278   Printf("ValidateTimeDrift..." );
279
280   TGraphErrors* gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
281   if(!gr) return kFALSE;
282   if(gr->GetN()<1) return kFALSE;
283
284   // check whether drift velocity corrections in the range
285   for(Int_t iPoint = 0; iPoint<gr->GetN(); iPoint++) 
286   {
287     if(TMath::Abs(gr->GetY()[iPoint]) > maxVDriftCorr)  
288       return kFALSE;
289   }
290
291 return kTRUE;
292 }
293
294 void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
295   //
296   //  update the OCDB entry for the nominal time0
297   //
298   //
299   //  AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
300   AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
301   TGraphErrors *grT =  (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
302   Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
303   Double_t deltaT   = deltaTcm/param->GetDriftV();
304   paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
305   paramNew->Update();
306
307   AliCDBMetaData *metaData= new AliCDBMetaData();
308   metaData->SetObjectClassName("TObjArray");
309   metaData->SetResponsible("Marian Ivanov");
310   metaData->SetBeamPeriod(1);
311   metaData->SetAliRootVersion("05-25-02"); //root version
312   metaData->SetComment("Updated calibration of nominal time 0");
313   AliCDBId* id1=NULL;
314   id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
315   AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
316   gStorage->Put(param, (*id1), metaData);
317
318 }
319
320
321 void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
322   //
323   // Print the names of the entries in array
324   //
325   Int_t entries = array->GetEntries();
326   for (Int_t i=0; i<entries; i++){
327     if (!array->At(i)) continue;
328     printf("%d\t %s\n", i,  array->At(i)->GetName());
329   }
330 }
331
332
333
334 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
335   // 2 filters:
336   //    1. filter graph - error cut errSigmaCut
337   //    2. filter graph - medianCutAbs around median
338   //
339   // errSigmaCut   - cut on error
340   // medianCutAbs  - cut on value around median
341   Double_t dummy=0;               //   
342   //
343   // 1. filter graph - error cut errSigmaCut
344   //              
345   TGraphErrors *graphF; 
346   graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
347   delete graph;
348   if (!graphF) return 0;
349   graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
350   delete graphF;
351   if (!graph) return 0;
352   //
353   // filter graph - kMedianCutAbs around median
354   // 
355   graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
356   delete graph;
357   if (!graphF) return 0;
358   graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
359   delete graphF;
360   if (!graph) return 0;
361   return graph;
362 }
363
364
365
366 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
367   //
368   // filter outlyer measurement
369   // Only points around median +- cut filtered 
370   //
371   if (!graph) return  0;
372   Int_t kMinPoints=2;
373   Int_t npoints0 = graph->GetN();
374   Int_t npoints=0;
375   Float_t  rmsY=0;
376   Double_t *outx=new Double_t[npoints0];
377   Double_t *outy=new Double_t[npoints0];
378   Double_t *errx=new Double_t[npoints0];
379   Double_t *erry=new Double_t[npoints0];
380   //
381   //
382   if (npoints0<kMinPoints) return 0;
383   for (Int_t iter=0; iter<3; iter++){
384     npoints=0;
385     for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
386       if (graph->GetY()[ipoint]==0) continue;
387       if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;  
388       outx[npoints]  = graph->GetX()[ipoint];
389       outy[npoints]  = graph->GetY()[ipoint];
390       errx[npoints]  = graph->GetErrorX(ipoint);
391       erry[npoints]  = graph->GetErrorY(ipoint);
392       npoints++;
393     }
394     if (npoints<=1) break;
395     medianY  =TMath::Median(npoints,outy);
396     rmsY   =TMath::RMS(npoints,outy);
397   }
398   TGraphErrors *graphOut=0;
399   if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry); 
400   delete []outx;
401   delete []outy;
402   delete []errx;
403   delete []erry;
404   return graphOut;
405 }
406
407
408 void AliTPCPreprocessorOffline::AddHistoGraphs(  TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
409   //
410   // Add graphs corresponding to the alignment
411   //
412   const Double_t kErrSigmaCut=5;      // error sigma cut - for filtering
413   const Double_t kMedianCutAbs=0.03;  // error sigma cut - for filtering
414   //
415   TObjArray * array=timeDrift->GetHistoDrift();
416   if (array){
417     THnSparse* hist=NULL;
418     // 2.a) cosmics with different triggers
419     for (Int_t i=0; i<array->GetEntriesFast();i++){
420       hist=(THnSparseF*)array->UncheckedAt(i);
421       if(!hist) continue;
422       if (hist->GetEntries()<minEntries) continue;
423       //hist->Print();
424       TString name=hist->GetName();
425       Int_t dim[4]={0,1,2,3};
426       THnSparse* newHist=hist->Projection(4,dim);
427       newHist->SetName(name);
428       TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
429       printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
430       Int_t pos=name.Index("_");
431       name=name(pos,name.Capacity()-pos);
432       TString graphName=graph->ClassName();
433       graphName+=name;
434       graphName.ToUpper();
435       //
436       graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
437       if (!graph) {
438         printf("Graph =%s filtered out\n", name.Data());
439         continue;
440       }
441       //
442       if (graph){
443         graph->SetMarkerStyle(i%8+20);
444         graph->SetMarkerColor(i%7);
445         graph->GetXaxis()->SetTitle("Time");
446         graph->GetYaxis()->SetTitle("v_{dcor}");
447         graph->SetName(graphName);
448         graph->SetTitle(graphName);
449         printf("Graph %d\t=\t%s\n", i, graphName.Data());
450         vdriftArray->Add(graph);
451       }
452     }
453   }
454 }
455
456
457
458
459 void AliTPCPreprocessorOffline::AddAlignmentGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
460   //
461   // Add graphs corresponding to alignment to the object array
462   //
463   TObjArray *arrayITS=0;
464   TObjArray *arrayTOF=0;
465   TObjArray *arrayTRD=0;
466   TMatrixD *mstatITS=0;
467   TMatrixD *mstatTOF=0;
468   TMatrixD *mstatTRD=0;
469   //
470   arrayITS=timeDrift->GetAlignITSTPC();
471   arrayTRD=timeDrift->GetAlignTRDTPC();
472   arrayTOF=timeDrift->GetAlignTOFTPC();
473
474   if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.9,50,0.025);
475   if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.9,1000,0.025);
476   if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.9,50,0.025);
477   //
478   TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
479   TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
480   TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
481   TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
482   TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
483   TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
484
485   TObjArray * arrayTRDP= 0x0;
486   TObjArray * arrayTRDM= 0x0;
487   TObjArray * arrayTRDB= 0x0;
488   arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
489   arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
490   arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
491   //
492   //
493   Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
494   TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
495                          arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
496                          arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
497   TString   grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
498                          "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
499                          "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
500   TString   grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
501                       "DELTAX", "DELTAY", "DELTAZ",
502                       "DRIFTVD", "T0", "VDGY"};
503
504   
505   TVectorD vX(entries);
506   TVectorD vY(entries);
507   TVectorD vEx(entries);
508   TVectorD vEy(entries);
509   TObjArray *arr=0;
510   for (Int_t iarray=0; iarray<12; iarray++){
511     arr = arrays[iarray];
512     if (arr==0) continue;
513     for (Int_t ipar=0; ipar<9; ipar++){      
514       Int_t counter=0;
515       for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
516         AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
517         if (!kalman) continue;
518         vX[counter]=kalman->GetTimeStamp();
519         vY[counter]=(*(kalman->GetState()))[ipar];
520         if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
521         vEx[counter]=0;
522         vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
523         counter++;
524       }
525     
526       TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
527                                           vY.GetMatrixArray(),
528                                           vEx.GetMatrixArray(),
529                                           vEy.GetMatrixArray());
530       TString grName=grnames[iarray];
531       grName+="_TPC_";
532       grName+=grpar[ipar];
533       graph->SetName(grName.Data());
534       vdriftArray->AddLast(graph);
535     }
536   }  
537 }
538
539
540
541
542 void AliTPCPreprocessorOffline::AddLaserGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
543   //
544   // add graphs for laser
545   //
546   const Double_t delayL0L1 = 0.071;  //this is hack for 1/2 weeks
547   THnSparse *hisN=0;
548   TGraphErrors *grLaser[6]={0,0,0,0,0,0};
549   hisN = timeDrift->GetHistVdriftLaserA(0);
550   if (timeDrift->GetHistVdriftLaserA(0)){
551     grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
552     grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
553     vdriftArray->AddLast(grLaser[0]);
554   }    
555   if (timeDrift->GetHistVdriftLaserA(1)){
556     grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
557     grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
558     vdriftArray->AddLast(grLaser[1]);
559   }    
560   if (timeDrift->GetHistVdriftLaserA(2)){
561     grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
562     grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
563     vdriftArray->AddLast(grLaser[2]);
564   }    
565   if (timeDrift->GetHistVdriftLaserC(0)){
566     grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
567     grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
568     vdriftArray->AddLast(grLaser[3]);
569   }    
570   if (timeDrift->GetHistVdriftLaserC(1)){
571     grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
572     grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
573     vdriftArray->AddLast(grLaser[4]);
574   }    
575   if (timeDrift->GetHistVdriftLaserC(2)){
576     grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
577     grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");    
578     vdriftArray->AddLast(grLaser[5]);
579   }    
580   for (Int_t i=0; i<6;i++){
581     if (grLaser[i]) {
582       SetDefaultGraphDrift(grLaser[i], 1,(i+20));
583       grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
584     }
585   }
586 }
587  
588  
589 TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
590   //
591   // Make graph with mean values and rms
592   //
593   hisN->GetAxis(itime)->SetRange(0,100000000);
594   hisN->GetAxis(ival)->SetRange(0,100000000);
595   TH1 * hisT      = hisN->Projection(itime);
596   TH1 * hisV      = hisN->Projection(ival);
597   //
598   Int_t firstBinA = hisT->FindFirstBinAbove(2);
599   Int_t lastBinA  = hisT->FindLastBinAbove(2);    
600   Int_t firstBinV = hisV->FindFirstBinAbove(0);
601   Int_t lastBinV  = hisV->FindLastBinAbove(0);    
602   hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
603   hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
604   Int_t entries=0;
605   for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
606     Double_t cont = hisT->GetBinContent(ibin);
607     if (cont<minEntries) continue;
608     entries++;
609   }
610   TVectorD vecTime(entries);
611   TVectorD vecMean0(entries);
612   TVectorD vecRMS0(entries);
613   TVectorD vecMean1(entries);
614   TVectorD vecRMS1(entries);
615   entries=0;
616   for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
617       Double_t cont = hisT->GetBinContent(ibin);
618       if (cont<minEntries) continue;
619       //hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
620       Int_t minBin = ibin-1;
621       Int_t maxBin = ibin+1;
622       if(minBin <= 0) minBin = 1;
623       if(maxBin >= hisN->GetAxis(itime)->GetNbins()) maxBin = hisN->GetAxis(itime)->GetNbins()-1;
624       hisN->GetAxis(itime)->SetRange(minBin,maxBin);
625       
626       Double_t time = hisT->GetBinCenter(ibin);
627       TH1 * his = hisN->Projection(ival);
628       Double_t nentries0= his->GetBinContent(his->FindBin(0));
629       if (cont-nentries0<minEntries) continue;
630       //
631       his->SetBinContent(his->FindBin(0),0);
632       vecTime[entries]=time;
633       vecMean0[entries]=his->GetMean()+offset;
634       vecMean1[entries]=his->GetMeanError();
635       vecRMS0[entries] =his->GetRMS();
636       vecRMS1[entries] =his->GetRMSError();
637       delete his;  
638       entries++;
639   }
640   delete hisT;
641   delete hisV;
642   TGraphErrors * graph =  new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(),                                    0, vecMean1.GetMatrixArray());
643   return graph;
644 }
645
646
647
648
649
650
651
652
653 void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
654   //
655   // Set default style for QA views
656   //
657   graph->GetXaxis()->SetTimeDisplay(kTRUE);
658   graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
659   graph->SetMaximum( 0.025);
660   graph->SetMinimum(-0.025);
661   graph->GetXaxis()->SetTitle("Time");
662   graph->GetYaxis()->SetTitle("v_{dcorr}");
663   //
664   graph->GetYaxis()->SetLabelSize(0.03);
665   graph->GetXaxis()->SetLabelSize(0.03);
666   //
667   graph->GetXaxis()->SetNdivisions(10,5,0);
668   graph->GetYaxis()->SetNdivisions(10,5,0);
669   //
670   graph->GetXaxis()->SetLabelOffset(0.02);
671   graph->GetYaxis()->SetLabelOffset(0.005);
672   //
673   graph->GetXaxis()->SetTitleOffset(1.3);
674   graph->GetYaxis()->SetTitleOffset(1.2);
675   //
676   graph->SetMarkerColor(color);
677   graph->SetLineColor(color);
678   graph->SetMarkerStyle(style);
679 }
680
681 void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
682   // 
683   // Set default pad style for QA
684   // 
685   pad->SetTicks(1,1);
686   pad->SetMargin(mx0,mx1,my0,my1);
687 }
688
689
690 void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray * /*picArray*/){
691   //
692   // 0. make a default QA plots
693   // 1. Store them in the array
694   //
695   //
696   Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
697   //
698   TGraphErrors* laserA       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
699   TGraphErrors* laserC       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
700   TGraphErrors* cosmic       =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
701   TGraphErrors* cross        =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
702   TGraphErrors* itstpcP       =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
703   TGraphErrors* itstpcM       =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
704   TGraphErrors* itstpcB       =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
705   //
706   if (laserA)  SetDefaultGraphDrift(laserA,2,25);
707   if (laserC)  SetDefaultGraphDrift(laserC,4,26);
708   if (cosmic)  SetDefaultGraphDrift(cosmic,3,27);
709   if (cross)   SetDefaultGraphDrift(cross,4,28);
710   if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
711   if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
712   if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
713   //
714   //
715   TPad *pad=0;
716   //
717   // Laser-Laser
718   //
719   if (laserA&&laserC){
720     pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
721     laserA->Draw("alp");
722     SetPadStyle(pad,mx0,mx1,my0,my1);
723     laserA->Draw("apl");
724     laserC->Draw("p");
725     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
726     legend->AddEntry(laserA,"Laser A side");
727     legend->AddEntry(laserC,"Laser C side");
728     legend->Draw();    
729     //picArray->AddLast(pad);
730   }
731
732   if (itstpcP&&itstpcM&&itstpcB){
733     pad = new TCanvas("ITSTPC","ITSTPC");
734     itstpcP->Draw("alp");
735     SetPadStyle(pad,mx0,mx1,my0,my1);    
736     itstpcP->Draw("alp");
737     gPad->Clear();
738     itstpcM->Draw("apl");
739     itstpcP->Draw("p");
740     itstpcB->Draw("p");
741     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
742     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
743     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
744     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
745     legend->Draw();    
746     //picArray->AddLast(pad);
747   }
748
749   if (itstpcB&&laserA&&itstpcP&&itstpcM){
750     pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
751     SetPadStyle(pad,mx0,mx1,my0,my1);    
752     laserA->Draw("alp");
753     itstpcP->Draw("p");
754     itstpcM->Draw("p");
755     itstpcB->Draw("p");
756     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
757     legend->AddEntry(laserA,"TPC laser");
758     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");   
759     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");   
760     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
761     legend->Draw();
762     //picArray->AddLast(pad);
763   }
764
765   if (itstpcP&&cross){ 
766     pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
767     SetPadStyle(pad,mx0,mx1,my0,my1);    
768     itstpcP->Draw("alp");
769     pad->Clear();
770     cross->Draw("ap");
771     itstpcP->Draw("p");
772     //
773     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
774
775     legend->AddEntry(cross,"TPC cross tracks");
776     legend->AddEntry(itstpcB,"ITS-TPC smooth");
777     legend->Draw();        
778     //picArray->AddLast(pad);
779   }
780   if (itstpcP&&cosmic){ 
781     pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
782     SetPadStyle(pad,mx0,mx1,my0,my1);    
783     itstpcP->Draw("alp");
784     pad->Clear();
785     cosmic->Draw("ap");
786     itstpcP->Draw("p");
787     //
788     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
789
790     legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
791     legend->AddEntry(itstpcB,"ITS-TPC smooth");
792     legend->Draw();        
793     //picArray->AddLast(pad);
794   }
795 }
796
797
798
799
800 void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber,  TString  pocdbStorage){
801   //
802   // Update OCDB gain
803   //
804   if (pocdbStorage.Length()==0) pocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
805
806   //
807   // 1. Read gain values
808   //
809   ReadGainGlobal(fileName);
810
811   //
812   // 2. Extract calibration values
813   //
814   AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
815   AnalyzeAttachment(startRunNumber,endRunNumber);
816   AnalyzePadRegionGain();
817   AnalyzeGainMultiplicity();
818   //
819   // 3. Make control plots
820   //
821   MakeQAPlot(1.43);  
822
823   //
824   // 4. validate OCDB entries
825   //
826   if(fSwitchOnValidation==kTRUE && ValidateTimeGain()==kFALSE) { 
827     Printf("TPC time gain OCDB parameters out of range!");
828     return;
829   }
830
831   //
832   // 5. Update OCDB
833   //
834   UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage.Data());
835 }
836
837 void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
838   //
839   // read calibration entries from file
840   // 
841   TFile fcalib(fileName);
842   TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
843   if (array){
844     fGainMIP    = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
845     fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
846     fGainMult   = ( AliTPCcalibGainMult *)array->FindObject("calibGainMult");
847   }else{
848     fGainMIP    = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
849     fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
850     fGainMult   = ( AliTPCcalibGainMult *)fcalib.Get("calibGainMult");
851   }
852   if (!fGainMult){
853     TFile fcalibMult("TPCMultObjects.root");
854     fGainMult   = ( AliTPCcalibGainMult *)fcalibMult.Get("calibGainMult");
855   }
856   TH1 * hisT=0;
857   Int_t firstBinA =0, lastBinA=0;
858
859   if (fGainCosmic){ 
860     hisT= fGainCosmic->GetHistGainTime()->Projection(1);
861     firstBinA = hisT->FindFirstBinAbove(2);
862     lastBinA  = hisT->FindLastBinAbove(2);    
863     fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
864     delete hisT;
865   }
866
867   if (fGainMIP){ 
868     hisT= fGainMIP->GetHistGainTime()->Projection(1);
869     firstBinA = hisT->FindFirstBinAbove(2);
870     lastBinA  = hisT->FindLastBinAbove(2);    
871     fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
872     delete hisT;
873   }
874
875 }
876
877
878
879 Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit,  Float_t FPtoMIPratio){
880   //
881   // Analyze gain - produce the calibration graphs
882   //
883
884   // 1.) try to create MIP spline
885   if (fGainMIP) 
886   {
887     fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
888     fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
889     fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
890     //
891     fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
892     if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
893     if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
894     if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
895     fGainArray->AddAt(fFitMIP,0);
896   } 
897
898   // 2.) try to create Cosmic spline
899   if (fGainCosmic)
900   {
901     fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
902     fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100);    // only Fermi-Plateau muons
903     //
904     fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
905     if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
906     //
907     if (fGraphCosmic) {
908       for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
909         fGraphCosmic->GetY()[i] /= FPtoMIPratio;
910         fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
911       }
912     }
913     //
914     if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
915     if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
916     fGainArray->AddAt(fFitCosmic,1);
917   }
918   // with naming convention and backward compatibility
919   fGainArray->AddAt(fGraphMIP,2);
920   fGainArray->AddAt(fGraphCosmic,3);
921   cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
922   return kTRUE;
923
924 }
925
926 Bool_t AliTPCPreprocessorOffline::AnalyzeAttachment(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesFit) {
927   //
928   // determine slope as a function of mean driftlength
929   //
930   if(!fGainMIP) return kFALSE;
931
932   fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
933   //
934   fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
935   fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
936   //
937   fGainMIP->GetHistGainTime()->GetAxis(3)->SetRangeUser(125,250);// only full tracking region (driftlength)
938   fGainMIP->GetHistGainTime()->GetAxis(0)->SetRangeUser(1.5,3.5);// only full tracking region (driftlength)
939   //
940   TH3D * hist = fGainMIP->GetHistGainTime()->Projection(1, 0, 3);
941   //
942   Double_t *xvec = new Double_t[hist->GetNbinsX()];
943   Double_t *yvec = new Double_t[hist->GetNbinsX()];
944   Double_t *xerr = new Double_t[hist->GetNbinsX()];
945   Double_t *yerr = new Double_t[hist->GetNbinsX()];
946   Int_t counter  = 0;
947   //
948   for(Int_t i=1; i < hist->GetNbinsX(); i++) {
949     Int_t nsum=0;
950     Int_t imin   =  i;
951     Int_t imax   =  i;    
952     for (Int_t idelta=0; idelta<5; idelta++){
953       //
954       imin   =  TMath::Max(i-idelta,1);
955       imax   =  TMath::Min(i+idelta,hist->GetNbinsX());
956       nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
957       //if (nsum==0) break;
958       if (nsum>minEntriesFit) break;
959     }
960     if (nsum<minEntriesFit) continue;
961     //
962     fGainMIP->GetHistGainTime()->GetAxis(1)->SetRangeUser(hist->GetXaxis()->GetBinCenter(imin-1),hist->GetXaxis()->GetBinCenter(imax+1)); // define time range
963     TH2D * histZdep = fGainMIP->GetHistGainTime()->Projection(0,3);
964     TObjArray arr;
965     histZdep->FitSlicesY(0,0,-1,0,"QNR",&arr);
966     TH1D * driftDep = (TH1D*)arr.At(1);
967     delete histZdep;
968     //TGraphErrors * driftDep = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,3,100,1,0.,1);
969     /*if (driftDep->GetN() < 4) {
970       delete driftDep;
971       continue;
972       }*/
973     //
974     //TObjArray arr;
975     //
976     TF1 pol1("polynom1","pol1",125,240);
977     //driftDep->Fit(&pol1,"QNRROB=0.8");
978     driftDep->Fit(&pol1,"QNR");
979     xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin-1)+hist->GetXaxis()->GetBinCenter(imax+1));
980     yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
981     xerr[counter] = hist->GetXaxis()->GetBinCenter(imax+1)-hist->GetXaxis()->GetBinCenter(imin-1);
982     yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
983     counter++;
984     //
985     //delete driftDep;
986   }
987   //
988   fGraphAttachmentMIP = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
989   if (fGraphAttachmentMIP) fGraphAttachmentMIP->SetName("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");// set proper names according to naming convention
990   fGainArray->AddLast(fGraphAttachmentMIP);
991   //
992   delete [] xvec;
993   delete [] yvec;
994   delete [] xerr;
995   delete [] yerr;
996   delete hist;
997   //
998   if (counter < 1) return kFALSE;
999   return kTRUE;
1000
1001 }
1002
1003
1004 Bool_t AliTPCPreprocessorOffline::AnalyzePadRegionGain(){
1005   //
1006   // Analyze gain for different pad regions - produce the calibration graphs 0,1,2
1007   //
1008   if (fGainMult) 
1009   {
1010     TH2D * histQmax = (TH2D*) fGainMult->GetHistPadEqual()->Projection(0,2);
1011     TH2D * histQtot = (TH2D*) fGainMult->GetHistPadEqual()->Projection(1,2);
1012     //
1013     TObjArray arr;
1014     histQmax->FitSlicesY(0,0,-1,0,"QNR",&arr);
1015     Double_t xMax[3] = {0,1,2};
1016     Double_t yMax[3]    = {((TH1D*)arr.At(1))->GetBinContent(1),
1017                            ((TH1D*)arr.At(1))->GetBinContent(2),
1018                            ((TH1D*)arr.At(1))->GetBinContent(3)};
1019     Double_t yMaxErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1020                            ((TH1D*)arr.At(1))->GetBinError(2),
1021                            ((TH1D*)arr.At(1))->GetBinError(3)};
1022     TGraphErrors * fitPadRegionQmax = new TGraphErrors(3, xMax, yMax, 0, yMaxErr);
1023     //
1024     histQtot->FitSlicesY(0,0,-1,0,"QNR",&arr);
1025     Double_t xTot[3] = {0,1,2};
1026     Double_t yTot[3]    = {((TH1D*)arr.At(1))->GetBinContent(1),
1027                            ((TH1D*)arr.At(1))->GetBinContent(2),
1028                            ((TH1D*)arr.At(1))->GetBinContent(3)};
1029     Double_t yTotErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1030                            ((TH1D*)arr.At(1))->GetBinError(2),
1031                            ((TH1D*)arr.At(1))->GetBinError(3)};
1032     TGraphErrors * fitPadRegionQtot = new TGraphErrors(3, xTot, yTot, 0, yTotErr);
1033     //
1034     fitPadRegionQtot->SetName("TGRAPHERRORS_MEANQTOT_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1035     fitPadRegionQmax->SetName("TGRAPHERRORS_MEANQMAX_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1036     //
1037     fGainArray->AddLast(fitPadRegionQtot);
1038     fGainArray->AddLast(fitPadRegionQmax);
1039     return kTRUE;
1040   } 
1041   return kFALSE;
1042
1043 }
1044
1045
1046 Bool_t AliTPCPreprocessorOffline::AnalyzeGainMultiplicity() {
1047   //
1048   // Analyze gain as a function of multiplicity and produce calibration graphs
1049   //
1050   if (!fGainMult) return kFALSE;
1051   fGainMult->GetHistGainMult()->GetAxis(3)->SetRangeUser(3,3);
1052   TH2D * histMultMax = fGainMult->GetHistGainMult()->Projection(0,4);
1053   TH2D * histMultTot = fGainMult->GetHistGainMult()->Projection(1,4);
1054   histMultMax->RebinX(4);
1055   histMultTot->RebinX(4);
1056   //
1057   TObjArray arrMax;
1058   TObjArray arrTot;
1059   histMultMax->FitSlicesY(0,0,-1,0,"QNR",&arrMax);
1060   histMultTot->FitSlicesY(0,0,-1,0,"QNR",&arrTot);
1061   //
1062   TH1D * meanMax = (TH1D*)arrMax.At(1);
1063   TH1D * meanTot = (TH1D*)arrTot.At(1);
1064   Float_t meanMult = histMultMax->GetMean();
1065   if(meanMax->GetBinContent(meanMax->FindBin(meanMult))) {
1066     meanMax->Scale(1./meanMax->GetBinContent(meanMax->FindBin(meanMult)));
1067   }
1068   else {
1069    return kFALSE;
1070   }
1071   if(meanTot->GetBinContent(meanTot->FindBin(meanMult))) {
1072     meanTot->Scale(1./meanTot->GetBinContent(meanTot->FindBin(meanMult)));
1073   }
1074   else {
1075    return kFALSE;
1076   }
1077   Float_t xMultMax[50];
1078   Float_t yMultMax[50];
1079   Float_t yMultErrMax[50];
1080   Float_t xMultTot[50];
1081   Float_t yMultTot[50];
1082   Float_t yMultErrTot[50];
1083   //
1084   Int_t nCountMax = 0;
1085   for(Int_t iBin = 1; iBin < meanMax->GetXaxis()->GetNbins(); iBin++) {
1086     Float_t yValMax = meanMax->GetBinContent(iBin);
1087     if (yValMax < 0.7) continue;
1088     if (yValMax > 1.3) continue;
1089     if (meanMax->GetBinError(iBin)/yValMax > 0.01) continue;
1090     xMultMax[nCountMax] = meanMax->GetXaxis()->GetBinCenter(iBin);
1091     yMultMax[nCountMax] = yValMax;
1092     yMultErrMax[nCountMax] = meanMax->GetBinError(iBin);
1093     nCountMax++;
1094   }
1095   //
1096   if (nCountMax < 10) return kFALSE;
1097   TGraphErrors * fitMultMax = new TGraphErrors(nCountMax, xMultMax, yMultMax, 0, yMultErrMax);
1098   fitMultMax->SetName("TGRAPHERRORS_MEANQMAX_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1099   //
1100   Int_t nCountTot = 0;
1101   for(Int_t iBin = 1; iBin < meanTot->GetXaxis()->GetNbins(); iBin++) {
1102     Float_t yValTot = meanTot->GetBinContent(iBin);
1103     if (yValTot < 0.1) continue;
1104     if (meanTot->GetBinError(iBin)/yValTot > 0.1) continue;
1105     xMultTot[nCountTot] = meanTot->GetXaxis()->GetBinCenter(iBin);
1106     yMultTot[nCountTot] = yValTot;
1107     yMultErrTot[nCountTot] = meanTot->GetBinError(iBin);
1108     nCountTot++;
1109   }
1110   //
1111   if (nCountTot < 10) return kFALSE;
1112   TGraphErrors *  fitMultTot = new TGraphErrors(nCountTot, xMultTot, yMultTot, 0, yMultErrTot);
1113   fitMultTot->SetName("TGRAPHERRORS_MEANQTOT_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1114   //
1115   fGainArray->AddLast(fitMultMax);
1116   fGainArray->AddLast(fitMultTot);
1117   //
1118   return kTRUE;
1119
1120 }
1121
1122
1123 void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
1124   //
1125   // Update OCDB entry
1126   //
1127   AliCDBMetaData *metaData= new AliCDBMetaData();
1128   metaData->SetObjectClassName("TObjArray");
1129   metaData->SetResponsible("Alexander Kalweit");
1130   metaData->SetBeamPeriod(1);
1131   metaData->SetAliRootVersion("05-24-00"); //root version
1132   metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
1133   AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
1134   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
1135   gStorage->Put(fGainArray, id1, metaData);    
1136 }
1137
1138 void AliTPCPreprocessorOffline::MakeQAPlot(Float_t  FPtoMIPratio) {
1139   //
1140   // Make QA plot to visualize results
1141   //
1142   //
1143   //
1144   if (fGraphCosmic) {
1145     TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
1146     canvasCosmic->cd();
1147     TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
1148     gainHistoCosmic->SetDirectory(0);
1149     gainHistoCosmic->SetName("GainHistoCosmic");
1150     gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
1151     gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1152     gainHistoCosmic->Draw("colz");
1153     fGraphCosmic->SetMarkerStyle(25);
1154     fGraphCosmic->Draw("lp");
1155     fGraphCosmic->SetMarkerStyle(25);
1156     TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
1157     if (grfFitCosmic) {
1158       for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
1159         grfFitCosmic->GetY()[i] *= FPtoMIPratio;        
1160       }
1161       for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
1162         fGraphCosmic->GetY()[i] *= FPtoMIPratio;        
1163       }
1164     }
1165     fGraphCosmic->Draw("lp");
1166     grfFitCosmic->SetLineColor(2);
1167     grfFitCosmic->Draw("lu");
1168     fGainArray->AddLast(gainHistoCosmic);
1169     //fGainArray->AddLast(canvasCosmic->Clone());
1170     delete canvasCosmic;    
1171   }
1172   if (fFitMIP) {
1173     TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
1174     canvasMIP->cd();
1175     TH2D * gainHistoMIP    = fGainMIP->GetHistGainTime()->Projection(0,1);
1176     gainHistoMIP->SetName("GainHistoCosmic");
1177     gainHistoMIP->SetDirectory(0);
1178     gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
1179     gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1180     gainHistoMIP->Draw("colz");
1181     fGraphMIP->SetMarkerStyle(25);
1182     fGraphMIP->Draw("lp");
1183     TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
1184     grfFitMIP->SetLineColor(2);
1185     grfFitMIP->Draw("lu");    
1186     fGainArray->AddLast(gainHistoMIP);
1187     //fGainArray->AddLast(canvasMIP->Clone());
1188     delete canvasMIP;    
1189   }  
1190 }
1191
1192 void AliTPCPreprocessorOffline::MakeFitTime(){
1193   //
1194   // mak aligment fit - store results in the file
1195   //
1196   const Int_t kMinEntries=1000;
1197   MakeChainTime();
1198   MakePrimitivesTime();
1199   if (!fAlignTree) return;
1200   if (fAlignTree->GetEntries()<kMinEntries) return;
1201   fAlignTree->SetAlias("ptype","type");
1202   fAlignTree->SetAlias("hasITS","(1+0)");
1203   fAlignTree->SetAlias("dITS","1-2*(refX<40)");
1204   fAlignTree->SetAlias("isITS","refX>10");
1205   fAlignTree->SetAlias("isVertex","refX<10");
1206   // 
1207   Int_t  npointsMax=30000000;
1208   TStatToolkit toolkit;
1209   Double_t chi2=0;
1210   Int_t    npoints=0;
1211   TVectorD param;
1212   TMatrixD covar;
1213
1214   TString fstringFast="";
1215   fstringFast+="FExBTwistX++";
1216   fstringFast+="FExBTwistY++";
1217   fstringFast+="FAlignRot0D++";
1218   fstringFast+="FAlignTrans0D++";
1219   fstringFast+="FAlignTrans1D++";
1220   //
1221   fstringFast+="hasITS*FAlignTrans0++";
1222   fstringFast+="hasITS*FAlignTrans1++";
1223   fstringFast+="hasITS*FAlignRot0++";
1224   fstringFast+="hasITS*FAlignRot1++";
1225   fstringFast+="hasITS*FAlignRot2++";
1226   //
1227   fstringFast+="dITS*FAlignTrans0++";
1228   fstringFast+="dITS*FAlignTrans1++";
1229   fstringFast+="dITS*FAlignRot0++";
1230   fstringFast+="dITS*FAlignRot1++";
1231   fstringFast+="dITS*FAlignRot2++";
1232   
1233   TCut cutFit="entries>10&&abs(mean)>0.00001";
1234   fAlignTree->SetAlias("err","rms");
1235
1236   TString *strDeltaITS = TStatToolkit::FitPlaneConstrain(fAlignTree,"mean:err", fstringFast.Data(),cutFit, chi2,npoints,param,covar,-1,0, npointsMax, 1);
1237   strDeltaITS->Tokenize("++")->Print();
1238   fAlignTree->SetAlias("fitYFast",strDeltaITS->Data());
1239   // 
1240   TVectorD paramC= param;
1241   TMatrixD covarC= covar;
1242   TStatToolkit::Constrain1D(fstringFast,"Trans0D",paramC,covarC,0, 0.1);
1243   TStatToolkit::Constrain1D(fstringFast,"Trans1D",paramC,covarC,0, 0.1);
1244   TStatToolkit::Constrain1D(fstringFast,"TwistX",paramC,covarC,0, 0.1);
1245   TStatToolkit::Constrain1D(fstringFast,"TwistY",paramC,covarC,0, 0.1);
1246   TString strFitConst=TStatToolkit::MakeFitString(fstringFast, paramC,covar);
1247   fAlignTree->SetAlias("fitYFastC",strFitConst.Data());
1248   CreateAlignTime(fstringFast,paramC);
1249
1250
1251 }
1252
1253
1254 void AliTPCPreprocessorOffline::MakeChainTime(){
1255   //
1256   TFile f("TPCTimeObjects.root");
1257   //  const char *cdtype[7]={"ITS","TRD","Vertex","TOF","TPC","TPC0","TPC1"};
1258   //const char *cptype[5]={"dy","dz","dsnp","dtheta","d1pt"}; 
1259   const char * hname[5]={"dy","dz","dsnp","dtheta","d1pt"};
1260   Int_t run=0;
1261   AliTPCcalibTime  *calibTime= (AliTPCcalibTime*) f.Get("calibTime");
1262   if (!calibTime) return;
1263   TTreeSRedirector *pcstream = new TTreeSRedirector("meanITSVertex.root");
1264   Int_t ihis=0;
1265   THnSparse *his = calibTime->GetResHistoTPCITS(ihis);
1266   if (his){
1267     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1268     his->GetAxis(2)->SetRange(0,1000000);
1269     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1270     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,5);
1271   }
1272   ihis=2;
1273   his = calibTime->GetResHistoTPCITS(ihis);
1274   if (his){
1275     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1276     his->GetAxis(2)->SetRange(0,1000000);
1277     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1278     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,5);
1279   }
1280   ihis=0;
1281   his = calibTime->GetResHistoTPCvertex(ihis);
1282   if (his){
1283     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1284     his->GetAxis(2)->SetRange(0,1000000);
1285     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1286     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,5);
1287   }
1288   ihis=2;
1289   his = calibTime->GetResHistoTPCvertex(ihis);
1290   if (his){
1291     his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1292     his->GetAxis(2)->SetRange(0,1000000);
1293     his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1294     AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,5);
1295
1296   }
1297   delete pcstream;
1298 }
1299
1300
1301 Double_t AliTPCPreprocessorOffline::EvalAt(Double_t phi, Double_t refX, Double_t theta, Int_t corr, Int_t ptype){
1302   //
1303   //
1304   //
1305   Double_t sector = 9*phi/TMath::Pi();
1306   if (sector<0) sector+=18;
1307   Double_t y85=AliTPCCorrection::GetCorrSector(sector,85,theta,1,corr);
1308   Double_t y245=AliTPCCorrection::GetCorrSector(sector,245,theta,1,corr);
1309   if (ptype==0) return y85+(y245-y85)*(refX-85.)/(245.-85.);
1310   if (ptype==2) return (y245-y85)/(245.-85.);
1311   return 0;
1312 }
1313
1314
1315
1316 void AliTPCPreprocessorOffline::MakePrimitivesTime(){
1317   //
1318   // Create primitive transformation to fit
1319   //
1320   fAlignTree=new TChain("fit","fit");
1321   fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdy");
1322   fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdsnp");
1323   fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdy");
1324   fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdsnp");
1325   // 
1326   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1327   Double_t bzField=AliTrackerBase::GetBz(); 
1328   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
1329   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
1330   Double_t wtP = -10.0 * (bzField) * vdrift /  ezField ; 
1331   AliTPCExBTwist *fitExBTwistX= new  AliTPCExBTwist;
1332   AliTPCExBTwist *fitExBTwistY= new  AliTPCExBTwist;
1333   AliTPCCalibGlobalMisalignment *trans0   =new  AliTPCCalibGlobalMisalignment;
1334   AliTPCCalibGlobalMisalignment *trans1   =new  AliTPCCalibGlobalMisalignment;
1335   AliTPCCalibGlobalMisalignment *trans0D  =new  AliTPCCalibGlobalMisalignment;
1336   AliTPCCalibGlobalMisalignment *trans1D  =new  AliTPCCalibGlobalMisalignment;
1337   AliTPCCalibGlobalMisalignment *rot0     =new  AliTPCCalibGlobalMisalignment;
1338   AliTPCCalibGlobalMisalignment *rot1     =new  AliTPCCalibGlobalMisalignment;
1339   AliTPCCalibGlobalMisalignment *rot2     =new  AliTPCCalibGlobalMisalignment;
1340   AliTPCCalibGlobalMisalignment *rot3     =new  AliTPCCalibGlobalMisalignment;
1341   //
1342   //
1343   fitExBTwistX->SetXTwist(0.001);
1344   fitExBTwistX->SetOmegaTauT1T2(wtP,1,1);  
1345   //
1346   fitExBTwistY->SetYTwist(0.001);
1347   fitExBTwistY->SetOmegaTauT1T2(wtP,1,1);  
1348   //
1349   TGeoHMatrix *matrixRot = new TGeoHMatrix; 
1350   TGeoHMatrix *matrixX = new TGeoHMatrix; 
1351   TGeoHMatrix *matrixY = new TGeoHMatrix; 
1352   matrixX->SetDx(0.1);
1353   matrixY->SetDy(0.1);
1354   Double_t rotAngles0[9]={0};
1355   Double_t rotAngles1[9]={0};
1356   Double_t rotAngles2[9]={0};
1357   //
1358   Double_t rotAngles3[9]={0};
1359
1360   rotAngles0[0]=1; rotAngles0[4]=1; rotAngles0[8]=1;
1361   rotAngles1[0]=1; rotAngles1[4]=1; rotAngles1[8]=1;
1362   rotAngles2[0]=1; rotAngles2[4]=1; rotAngles2[8]=1;
1363   rotAngles3[0]=1; rotAngles3[4]=1; rotAngles3[8]=1;
1364
1365   rotAngles0[1]=-0.001;rotAngles0[3]=0.001;
1366   rotAngles1[5]=-0.001;rotAngles1[7]=0.001;
1367   rotAngles2[2]=0.001;rotAngles2[6]=-0.001;
1368   rotAngles3[1]=0.001;rotAngles3[3]=-0.001;
1369   matrixRot->SetRotation(rotAngles0);
1370   rot0->SetAlignGlobal(matrixRot);
1371   matrixRot->SetRotation(rotAngles1);
1372   rot1->SetAlignGlobal(matrixRot);
1373   matrixRot->SetRotation(rotAngles2);
1374   rot2->SetAlignGlobal(matrixRot); 
1375   matrixRot->SetRotation(rotAngles3);
1376   rot3->SetAlignGlobalDelta(matrixRot); 
1377   //
1378   trans0->SetAlignGlobal(matrixX);
1379   trans1->SetAlignGlobal(matrixY);
1380   trans0D->SetAlignGlobalDelta(matrixX);
1381   trans1D->SetAlignGlobalDelta(matrixY);
1382   fitExBTwistX->Init();
1383   fitExBTwistY->Init();
1384   //
1385   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistX->Clone()),100);
1386   fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistY->Clone()),101);
1387   //
1388   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot0->Clone()),102);
1389   fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(rot1->Clone()),103);
1390   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot2->Clone()),104);
1391   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot3->Clone()),105);
1392
1393   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0->Clone()),106);
1394   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1->Clone()),107);
1395   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0D->Clone()),108);
1396   fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1D->Clone()),109);
1397   //
1398   fAlignTree->SetAlias("FExBTwistX", "AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,100,ptype)+0");
1399   fAlignTree->SetAlias("FExBTwistY","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,101,ptype)+0");
1400   fAlignTree->SetAlias("FAlignRot0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,102,ptype)+0");
1401   fAlignTree->SetAlias("FAlignRot0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,105,ptype)+0");
1402   fAlignTree->SetAlias("FAlignRot1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,103,ptype)+0");
1403   fAlignTree->SetAlias("FAlignRot2","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,104,ptype)+0");
1404   fAlignTree->SetAlias("FAlignTrans0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,106,ptype)+0");
1405   fAlignTree->SetAlias("FAlignTrans1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,107,ptype)+0");
1406   fAlignTree->SetAlias("FAlignTrans0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,108,ptype)+0");
1407   fAlignTree->SetAlias("FAlignTrans1D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,109,ptype)+0");
1408   //
1409   // test fast function
1410   //
1411 //   fAlignTree->Draw("FExBTwistX:ExBTwistX","isITS&&ptype==0&&abs(snp)<0.05","");
1412 //   fAlignTree->Draw("FExBTwistY:ExBTwistY","isITS&&ptype==0&&abs(snp)<0.05","");
1413 //   fAlignTree->Draw("FAlignRot0:alignRot0","isITS&&ptype==0&&abs(snp)<0.05","");
1414 //   fAlignTree->Draw("FAlignRot1:alignRot1","isITS&&ptype==0&&abs(snp)<0.05","");
1415 //   fAlignTree->Draw("FAlignRot2:alignRot2","isITS&&ptype==0&&abs(snp)<0.05","");
1416 //   //
1417 //   fAlignTree->Draw("FAlignTrans0:alignTrans0","isITS&&ptype==0&&abs(snp)<0.05","");
1418 //   fAlignTree->Draw("FAlignTrans1:alignTrans1","isITS&&ptype==0&&abs(snp)<0.05","");
1419
1420
1421
1422
1423 void AliTPCPreprocessorOffline::CreateAlignTime(TString fstring, TVectorD paramC){
1424   //
1425   //
1426   //
1427   //
1428   TGeoHMatrix *matrixDelta     = new TGeoHMatrix; 
1429   TGeoHMatrix *matrixGlobal    = new TGeoHMatrix; 
1430   Double_t rAngles[9];
1431   Int_t index=0;
1432   //
1433   index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans0D");
1434   if (index>=0) matrixDelta->SetDx(paramC[index+1]*0.1);
1435   index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans1D");
1436   if (index>=0) matrixDelta->SetDy(paramC[index+1]*0.1);
1437   rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1438   index=TStatToolkit::GetFitIndex(fstring,"FAlignRot0D");
1439   rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1440   rAngles[5]=0; rAngles[7] =0;
1441   rAngles[2]=0; rAngles[6] =0;
1442   matrixDelta->SetRotation(rAngles);
1443   //
1444   //
1445   //
1446   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans0");
1447   if (index>=0) matrixGlobal->SetDx(paramC[index+1]*0.1);
1448   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans1");
1449   if (index>=0) matrixGlobal->SetDy(paramC[index+1]*0.1);
1450   rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1451   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot0");
1452   rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1453   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot1");  
1454   rAngles[5]=-paramC[index+1]*0.001; rAngles[7]=paramC[index+1]*0.001;
1455   index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot2");  
1456   rAngles[2]=paramC[index+1]*0.001; rAngles[6] =-paramC[index+1]*0.001;
1457   matrixGlobal->SetRotation(rAngles);
1458   //
1459   AliTPCCalibGlobalMisalignment *fitAlignTime  =new  AliTPCCalibGlobalMisalignment;
1460   fitAlignTime  =new  AliTPCCalibGlobalMisalignment;
1461   fitAlignTime->SetName("FitAlignTime");
1462   fitAlignTime->SetTitle("FitAlignTime");
1463   fitAlignTime->SetAlignGlobalDelta(matrixDelta);
1464   fitAlignTime->SetAlignGlobal(matrixGlobal);
1465   //
1466   AliTPCExBTwist * fitExBTwist= new  AliTPCExBTwist;
1467   Int_t indexX=TStatToolkit::GetFitIndex(fstring,"ExBTwistX");
1468   Int_t indexY=TStatToolkit::GetFitIndex(fstring,"ExBTwistY");  
1469   fitExBTwist->SetXTwist(0.001*paramC[indexX+1]);  // 1 mrad twist in x
1470   fitExBTwist->SetYTwist(0.001*paramC[indexY+1]);  // 1 mrad twist in x
1471   fitExBTwist->SetName("FitExBTwistTime");
1472   fitExBTwist->SetTitle("FitExBTwistTime"); 
1473   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1474   Double_t bzField=AliTrackerBase::GetBz();
1475   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
1476
1477   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
1478   Double_t wt = -10.0 * (bzField) * vdrift /  ezField ; 
1479   //
1480   fitExBTwist->SetOmegaTauT1T2(wt,1,1);  
1481   fitExBTwist->Init();  
1482
1483   AliTPCComposedCorrection *corrTime =  new AliTPCComposedCorrection;
1484   TObjArray *arr = new TObjArray;
1485   corrTime->SetCorrections(arr);
1486   
1487   corrTime->GetCorrections()->Add(fitExBTwist);
1488   corrTime->GetCorrections()->Add(fitAlignTime);
1489   corrTime->SetName("FitCorrectionTime");
1490   corrTime->SetTitle("FitCorrectionTime");
1491
1492   fitExBTwist->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwist->Clone()),1001);
1493   fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(fitAlignTime->Clone()),1002);
1494   fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(corrTime->Clone()),1003);
1495   
1496   
1497   fAlignTree->SetAlias("ExBTwistTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1001,ptype)+0");
1498   fAlignTree->SetAlias("AlignTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1002,ptype)+0");
1499   fAlignTree->SetAlias("FitCorrectionTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1003,ptype)+0");
1500
1501
1502   TFile *f = new TFile("fitITSVertex.root","update");
1503   corrTime->Write("FitCorrectionTime");
1504   f->Close();
1505 }
1506
1507
1508
1509