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