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