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