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