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