1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 Responsible: marian.ivanov@cern.ch
20 Code to analyze the TPC calibration and to produce OCDB entries
24 gSystem->Load("libANALYSIS");
25 gSystem->Load("libTPCcalib");
27 AliTPCPreprocessorOffline proces;
28 TString ocdbPath="local:////"
29 ocdbPath+=gSystem->GetFromPipe("pwd");
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
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")
48 #include "Riostream.h"
51 #include "TGraphErrors.h"
52 #include "AliExternalTrackParam.h"
55 #include "TMultiGraph.h"
57 #include "THnSparse.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"
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"
86 #include "AliTrackerBase.h"
87 #include "AliTracker.h"
88 #include "AliTPCPreprocessorOffline.h"
89 #include "AliTPCCorrectionFit.h"
93 ClassImp(AliTPCPreprocessorOffline)
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),
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
114 fAlignTree(0), // alignment tree
115 fSwitchOnValidation(kFALSE), // flag to switch on validation of OCDB parameters
121 // default constructor
125 AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
134 void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const timeDrift){
136 // find the fist and last run
138 TObjArray *hisArray =timeDrift->GetHistoDrift();
139 {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
140 THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
141 if (!addHist) continue;
142 if (addHist->GetEntries()<fMinEntries) continue;
143 TH1D* histo =addHist->Projection(3);
144 TH1D* histoTime=addHist->Projection(0);
145 printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
148 startRun=histo->FindFirstBinAbove(0);
149 endRun =histo->FindLastBinAbove(0);
151 startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun);
152 endRun =TMath::Max(histo->FindLastBinAbove(0),endRun);
155 startTime=histoTime->FindFirstBinAbove(0);
156 endTime =histoTime->FindLastBinAbove(0);
158 startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime);
159 endTime =TMath::Max(histoTime->FindLastBinAbove(0),endTime);
164 if (startRun<0) startRun=0;
165 if (endRun<0) endRun=100000000;
166 printf("Run range :\t%d-%d\n", startRun, endRun);
167 printf("Time range :\t%d-%d\n", startTime, endTime);
173 void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, TString pocdbStorage){
175 // make calibration of the drift velocity
177 // file - the location of input file
178 // ustartRun, uendrun - run validity period
179 // pocdbStorage - path to hte OCDB storage
180 // - if empty - local storage 'pwd' uesed
181 if (pocdbStorage.Length()>0) ocdbStorage=pocdbStorage;
183 ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
185 // 1. Initialization and run range setting
187 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
189 fTimeDrift = (AliTPCcalibTime *)array->FindObject("calibTime");
191 fTimeDrift = (AliTPCcalibTime*)fcalib.Get("calibTime");
193 if(!fTimeDrift) return;
197 TObjArray *hisArray =fTimeDrift->GetHistoDrift();
198 GetRunRange(fTimeDrift);
199 for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
200 THnSparse* addHist=(THnSparse*)hisArray->At(i);
201 if (!addHist) continue;
202 if (startTime<endTime) addHist->GetAxis(0)->SetRange(startTime-1,endTime+1);
203 if (startRun<endRun) addHist->GetAxis(3)->SetRange(startRun-1,endRun+1);
207 // 2. extraction of the information
209 fVdriftArray = new TObjArray();
210 AddAlignmentGraphs(fVdriftArray,fTimeDrift);
211 AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
212 AddLaserGraphs(fVdriftArray,fTimeDrift);
214 // 3. Append QA plots
216 MakeDefaultPlots(fVdriftArray,fVdriftArray);
219 // 4. validate OCDB entries
221 if(fSwitchOnValidation==kTRUE && ValidateTimeDrift()==kFALSE) {
222 Printf("TPC time drift OCDB parameters out of range!");
229 TFile * ftime= TFile::Open("fitITSVertex.root");
231 TObject * alignmentTime=ftime->Get("FitCorrectionTime");
232 if (alignmentTime) fVdriftArray->AddLast(alignmentTime);
239 UpdateOCDBDrift(ustartRun,uendRun,ocdbStorage);
242 void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun, const char* storagePath ){
246 AliCDBMetaData *metaData= new AliCDBMetaData();
247 metaData->SetObjectClassName("TObjArray");
248 metaData->SetResponsible("Marian Ivanov");
249 metaData->SetBeamPeriod(1);
250 metaData->SetAliRootVersion("05-25-01"); //root version
251 metaData->SetComment("Calibration of the time dependence of the drift velocity");
253 id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
254 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
255 gStorage->Put(fVdriftArray, (*id1), metaData);
258 Bool_t AliTPCPreprocessorOffline::ValidateTimeGain()
261 // Validate time gain corrections
263 Printf("ValidateTimeGain..." );
264 Float_t minGain = fMinGain;
265 Float_t maxGain = fMaxGain;
267 TGraphErrors *gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
269 gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
270 if (!gr) return kFALSE;
271 Printf("Assuming given run is a cosmic run. Using gain calibration from Fermi-plateau muons.");
273 if(gr->GetN()<1) return kFALSE;
275 // check whether gain in the range
276 for(Int_t iPoint=0; iPoint<gr->GetN(); iPoint++)
278 if(gr->GetY()[iPoint] < minGain || gr->GetY()[iPoint] > maxGain)
286 Bool_t AliTPCPreprocessorOffline::ValidateTimeDrift()
289 // Validate time drift velocity corrections
291 Printf("ValidateTimeDrift..." );
293 Float_t maxVDriftCorr = fMaxVdriftCorr;
295 TGraphErrors* gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
296 Printf("ALIGN_ITSB_TPC_DRIFTVD graph = %p",gr);
298 if(!gr) return kFALSE;
300 Printf("ALIGN_ITSB_TPC_DRIFTVD number of points = %d",gr->GetN());
304 // check whether drift velocity corrections in the range
305 for(Int_t iPoint = 0; iPoint<gr->GetN(); iPoint++)
307 Printf("Y value from the graph: %f",TMath::Abs(gr->GetY()[iPoint]));
308 if(TMath::Abs(gr->GetY()[iPoint]) > maxVDriftCorr)
315 void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
317 // update the OCDB entry for the nominal time0
320 // AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
321 AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
322 TGraphErrors *grT = (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
323 Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
324 Double_t deltaT = deltaTcm/param->GetDriftV();
325 paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
328 AliCDBMetaData *metaData= new AliCDBMetaData();
329 metaData->SetObjectClassName("TObjArray");
330 metaData->SetResponsible("Marian Ivanov");
331 metaData->SetBeamPeriod(1);
332 metaData->SetAliRootVersion("05-25-02"); //root version
333 metaData->SetComment("Updated calibration of nominal time 0");
335 id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
336 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
337 gStorage->Put(param, (*id1), metaData);
342 void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
344 // Print the names of the entries in array
346 Int_t entries = array->GetEntries();
347 for (Int_t i=0; i<entries; i++){
348 if (!array->At(i)) continue;
349 printf("%d\t %s\n", i, array->At(i)->GetName());
355 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
357 // 1. filter graph - error cut errSigmaCut
358 // 2. filter graph - medianCutAbs around median
360 // errSigmaCut - cut on error
361 // medianCutAbs - cut on value around median
364 // 1. filter graph - error cut errSigmaCut
366 TGraphErrors *graphF;
367 graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
369 if (!graphF) return 0;
370 graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
372 if (!graph) return 0;
374 // filter graph - kMedianCutAbs around median
376 graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
378 if (!graphF) return 0;
379 graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
381 if (!graph) return 0;
387 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
389 // filter outlyer measurement
390 // Only points around median +- cut filtered
392 if (!graph) return 0;
394 Int_t npoints0 = graph->GetN();
397 Double_t *outx=new Double_t[npoints0];
398 Double_t *outy=new Double_t[npoints0];
399 Double_t *errx=new Double_t[npoints0];
400 Double_t *erry=new Double_t[npoints0];
403 if (npoints0<kMinPoints) {
410 for (Int_t iter=0; iter<3; iter++){
412 for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
413 if (graph->GetY()[ipoint]==0) continue;
414 if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;
415 outx[npoints] = graph->GetX()[ipoint];
416 outy[npoints] = graph->GetY()[ipoint];
417 errx[npoints] = graph->GetErrorX(ipoint);
418 erry[npoints] = graph->GetErrorY(ipoint);
421 if (npoints<=1) break;
422 medianY =TMath::Median(npoints,outy);
423 rmsY =TMath::RMS(npoints,outy);
425 TGraphErrors *graphOut=0;
426 if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry);
435 void AliTPCPreprocessorOffline::AddHistoGraphs( TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
437 // Add graphs corresponding to the alignment
439 const Double_t kErrSigmaCut=5; // error sigma cut - for filtering
440 const Double_t kMedianCutAbs=0.03; // error sigma cut - for filtering
442 TObjArray * array=timeDrift->GetHistoDrift();
444 THnSparse* hist=NULL;
445 // 2.a) cosmics with different triggers
446 for (Int_t i=0; i<array->GetEntriesFast();i++){
447 hist=(THnSparseF*)array->UncheckedAt(i);
449 if (hist->GetEntries()<minEntries) continue;
451 TString name=hist->GetName();
452 Int_t dim[4]={0,1,2,3};
453 THnSparse* newHist=hist->Projection(4,dim);
454 newHist->SetName(name);
455 TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
457 printf("Graph =%s filtered out\n", name.Data());
460 printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
461 Int_t pos=name.Index("_");
462 name=name(pos,name.Capacity()-pos);
463 TString graphName=graph->ClassName();
467 graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
470 graph->SetMarkerStyle(i%8+20);
471 graph->SetMarkerColor(i%7);
472 graph->GetXaxis()->SetTitle("Time");
473 graph->GetYaxis()->SetTitle("v_{dcor}");
474 graph->SetName(graphName);
475 graph->SetTitle(graphName);
476 printf("Graph %d\t=\t%s\n", i, graphName.Data());
477 vdriftArray->Add(graph);
486 void AliTPCPreprocessorOffline::AddAlignmentGraphs( TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
488 // Add graphs corresponding to alignment to the object array
490 TObjArray *arrayITS=0;
491 TObjArray *arrayTOF=0;
492 TObjArray *arrayTRD=0;
493 TMatrixD *mstatITS=0;
494 TMatrixD *mstatTOF=0;
495 TMatrixD *mstatTRD=0;
497 arrayITS=timeDrift->GetAlignITSTPC();
498 arrayTRD=timeDrift->GetAlignTRDTPC();
499 arrayTOF=timeDrift->GetAlignTOFTPC();
501 if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.7,50,fMaxVdriftCorr);
502 if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.7,1000,fMaxVdriftCorr);
503 if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.7,50,fMaxVdriftCorr);
505 TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
506 TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
507 TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
508 TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
509 TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
510 TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
512 TObjArray * arrayTRDP= 0x0;
513 TObjArray * arrayTRDM= 0x0;
514 TObjArray * arrayTRDB= 0x0;
515 arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
516 arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
517 arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
520 Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
521 TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
522 arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
523 arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
524 TString grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
525 "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
526 "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
527 TString grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
528 "DELTAX", "DELTAY", "DELTAZ",
529 "DRIFTVD", "T0", "VDGY"};
532 TVectorD vX(entries);
533 TVectorD vY(entries);
534 TVectorD vEx(entries);
535 TVectorD vEy(entries);
537 for (Int_t iarray=0; iarray<12; iarray++){
538 arr = arrays[iarray];
539 if (arr==0) continue;
540 for (Int_t ipar=0; ipar<9; ipar++){
542 for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
543 AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
544 if (!kalman) continue;
545 vX[counter]=kalman->GetTimeStamp();
546 vY[counter]=(*(kalman->GetState()))[ipar];
547 if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
549 vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
553 TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
555 vEx.GetMatrixArray(),
556 vEy.GetMatrixArray());
557 TString grName=grnames[iarray];
560 graph->SetName(grName.Data());
561 vdriftArray->AddLast(graph);
569 void AliTPCPreprocessorOffline::AddLaserGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
571 // add graphs for laser
573 const Double_t delayL0L1 = 0.071; //this is hack for 1/2 weeks
575 TGraphErrors *grLaser[6]={0,0,0,0,0,0};
576 //hisN = timeDrift->GetHistVdriftLaserA(0);
577 if (timeDrift->GetHistVdriftLaserA(0)){
578 grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
579 grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
580 vdriftArray->AddLast(grLaser[0]);
582 if (timeDrift->GetHistVdriftLaserA(1)){
583 grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
584 grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
585 vdriftArray->AddLast(grLaser[1]);
587 if (timeDrift->GetHistVdriftLaserA(2)){
588 grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
589 grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
590 vdriftArray->AddLast(grLaser[2]);
592 if (timeDrift->GetHistVdriftLaserC(0)){
593 grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
594 grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
595 vdriftArray->AddLast(grLaser[3]);
597 if (timeDrift->GetHistVdriftLaserC(1)){
598 grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
599 grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
600 vdriftArray->AddLast(grLaser[4]);
602 if (timeDrift->GetHistVdriftLaserC(2)){
603 grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
604 grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
605 vdriftArray->AddLast(grLaser[5]);
607 for (Int_t i=0; i<6;i++){
609 SetDefaultGraphDrift(grLaser[i], 1,(i+20));
610 grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
616 TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
618 // Make graph with mean values and rms
620 hisN->GetAxis(itime)->SetRange(0,100000000);
621 hisN->GetAxis(ival)->SetRange(0,100000000);
622 TH1 * hisT = hisN->Projection(itime);
623 TH1 * hisV = hisN->Projection(ival);
625 Int_t firstBinA = hisT->FindFirstBinAbove(2);
626 Int_t lastBinA = hisT->FindLastBinAbove(2);
627 Int_t firstBinV = hisV->FindFirstBinAbove(0);
628 Int_t lastBinV = hisV->FindLastBinAbove(0);
629 hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
630 hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
632 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
633 Double_t cont = hisT->GetBinContent(ibin);
634 if (cont<minEntries) continue;
637 TVectorD vecTime(entries);
638 TVectorD vecMean0(entries);
639 TVectorD vecRMS0(entries);
640 TVectorD vecMean1(entries);
641 TVectorD vecRMS1(entries);
643 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
644 Double_t cont = hisT->GetBinContent(ibin);
645 if (cont<minEntries) continue;
646 //hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
647 Int_t minBin = ibin-1;
648 Int_t maxBin = ibin+1;
649 if(minBin <= 0) minBin = 1;
650 if(maxBin >= hisN->GetAxis(itime)->GetNbins()) maxBin = hisN->GetAxis(itime)->GetNbins()-1;
651 hisN->GetAxis(itime)->SetRange(minBin,maxBin);
653 Double_t time = hisT->GetBinCenter(ibin);
654 TH1 * his = hisN->Projection(ival);
655 Double_t nentries0= his->GetBinContent(his->FindBin(0));
656 if (cont-nentries0<minEntries) continue;
658 his->SetBinContent(his->FindBin(0),0);
659 vecTime[entries]=time;
660 vecMean0[entries]=his->GetMean()+offset;
661 vecMean1[entries]=his->GetMeanError();
662 vecRMS0[entries] =his->GetRMS();
663 vecRMS1[entries] =his->GetRMSError();
669 TGraphErrors * graph = new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(), 0, vecMean1.GetMatrixArray());
680 void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
682 // Set default style for QA views
684 graph->GetXaxis()->SetTimeDisplay(kTRUE);
685 graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
686 graph->SetMaximum( 0.025);
687 graph->SetMinimum(-0.025);
688 graph->GetXaxis()->SetTitle("Time");
689 graph->GetYaxis()->SetTitle("v_{dcorr}");
691 graph->GetYaxis()->SetLabelSize(0.03);
692 graph->GetXaxis()->SetLabelSize(0.03);
694 graph->GetXaxis()->SetNdivisions(10,5,0);
695 graph->GetYaxis()->SetNdivisions(10,5,0);
697 graph->GetXaxis()->SetLabelOffset(0.02);
698 graph->GetYaxis()->SetLabelOffset(0.005);
700 graph->GetXaxis()->SetTitleOffset(1.3);
701 graph->GetYaxis()->SetTitleOffset(1.2);
703 graph->SetMarkerColor(color);
704 graph->SetLineColor(color);
705 graph->SetMarkerStyle(style);
708 void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
710 // Set default pad style for QA
713 pad->SetMargin(mx0,mx1,my0,my1);
717 void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray * /*picArray*/){
719 // 0. make a default QA plots
720 // 1. Store them in the array
723 Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
725 TGraphErrors* laserA =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
726 TGraphErrors* laserC =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
727 TGraphErrors* cosmic =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
728 TGraphErrors* cross =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
729 TGraphErrors* itstpcP =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
730 TGraphErrors* itstpcM =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
731 TGraphErrors* itstpcB =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
733 if (laserA) SetDefaultGraphDrift(laserA,2,25);
734 if (laserC) SetDefaultGraphDrift(laserC,4,26);
735 if (cosmic) SetDefaultGraphDrift(cosmic,3,27);
736 if (cross) SetDefaultGraphDrift(cross,4,28);
737 if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
738 if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
739 if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
747 pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
749 SetPadStyle(pad,mx0,mx1,my0,my1);
752 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
753 legend->AddEntry(laserA,"Laser A side");
754 legend->AddEntry(laserC,"Laser C side");
756 //picArray->AddLast(pad);
759 if (itstpcP&&itstpcM&&itstpcB){
760 pad = new TCanvas("ITSTPC","ITSTPC");
761 itstpcP->Draw("alp");
762 SetPadStyle(pad,mx0,mx1,my0,my1);
763 itstpcP->Draw("alp");
765 itstpcM->Draw("apl");
768 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
769 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
770 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
771 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
773 //picArray->AddLast(pad);
776 if (itstpcB&&laserA&&itstpcP&&itstpcM){
777 pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
778 SetPadStyle(pad,mx0,mx1,my0,my1);
783 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
784 legend->AddEntry(laserA,"TPC laser");
785 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
786 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
787 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
789 //picArray->AddLast(pad);
793 pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
794 SetPadStyle(pad,mx0,mx1,my0,my1);
795 itstpcP->Draw("alp");
800 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
802 legend->AddEntry(cross,"TPC cross tracks");
803 legend->AddEntry(itstpcB,"ITS-TPC smooth");
805 //picArray->AddLast(pad);
807 if (itstpcP&&cosmic){
808 pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
809 SetPadStyle(pad,mx0,mx1,my0,my1);
810 itstpcP->Draw("alp");
815 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
817 legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
818 legend->AddEntry(itstpcB,"ITS-TPC smooth");
820 //picArray->AddLast(pad);
827 void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber, TString pocdbStorage){
831 if (pocdbStorage.Length()==0) pocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
834 // 1. Read gain values
836 ReadGainGlobal(fileName);
839 // 2. Extract calibration values
841 AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
842 AnalyzeAttachment(startRunNumber,endRunNumber);
843 AnalyzePadRegionGain();
844 AnalyzeGainMultiplicity();
845 AnalyzeGainChamberByChamber();
847 // 3. Make control plots
852 // 4. validate OCDB entries
854 if(fSwitchOnValidation==kTRUE && ValidateTimeGain()==kFALSE) {
855 Printf("TPC time gain OCDB parameters out of range!");
862 UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage.Data());
865 void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
867 // read calibration entries from file
869 TFile fcalib(fileName);
870 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
872 fGainMIP = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
873 fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
874 fGainMult = ( AliTPCcalibGainMult *)array->FindObject("calibGainMult");
876 fGainMIP = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
877 fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
878 fGainMult = ( AliTPCcalibGainMult *)fcalib.Get("calibGainMult");
881 TFile fcalibMult("TPCMultObjects.root");
882 fGainMult = ( AliTPCcalibGainMult *)fcalibMult.Get("calibGainMult");
885 Int_t firstBinA =0, lastBinA=0;
888 hisT= fGainCosmic->GetHistGainTime()->Projection(1);
889 firstBinA = hisT->FindFirstBinAbove(2);
890 lastBinA = hisT->FindLastBinAbove(2);
891 fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
896 hisT= fGainMIP->GetHistGainTime()->Projection(1);
897 firstBinA = hisT->FindFirstBinAbove(2);
898 lastBinA = hisT->FindLastBinAbove(2);
899 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
907 Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit, Float_t FPtoMIPratio){
909 // Analyze gain - produce the calibration graphs
912 // 1.) try to create MIP spline
915 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
916 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
917 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
919 fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
920 if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
921 if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
922 if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
923 fGainArray->AddAt(fFitMIP,0);
926 // 2.) try to create Cosmic spline
929 fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
930 fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
932 fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
933 if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
936 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
937 fGraphCosmic->GetY()[i] /= FPtoMIPratio;
938 fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
942 if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
943 if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
944 fGainArray->AddAt(fFitCosmic,1);
946 // with naming convention and backward compatibility
947 fGainArray->AddAt(fGraphMIP,2);
948 fGainArray->AddAt(fGraphCosmic,3);
949 cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
954 Bool_t AliTPCPreprocessorOffline::AnalyzeAttachment(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesFit) {
956 // determine slope as a function of mean driftlength
958 if(!fGainMIP) return kFALSE;
960 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
962 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
963 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
965 fGainMIP->GetHistGainTime()->GetAxis(3)->SetRangeUser(125,250);// only full tracking region (driftlength)
966 fGainMIP->GetHistGainTime()->GetAxis(0)->SetRangeUser(1.5,3.5);// only full tracking region (driftlength)
968 TH3D * hist = fGainMIP->GetHistGainTime()->Projection(1, 0, 3);
970 Double_t *xvec = new Double_t[hist->GetNbinsX()];
971 Double_t *yvec = new Double_t[hist->GetNbinsX()];
972 Double_t *xerr = new Double_t[hist->GetNbinsX()];
973 Double_t *yerr = new Double_t[hist->GetNbinsX()];
976 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
980 for (Int_t idelta=0; idelta<5; idelta++){
982 imin = TMath::Max(i-idelta,1);
983 imax = TMath::Min(i+idelta,hist->GetNbinsX());
984 nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
985 //if (nsum==0) break;
986 if (nsum>minEntriesFit) break;
988 if (nsum<minEntriesFit) continue;
990 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRangeUser(hist->GetXaxis()->GetBinCenter(imin-1),hist->GetXaxis()->GetBinCenter(imax+1)); // define time range
991 TH2D * histZdep = fGainMIP->GetHistGainTime()->Projection(0,3);
993 histZdep->FitSlicesY(0,0,-1,0,"QNR",&arr);
994 TH1D * driftDep = (TH1D*)arr.At(1);
996 //TGraphErrors * driftDep = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,3,100,1,0.,1);
997 /*if (driftDep->GetN() < 4) {
1004 TF1 pol1("polynom1","pol1",125,240);
1005 //driftDep->Fit(&pol1,"QNRROB=0.8");
1006 driftDep->Fit(&pol1,"QNR");
1007 xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin-1)+hist->GetXaxis()->GetBinCenter(imax+1));
1008 yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
1009 xerr[counter] = hist->GetXaxis()->GetBinCenter(imax+1)-hist->GetXaxis()->GetBinCenter(imin-1);
1010 yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
1016 fGraphAttachmentMIP = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
1017 if (fGraphAttachmentMIP) fGraphAttachmentMIP->SetName("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");// set proper names according to naming convention
1018 fGainArray->AddLast(fGraphAttachmentMIP);
1026 if (counter < 1) return kFALSE;
1032 Bool_t AliTPCPreprocessorOffline::AnalyzePadRegionGain(){
1034 // Analyze gain for different pad regions - produce the calibration graphs 0,1,2
1038 TH2D * histQmax = (TH2D*) fGainMult->GetHistPadEqual()->Projection(0,2);
1039 TH2D * histQtot = (TH2D*) fGainMult->GetHistPadEqual()->Projection(1,2);
1042 histQmax->FitSlicesY(0,0,-1,0,"QNR",&arr);
1043 Double_t xMax[3] = {0,1,2};
1044 Double_t yMax[3] = {((TH1D*)arr.At(1))->GetBinContent(1),
1045 ((TH1D*)arr.At(1))->GetBinContent(2),
1046 ((TH1D*)arr.At(1))->GetBinContent(3)};
1047 Double_t yMaxErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1048 ((TH1D*)arr.At(1))->GetBinError(2),
1049 ((TH1D*)arr.At(1))->GetBinError(3)};
1050 TGraphErrors * fitPadRegionQmax = new TGraphErrors(3, xMax, yMax, 0, yMaxErr);
1052 histQtot->FitSlicesY(0,0,-1,0,"QNR",&arr);
1053 Double_t xTot[3] = {0,1,2};
1054 Double_t yTot[3] = {((TH1D*)arr.At(1))->GetBinContent(1),
1055 ((TH1D*)arr.At(1))->GetBinContent(2),
1056 ((TH1D*)arr.At(1))->GetBinContent(3)};
1057 Double_t yTotErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1058 ((TH1D*)arr.At(1))->GetBinError(2),
1059 ((TH1D*)arr.At(1))->GetBinError(3)};
1060 TGraphErrors * fitPadRegionQtot = new TGraphErrors(3, xTot, yTot, 0, yTotErr);
1062 fitPadRegionQtot->SetName("TGRAPHERRORS_MEANQTOT_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1063 fitPadRegionQmax->SetName("TGRAPHERRORS_MEANQMAX_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1065 fGainArray->AddLast(fitPadRegionQtot);
1066 fGainArray->AddLast(fitPadRegionQmax);
1074 Bool_t AliTPCPreprocessorOffline::AnalyzeGainMultiplicity() {
1076 // Analyze gain as a function of multiplicity and produce calibration graphs
1078 if (!fGainMult) return kFALSE;
1079 fGainMult->GetHistGainMult()->GetAxis(3)->SetRangeUser(3,3);
1080 TH2D * histMultMax = fGainMult->GetHistGainMult()->Projection(0,4);
1081 TH2D * histMultTot = fGainMult->GetHistGainMult()->Projection(1,4);
1082 histMultMax->RebinX(4);
1083 histMultTot->RebinX(4);
1087 histMultMax->FitSlicesY(0,0,-1,0,"QNR",&arrMax);
1088 histMultTot->FitSlicesY(0,0,-1,0,"QNR",&arrTot);
1090 TH1D * meanMax = (TH1D*)arrMax.At(1);
1091 TH1D * meanTot = (TH1D*)arrTot.At(1);
1092 Float_t meanMult = histMultMax->GetMean();
1093 if(meanMax->GetBinContent(meanMax->FindBin(meanMult))) {
1094 meanMax->Scale(1./meanMax->GetBinContent(meanMax->FindBin(meanMult)));
1099 if(meanTot->GetBinContent(meanTot->FindBin(meanMult))) {
1100 meanTot->Scale(1./meanTot->GetBinContent(meanTot->FindBin(meanMult)));
1105 Float_t xMultMax[50];
1106 Float_t yMultMax[50];
1107 Float_t yMultErrMax[50];
1108 Float_t xMultTot[50];
1109 Float_t yMultTot[50];
1110 Float_t yMultErrTot[50];
1112 Int_t nCountMax = 0;
1113 for(Int_t iBin = 1; iBin < meanMax->GetXaxis()->GetNbins(); iBin++) {
1114 Float_t yValMax = meanMax->GetBinContent(iBin);
1115 if (yValMax < 0.7) continue;
1116 if (yValMax > 1.3) continue;
1117 if (meanMax->GetBinError(iBin)/yValMax > 0.01) continue;
1118 xMultMax[nCountMax] = meanMax->GetXaxis()->GetBinCenter(iBin);
1119 yMultMax[nCountMax] = yValMax;
1120 yMultErrMax[nCountMax] = meanMax->GetBinError(iBin);
1124 if (nCountMax < 10) return kFALSE;
1125 TGraphErrors * fitMultMax = new TGraphErrors(nCountMax, xMultMax, yMultMax, 0, yMultErrMax);
1126 fitMultMax->SetName("TGRAPHERRORS_MEANQMAX_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1128 Int_t nCountTot = 0;
1129 for(Int_t iBin = 1; iBin < meanTot->GetXaxis()->GetNbins(); iBin++) {
1130 Float_t yValTot = meanTot->GetBinContent(iBin);
1131 if (yValTot < 0.7) continue;
1132 if (yValTot > 1.3) continue;
1133 if (meanTot->GetBinError(iBin)/yValTot > 0.1) continue;
1134 xMultTot[nCountTot] = meanTot->GetXaxis()->GetBinCenter(iBin);
1135 yMultTot[nCountTot] = yValTot;
1136 yMultErrTot[nCountTot] = meanTot->GetBinError(iBin);
1140 if (nCountTot < 10) return kFALSE;
1141 TGraphErrors * fitMultTot = new TGraphErrors(nCountTot, xMultTot, yMultTot, 0, yMultErrTot);
1142 fitMultTot->SetName("TGRAPHERRORS_MEANQTOT_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1144 fGainArray->AddLast(fitMultMax);
1145 fGainArray->AddLast(fitMultTot);
1151 Bool_t AliTPCPreprocessorOffline::AnalyzeGainChamberByChamber(){
1153 // get chamber by chamber gain
1155 TGraphErrors *grShort = fGainMult->GetGainPerChamber(0);
1156 TGraphErrors *grMedium = fGainMult->GetGainPerChamber(1);
1157 TGraphErrors *grLong = fGainMult->GetGainPerChamber(2);
1158 if (grShort==0x0 || grMedium==0x0 || grLong==0x0) {
1165 fGainArray->AddLast(grShort);
1166 fGainArray->AddLast(grMedium);
1167 fGainArray->AddLast(grLong);
1172 void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
1174 // Update OCDB entry
1176 AliCDBMetaData *metaData= new AliCDBMetaData();
1177 metaData->SetObjectClassName("TObjArray");
1178 metaData->SetResponsible("Alexander Kalweit");
1179 metaData->SetBeamPeriod(1);
1180 metaData->SetAliRootVersion("05-24-00"); //root version
1181 metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
1182 AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
1183 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
1184 gStorage->Put(fGainArray, id1, metaData);
1187 void AliTPCPreprocessorOffline::MakeQAPlot(Float_t FPtoMIPratio) {
1189 // Make QA plot to visualize results
1194 TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
1196 TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
1197 gainHistoCosmic->SetDirectory(0);
1198 gainHistoCosmic->SetName("GainHistoCosmic");
1199 gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
1200 gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1201 gainHistoCosmic->Draw("colz");
1202 fGraphCosmic->SetMarkerStyle(25);
1203 fGraphCosmic->Draw("lp");
1204 fGraphCosmic->SetMarkerStyle(25);
1205 TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
1207 for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
1208 grfFitCosmic->GetY()[i] *= FPtoMIPratio;
1210 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
1211 fGraphCosmic->GetY()[i] *= FPtoMIPratio;
1214 fGraphCosmic->Draw("lp");
1216 grfFitCosmic->SetLineColor(2);
1217 grfFitCosmic->Draw("lu");
1219 fGainArray->AddLast(gainHistoCosmic);
1220 //fGainArray->AddLast(canvasCosmic->Clone());
1221 delete canvasCosmic;
1224 TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
1226 TH2D * gainHistoMIP = fGainMIP->GetHistGainTime()->Projection(0,1);
1227 gainHistoMIP->SetName("GainHistoCosmic");
1228 gainHistoMIP->SetDirectory(0);
1229 gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
1230 gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1231 gainHistoMIP->Draw("colz");
1232 fGraphMIP->SetMarkerStyle(25);
1233 fGraphMIP->Draw("lp");
1234 TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
1235 grfFitMIP->SetLineColor(2);
1236 grfFitMIP->Draw("lu");
1237 fGainArray->AddLast(gainHistoMIP);
1238 //fGainArray->AddLast(canvasMIP->Clone());
1243 void AliTPCPreprocessorOffline::MakeFitTime(){
1245 // make aligment fit - store results in the file
1247 const Int_t kMinEntries=1000;
1249 MakePrimitivesTime();
1250 if (!fAlignTree) return;
1251 if (fAlignTree->GetEntries()<kMinEntries) return;
1252 fAlignTree->SetAlias("ptype","type");
1253 fAlignTree->SetAlias("hasITS","(1+0)");
1254 fAlignTree->SetAlias("dITS","1-2*(refX<40)");
1255 fAlignTree->SetAlias("isITS","refX>10");
1256 fAlignTree->SetAlias("isVertex","refX<10");
1258 Int_t npointsMax=30000000;
1259 TStatToolkit toolkit;
1265 TString fstringFast="";
1266 fstringFast+="FExBTwistX++";
1267 fstringFast+="FExBTwistY++";
1268 fstringFast+="FAlignRot0D++";
1269 fstringFast+="FAlignTrans0D++";
1270 fstringFast+="FAlignTrans1D++";
1272 fstringFast+="hasITS*FAlignTrans0++";
1273 fstringFast+="hasITS*FAlignTrans1++";
1274 fstringFast+="hasITS*FAlignRot0++";
1275 fstringFast+="hasITS*FAlignRot1++";
1276 fstringFast+="hasITS*FAlignRot2++";
1278 fstringFast+="dITS*FAlignTrans0++";
1279 fstringFast+="dITS*FAlignTrans1++";
1280 fstringFast+="dITS*FAlignRot0++";
1281 fstringFast+="dITS*FAlignRot1++";
1282 fstringFast+="dITS*FAlignRot2++";
1284 TCut cutFit="entries>10&&abs(mean)>0.00001&&rms>0";
1285 fAlignTree->SetAlias("err","rms");
1287 TString *strDeltaITS = TStatToolkit::FitPlaneConstrain(fAlignTree,"mean:err", fstringFast.Data(),cutFit, chi2,npoints,param,covar,-1,0, npointsMax, 1);
1288 strDeltaITS->Tokenize("++")->Print();
1289 fAlignTree->SetAlias("fitYFast",strDeltaITS->Data());
1291 TVectorD paramC= param;
1292 TMatrixD covarC= covar;
1293 TStatToolkit::Constrain1D(fstringFast,"Trans0D",paramC,covarC,0, 0.1);
1294 TStatToolkit::Constrain1D(fstringFast,"Trans1D",paramC,covarC,0, 0.1);
1295 TStatToolkit::Constrain1D(fstringFast,"TwistX",paramC,covarC,0, 0.1);
1296 TStatToolkit::Constrain1D(fstringFast,"TwistY",paramC,covarC,0, 0.1);
1297 TString strFitConst=TStatToolkit::MakeFitString(fstringFast, paramC,covar);
1298 fAlignTree->SetAlias("fitYFastC",strFitConst.Data());
1299 CreateAlignTime(fstringFast,paramC);
1305 void AliTPCPreprocessorOffline::MakeChainTime(){
1309 TFile f("CalibObjects.root");
1311 // const char *cdtype[7]={"ITS","TRD","Vertex","TOF","TPC","TPC0","TPC1"};
1312 //const char *cptype[5]={"dy","dz","dsnp","dtheta","d1pt"};
1313 const char * hname[5]={"dy","dz","dsnp","dtheta","d1pt"};
1315 AliTPCcalibTime *calibTime = 0;
1316 TObjArray * array = (TObjArray*)f.Get("TPCCalib");
1318 calibTime = (AliTPCcalibTime *)array->FindObject("calibTime");
1320 calibTime = (AliTPCcalibTime*)f.Get("calibTime");
1322 if (!calibTime) return;
1323 AliTPCCorrectionFit::CreateAlignMaps(AliTracker::GetBz(), run);
1324 TTreeSRedirector *pcstream = new TTreeSRedirector("meanITSVertex.root");
1327 THnSparse *his = calibTime->GetResHistoTPCITS(ihis);
1329 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1330 his->GetAxis(2)->SetRange(0,1000000);
1331 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1332 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
1335 his = calibTime->GetResHistoTPCITS(ihis);
1337 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1338 his->GetAxis(2)->SetRange(0,1000000);
1339 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1340 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
1343 his = calibTime->GetResHistoTPCITS(ihis);
1345 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1346 his->GetAxis(2)->SetRange(0,1000000);
1347 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1348 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
1351 his = calibTime->GetResHistoTPCvertex(ihis);
1353 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1354 his->GetAxis(2)->SetRange(0,1000000);
1355 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1356 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
1359 his = calibTime->GetResHistoTPCvertex(ihis);
1361 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1362 his->GetAxis(2)->SetRange(0,1000000);
1363 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1364 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
1368 his = calibTime->GetResHistoTPCvertex(ihis);
1370 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1371 his->GetAxis(2)->SetRange(0,1000000);
1372 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1373 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
1377 his = calibTime->GetResHistoTPCTOF(ihis);
1379 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1380 his->GetAxis(2)->SetRange(0,1000000);
1381 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1382 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("TOF%s",hname[ihis]),run,0.,ihis,3);
1386 his = calibTime->GetResHistoTPCTRD(ihis);
1388 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1389 his->GetAxis(2)->SetRange(0,1000000);
1390 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1391 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("TRD%s",hname[ihis]),run,0.,ihis,3);
1398 Double_t AliTPCPreprocessorOffline::EvalAt(Double_t phi, Double_t refX, Double_t theta, Int_t corr, Int_t ptype){
1402 Double_t sector = 9*phi/TMath::Pi();
1403 if (sector<0) sector+=18;
1404 Double_t y85=AliTPCCorrection::GetCorrSector(sector,85,theta,1,corr);
1405 Double_t y245=AliTPCCorrection::GetCorrSector(sector,245,theta,1,corr);
1406 if (ptype==0) return y85+(y245-y85)*(refX-85.)/(245.-85.);
1407 if (ptype==2) return (y245-y85)/(245.-85.);
1412 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){
1414 // Fit the distortion along the line with the parabolic model
1416 // phi0 - phi at the entrance of the TPC
1417 // snp - local inclination angle at the entrance of the TPC
1418 // refX - ref X where the distortion is evanluated
1421 static TLinearFitter fitter(3,"pol2");
1422 fitter.ClearPoints();
1423 if (nsteps<3) nsteps=3;
1424 Double_t deltaX=(245-85)/(nsteps);
1425 for (Int_t istep=0; istep<(nsteps+1); istep++){
1427 Double_t localX =85.+deltaX*istep;
1428 Double_t localPhi=phi0+deltaX*snp*istep;
1429 Double_t sector = 9*localPhi/TMath::Pi();
1430 if (sector<0) sector+=18;
1431 Double_t y=AliTPCCorrection::GetCorrSector(sector,localX,theta,1,corr);
1432 Double_t dlocalX=AliTPCCorrection::GetCorrSector(sector,localX,theta,0,corr);
1433 Double_t x[1]={localX-dlocalX};
1434 fitter.AddPoint(x,y);
1438 par[0]=fitter.GetParameter(0);
1439 par[1]=fitter.GetParameter(1);
1440 par[2]=fitter.GetParameter(2);
1442 if (ptype==0) return par[0]+par[1]*refX+par[2]*refX*refX;
1443 if (ptype==2) return par[1]+2*par[2]*refX;
1444 if (ptype==4) return par[2];
1453 void AliTPCPreprocessorOffline::MakePrimitivesTime(){
1455 // Create primitive transformation to fit
1457 fAlignTree=new TChain("fit","fit");
1458 fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdy");
1459 fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdsnp");
1460 fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdy");
1461 fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdsnp");
1463 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1464 Double_t bzField=AliTrackerBase::GetBz();
1465 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
1466 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
1467 Double_t wtP = -10.0 * (bzField) * vdrift / ezField ;
1468 AliTPCExBTwist *fitExBTwistX= new AliTPCExBTwist;
1469 AliTPCExBTwist *fitExBTwistY= new AliTPCExBTwist;
1470 AliTPCCalibGlobalMisalignment *trans0 =new AliTPCCalibGlobalMisalignment;
1471 AliTPCCalibGlobalMisalignment *trans1 =new AliTPCCalibGlobalMisalignment;
1472 AliTPCCalibGlobalMisalignment *trans0D =new AliTPCCalibGlobalMisalignment;
1473 AliTPCCalibGlobalMisalignment *trans1D =new AliTPCCalibGlobalMisalignment;
1474 AliTPCCalibGlobalMisalignment *rot0 =new AliTPCCalibGlobalMisalignment;
1475 AliTPCCalibGlobalMisalignment *rot1 =new AliTPCCalibGlobalMisalignment;
1476 AliTPCCalibGlobalMisalignment *rot2 =new AliTPCCalibGlobalMisalignment;
1477 AliTPCCalibGlobalMisalignment *rot3 =new AliTPCCalibGlobalMisalignment;
1480 fitExBTwistX->SetXTwist(0.001);
1481 fitExBTwistX->SetOmegaTauT1T2(wtP,1,1);
1483 fitExBTwistY->SetYTwist(0.001);
1484 fitExBTwistY->SetOmegaTauT1T2(wtP,1,1);
1486 TGeoHMatrix *matrixRot = new TGeoHMatrix;
1487 TGeoHMatrix *matrixX = new TGeoHMatrix;
1488 TGeoHMatrix *matrixY = new TGeoHMatrix;
1489 matrixX->SetDx(0.1);
1490 matrixY->SetDy(0.1);
1491 Double_t rotAngles0[9]={0};
1492 Double_t rotAngles1[9]={0};
1493 Double_t rotAngles2[9]={0};
1495 Double_t rotAngles3[9]={0};
1497 rotAngles0[0]=1; rotAngles0[4]=1; rotAngles0[8]=1;
1498 rotAngles1[0]=1; rotAngles1[4]=1; rotAngles1[8]=1;
1499 rotAngles2[0]=1; rotAngles2[4]=1; rotAngles2[8]=1;
1500 rotAngles3[0]=1; rotAngles3[4]=1; rotAngles3[8]=1;
1502 rotAngles0[1]=-0.001;rotAngles0[3]=0.001;
1503 rotAngles1[5]=-0.001;rotAngles1[7]=0.001;
1504 rotAngles2[2]=0.001;rotAngles2[6]=-0.001;
1505 rotAngles3[1]=0.001;rotAngles3[3]=-0.001;
1506 matrixRot->SetRotation(rotAngles0);
1507 rot0->SetAlignGlobal(matrixRot);
1508 matrixRot->SetRotation(rotAngles1);
1509 rot1->SetAlignGlobal(matrixRot);
1510 matrixRot->SetRotation(rotAngles2);
1511 rot2->SetAlignGlobal(matrixRot);
1512 matrixRot->SetRotation(rotAngles3);
1513 rot3->SetAlignGlobalDelta(matrixRot);
1515 trans0->SetAlignGlobal(matrixX);
1516 trans1->SetAlignGlobal(matrixY);
1517 trans0D->SetAlignGlobalDelta(matrixX);
1518 trans1D->SetAlignGlobalDelta(matrixY);
1519 fitExBTwistX->Init();
1520 fitExBTwistY->Init();
1522 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistX->Clone()),100);
1523 fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistY->Clone()),101);
1525 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot0->Clone()),102);
1526 fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(rot1->Clone()),103);
1527 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot2->Clone()),104);
1528 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot3->Clone()),105);
1530 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0->Clone()),106);
1531 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1->Clone()),107);
1532 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0D->Clone()),108);
1533 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1D->Clone()),109);
1535 fAlignTree->SetAlias("FExBTwistX", "AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,100,ptype)+0");
1536 fAlignTree->SetAlias("FExBTwistY","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,101,ptype)+0");
1537 fAlignTree->SetAlias("FAlignRot0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,102,ptype)+0");
1538 fAlignTree->SetAlias("FAlignRot0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,105,ptype)+0");
1539 fAlignTree->SetAlias("FAlignRot1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,103,ptype)+0");
1540 fAlignTree->SetAlias("FAlignRot2","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,104,ptype)+0");
1541 fAlignTree->SetAlias("FAlignTrans0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,106,ptype)+0");
1542 fAlignTree->SetAlias("FAlignTrans1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,107,ptype)+0");
1543 fAlignTree->SetAlias("FAlignTrans0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,108,ptype)+0");
1544 fAlignTree->SetAlias("FAlignTrans1D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,109,ptype)+0");
1546 // test fast function
1548 // fAlignTree->Draw("FExBTwistX:ExBTwistX","isITS&&ptype==0&&abs(snp)<0.05","");
1549 // fAlignTree->Draw("FExBTwistY:ExBTwistY","isITS&&ptype==0&&abs(snp)<0.05","");
1550 // fAlignTree->Draw("FAlignRot0:alignRot0","isITS&&ptype==0&&abs(snp)<0.05","");
1551 // fAlignTree->Draw("FAlignRot1:alignRot1","isITS&&ptype==0&&abs(snp)<0.05","");
1552 // fAlignTree->Draw("FAlignRot2:alignRot2","isITS&&ptype==0&&abs(snp)<0.05","");
1554 // fAlignTree->Draw("FAlignTrans0:alignTrans0","isITS&&ptype==0&&abs(snp)<0.05","");
1555 // fAlignTree->Draw("FAlignTrans1:alignTrans1","isITS&&ptype==0&&abs(snp)<0.05","");
1560 void AliTPCPreprocessorOffline::CreateAlignTime(TString fstring, TVectorD paramC){
1565 TGeoHMatrix *matrixDelta = new TGeoHMatrix;
1566 TGeoHMatrix *matrixGlobal = new TGeoHMatrix;
1567 Double_t rAngles[9];
1570 index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans0D");
1571 if (index>=0) matrixDelta->SetDx(paramC[index+1]*0.1);
1572 index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans1D");
1573 if (index>=0) matrixDelta->SetDy(paramC[index+1]*0.1);
1574 rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1575 index=TStatToolkit::GetFitIndex(fstring,"FAlignRot0D");
1576 rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1577 rAngles[5]=0; rAngles[7] =0;
1578 rAngles[2]=0; rAngles[6] =0;
1579 matrixDelta->SetRotation(rAngles);
1583 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans0");
1584 if (index>=0) matrixGlobal->SetDx(paramC[index+1]*0.1);
1585 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans1");
1586 if (index>=0) matrixGlobal->SetDy(paramC[index+1]*0.1);
1587 rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1588 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot0");
1589 rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1590 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot1");
1591 rAngles[5]=-paramC[index+1]*0.001; rAngles[7]=paramC[index+1]*0.001;
1592 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot2");
1593 rAngles[2]=paramC[index+1]*0.001; rAngles[6] =-paramC[index+1]*0.001;
1594 matrixGlobal->SetRotation(rAngles);
1596 AliTPCCalibGlobalMisalignment *fitAlignTime =0;
1597 fitAlignTime =new AliTPCCalibGlobalMisalignment;
1598 fitAlignTime->SetName("FitAlignTime");
1599 fitAlignTime->SetTitle("FitAlignTime");
1600 fitAlignTime->SetAlignGlobalDelta(matrixDelta);
1601 fitAlignTime->SetAlignGlobal(matrixGlobal);
1603 AliTPCExBTwist * fitExBTwist= new AliTPCExBTwist;
1604 Int_t indexX=TStatToolkit::GetFitIndex(fstring,"ExBTwistX");
1605 Int_t indexY=TStatToolkit::GetFitIndex(fstring,"ExBTwistY");
1606 fitExBTwist->SetXTwist(0.001*paramC[indexX+1]); // 1 mrad twist in x
1607 fitExBTwist->SetYTwist(0.001*paramC[indexY+1]); // 1 mrad twist in x
1608 fitExBTwist->SetName("FitExBTwistTime");
1609 fitExBTwist->SetTitle("FitExBTwistTime");
1610 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1611 Double_t bzField=AliTrackerBase::GetBz();
1612 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
1614 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
1615 Double_t wt = -10.0 * (bzField) * vdrift / ezField ;
1617 fitExBTwist->SetOmegaTauT1T2(wt,1,1);
1618 fitExBTwist->Init();
1620 AliTPCComposedCorrection *corrTime = new AliTPCComposedCorrection;
1621 TObjArray *arr = new TObjArray;
1622 corrTime->SetCorrections(arr);
1624 corrTime->GetCorrections()->Add(fitExBTwist);
1625 corrTime->GetCorrections()->Add(fitAlignTime);
1626 corrTime->SetName("FitCorrectionTime");
1627 corrTime->SetTitle("FitCorrectionTime");
1629 fitExBTwist->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwist->Clone()),1001);
1630 fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(fitAlignTime->Clone()),1002);
1631 fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(corrTime->Clone()),1003);
1634 fAlignTree->SetAlias("ExBTwistTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1001,ptype)+0");
1635 fAlignTree->SetAlias("AlignTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1002,ptype)+0");
1636 fAlignTree->SetAlias("FitCorrectionTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1003,ptype)+0");
1639 TFile *f = new TFile("fitITSVertex.root","update");
1640 corrTime->Write("FitCorrectionTime");