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