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