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 "AliSplineFit.h"
79 #include "AliTPCComposedCorrection.h"
80 #include "AliTPCExBTwist.h"
81 #include "AliTPCCalibGlobalMisalignment.h"
82 #include "TStatToolkit.h"
85 #include "AliTrackerBase.h"
86 #include "AliTPCPreprocessorOffline.h"
89 ClassImp(AliTPCPreprocessorOffline)
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),
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
110 fAlignTree(0), // alignment tree
111 fSwitchOnValidation(kFALSE), // flag to switch on validation of OCDB parameters
117 // default constructor
121 AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
130 void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const timeDrift){
132 // find the fist and last run
134 TObjArray *hisArray =timeDrift->GetHistoDrift();
135 {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
136 THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
137 if (!addHist) continue;
138 if (addHist->GetEntries()<fMinEntries) continue;
139 TH1D* histo =addHist->Projection(3);
140 TH1D* histoTime=addHist->Projection(0);
141 printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
144 startRun=histo->FindFirstBinAbove(0);
145 endRun =histo->FindLastBinAbove(0);
147 startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun);
148 endRun =TMath::Max(histo->FindLastBinAbove(0),endRun);
151 startTime=histoTime->FindFirstBinAbove(0);
152 endTime =histoTime->FindLastBinAbove(0);
154 startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime);
155 endTime =TMath::Max(histoTime->FindLastBinAbove(0),endTime);
160 if (startRun<0) startRun=0;
161 if (endRun<0) endRun=100000000;
162 printf("Run range :\t%d-%d\n", startRun, endRun);
163 printf("Time range :\t%d-%d\n", startTime, endTime);
169 void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, TString pocdbStorage){
171 // make calibration of the drift velocity
173 // file - the location of input file
174 // ustartRun, uendrun - run validity period
175 // pocdbStorage - path to hte OCDB storage
176 // - if empty - local storage 'pwd' uesed
177 if (pocdbStorage.Length()>0) ocdbStorage=pocdbStorage;
179 ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
181 // 1. Initialization and run range setting
183 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
185 fTimeDrift = (AliTPCcalibTime *)array->FindObject("calibTime");
187 fTimeDrift = (AliTPCcalibTime*)fcalib.Get("calibTime");
189 if(!fTimeDrift) return;
193 TObjArray *hisArray =fTimeDrift->GetHistoDrift();
194 GetRunRange(fTimeDrift);
195 for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
196 THnSparse* addHist=(THnSparse*)hisArray->At(i);
197 if (!addHist) continue;
198 if (startTime<endTime) addHist->GetAxis(0)->SetRange(startTime-1,endTime+1);
199 if (startRun<endRun) addHist->GetAxis(3)->SetRange(startRun-1,endRun+1);
203 // 2. extraction of the information
205 fVdriftArray = new TObjArray();
206 AddAlignmentGraphs(fVdriftArray,fTimeDrift);
207 AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
208 AddLaserGraphs(fVdriftArray,fTimeDrift);
210 // 3. Append QA plots
212 MakeDefaultPlots(fVdriftArray,fVdriftArray);
215 // 4. validate OCDB entries
217 if(fSwitchOnValidation==kTRUE && ValidateTimeDrift()==kFALSE) {
218 Printf("TPC time drift OCDB parameters out of range!");
225 TFile * ftime= TFile::Open("fitITSVertex.root");
227 TObject * alignmentTime=ftime->Get("FitCorrectionTime");
228 if (alignmentTime) fVdriftArray->AddLast(alignmentTime);
235 UpdateOCDBDrift(ustartRun,uendRun,ocdbStorage);
238 void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun, const char* storagePath ){
242 AliCDBMetaData *metaData= new AliCDBMetaData();
243 metaData->SetObjectClassName("TObjArray");
244 metaData->SetResponsible("Marian Ivanov");
245 metaData->SetBeamPeriod(1);
246 metaData->SetAliRootVersion("05-25-01"); //root version
247 metaData->SetComment("Calibration of the time dependence of the drift velocity");
249 id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
250 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
251 gStorage->Put(fVdriftArray, (*id1), metaData);
254 Bool_t AliTPCPreprocessorOffline::ValidateTimeGain()
257 // Validate time gain corrections
259 Printf("ValidateTimeGain..." );
260 Float_t minGain = fMinGain;
261 Float_t maxGain = fMaxGain;
263 TGraphErrors *gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
265 gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
266 if (!gr) return kFALSE;
267 Printf("Assuming given run is a cosmic run. Using gain calibration from Fermi-plateau muons.");
269 if(gr->GetN()<1) return kFALSE;
271 // check whether gain in the range
272 for(Int_t iPoint=0; iPoint<gr->GetN(); iPoint++)
274 if(gr->GetY()[iPoint] < minGain || gr->GetY()[iPoint] > maxGain)
282 Bool_t AliTPCPreprocessorOffline::ValidateTimeDrift()
285 // Validate time drift velocity corrections
287 Printf("ValidateTimeDrift..." );
289 Float_t maxVDriftCorr = fMaxVdriftCorr;
291 TGraphErrors* gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
292 Printf("ALIGN_ITSB_TPC_DRIFTVD graph = %p",gr);
294 if(!gr) return kFALSE;
296 Printf("ALIGN_ITSB_TPC_DRIFTVD number of points = %d",gr->GetN());
300 // check whether drift velocity corrections in the range
301 for(Int_t iPoint = 0; iPoint<gr->GetN(); iPoint++)
303 Printf("Y value from the graph: %f",TMath::Abs(gr->GetY()[iPoint]));
304 if(TMath::Abs(gr->GetY()[iPoint]) > maxVDriftCorr)
311 void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
313 // update the OCDB entry for the nominal time0
316 // AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
317 AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
318 TGraphErrors *grT = (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
319 Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
320 Double_t deltaT = deltaTcm/param->GetDriftV();
321 paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
324 AliCDBMetaData *metaData= new AliCDBMetaData();
325 metaData->SetObjectClassName("TObjArray");
326 metaData->SetResponsible("Marian Ivanov");
327 metaData->SetBeamPeriod(1);
328 metaData->SetAliRootVersion("05-25-02"); //root version
329 metaData->SetComment("Updated calibration of nominal time 0");
331 id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
332 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
333 gStorage->Put(param, (*id1), metaData);
338 void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
340 // Print the names of the entries in array
342 Int_t entries = array->GetEntries();
343 for (Int_t i=0; i<entries; i++){
344 if (!array->At(i)) continue;
345 printf("%d\t %s\n", i, array->At(i)->GetName());
351 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
353 // 1. filter graph - error cut errSigmaCut
354 // 2. filter graph - medianCutAbs around median
356 // errSigmaCut - cut on error
357 // medianCutAbs - cut on value around median
360 // 1. filter graph - error cut errSigmaCut
362 TGraphErrors *graphF;
363 graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
365 if (!graphF) return 0;
366 graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
368 if (!graph) return 0;
370 // filter graph - kMedianCutAbs around median
372 graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
374 if (!graphF) return 0;
375 graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
377 if (!graph) return 0;
383 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
385 // filter outlyer measurement
386 // Only points around median +- cut filtered
388 if (!graph) return 0;
390 Int_t npoints0 = graph->GetN();
393 Double_t *outx=new Double_t[npoints0];
394 Double_t *outy=new Double_t[npoints0];
395 Double_t *errx=new Double_t[npoints0];
396 Double_t *erry=new Double_t[npoints0];
399 if (npoints0<kMinPoints) {
406 for (Int_t iter=0; iter<3; iter++){
408 for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
409 if (graph->GetY()[ipoint]==0) continue;
410 if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;
411 outx[npoints] = graph->GetX()[ipoint];
412 outy[npoints] = graph->GetY()[ipoint];
413 errx[npoints] = graph->GetErrorX(ipoint);
414 erry[npoints] = graph->GetErrorY(ipoint);
417 if (npoints<=1) break;
418 medianY =TMath::Median(npoints,outy);
419 rmsY =TMath::RMS(npoints,outy);
421 TGraphErrors *graphOut=0;
422 if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry);
431 void AliTPCPreprocessorOffline::AddHistoGraphs( TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
433 // Add graphs corresponding to the alignment
435 const Double_t kErrSigmaCut=5; // error sigma cut - for filtering
436 const Double_t kMedianCutAbs=0.03; // error sigma cut - for filtering
438 TObjArray * array=timeDrift->GetHistoDrift();
440 THnSparse* hist=NULL;
441 // 2.a) cosmics with different triggers
442 for (Int_t i=0; i<array->GetEntriesFast();i++){
443 hist=(THnSparseF*)array->UncheckedAt(i);
445 if (hist->GetEntries()<minEntries) continue;
447 TString name=hist->GetName();
448 Int_t dim[4]={0,1,2,3};
449 THnSparse* newHist=hist->Projection(4,dim);
450 newHist->SetName(name);
451 TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
453 printf("Graph =%s filtered out\n", name.Data());
456 printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
457 Int_t pos=name.Index("_");
458 name=name(pos,name.Capacity()-pos);
459 TString graphName=graph->ClassName();
463 graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
466 graph->SetMarkerStyle(i%8+20);
467 graph->SetMarkerColor(i%7);
468 graph->GetXaxis()->SetTitle("Time");
469 graph->GetYaxis()->SetTitle("v_{dcor}");
470 graph->SetName(graphName);
471 graph->SetTitle(graphName);
472 printf("Graph %d\t=\t%s\n", i, graphName.Data());
473 vdriftArray->Add(graph);
482 void AliTPCPreprocessorOffline::AddAlignmentGraphs( TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
484 // Add graphs corresponding to alignment to the object array
486 TObjArray *arrayITS=0;
487 TObjArray *arrayTOF=0;
488 TObjArray *arrayTRD=0;
489 TMatrixD *mstatITS=0;
490 TMatrixD *mstatTOF=0;
491 TMatrixD *mstatTRD=0;
493 arrayITS=timeDrift->GetAlignITSTPC();
494 arrayTRD=timeDrift->GetAlignTRDTPC();
495 arrayTOF=timeDrift->GetAlignTOFTPC();
497 if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.9,50,fMaxVdriftCorr);
498 if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.9,1000,fMaxVdriftCorr);
499 if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.9,50,fMaxVdriftCorr);
501 TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
502 TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
503 TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
504 TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
505 TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
506 TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
508 TObjArray * arrayTRDP= 0x0;
509 TObjArray * arrayTRDM= 0x0;
510 TObjArray * arrayTRDB= 0x0;
511 arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
512 arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
513 arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
516 Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
517 TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
518 arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
519 arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
520 TString grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
521 "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
522 "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
523 TString grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
524 "DELTAX", "DELTAY", "DELTAZ",
525 "DRIFTVD", "T0", "VDGY"};
528 TVectorD vX(entries);
529 TVectorD vY(entries);
530 TVectorD vEx(entries);
531 TVectorD vEy(entries);
533 for (Int_t iarray=0; iarray<12; iarray++){
534 arr = arrays[iarray];
535 if (arr==0) continue;
536 for (Int_t ipar=0; ipar<9; ipar++){
538 for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
539 AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
540 if (!kalman) continue;
541 vX[counter]=kalman->GetTimeStamp();
542 vY[counter]=(*(kalman->GetState()))[ipar];
543 if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
545 vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
549 TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
551 vEx.GetMatrixArray(),
552 vEy.GetMatrixArray());
553 TString grName=grnames[iarray];
556 graph->SetName(grName.Data());
557 vdriftArray->AddLast(graph);
565 void AliTPCPreprocessorOffline::AddLaserGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
567 // add graphs for laser
569 const Double_t delayL0L1 = 0.071; //this is hack for 1/2 weeks
571 TGraphErrors *grLaser[6]={0,0,0,0,0,0};
572 //hisN = timeDrift->GetHistVdriftLaserA(0);
573 if (timeDrift->GetHistVdriftLaserA(0)){
574 grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
575 grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
576 vdriftArray->AddLast(grLaser[0]);
578 if (timeDrift->GetHistVdriftLaserA(1)){
579 grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
580 grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
581 vdriftArray->AddLast(grLaser[1]);
583 if (timeDrift->GetHistVdriftLaserA(2)){
584 grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
585 grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
586 vdriftArray->AddLast(grLaser[2]);
588 if (timeDrift->GetHistVdriftLaserC(0)){
589 grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
590 grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
591 vdriftArray->AddLast(grLaser[3]);
593 if (timeDrift->GetHistVdriftLaserC(1)){
594 grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
595 grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
596 vdriftArray->AddLast(grLaser[4]);
598 if (timeDrift->GetHistVdriftLaserC(2)){
599 grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
600 grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
601 vdriftArray->AddLast(grLaser[5]);
603 for (Int_t i=0; i<6;i++){
605 SetDefaultGraphDrift(grLaser[i], 1,(i+20));
606 grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
612 TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
614 // Make graph with mean values and rms
616 hisN->GetAxis(itime)->SetRange(0,100000000);
617 hisN->GetAxis(ival)->SetRange(0,100000000);
618 TH1 * hisT = hisN->Projection(itime);
619 TH1 * hisV = hisN->Projection(ival);
621 Int_t firstBinA = hisT->FindFirstBinAbove(2);
622 Int_t lastBinA = hisT->FindLastBinAbove(2);
623 Int_t firstBinV = hisV->FindFirstBinAbove(0);
624 Int_t lastBinV = hisV->FindLastBinAbove(0);
625 hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
626 hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
628 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
629 Double_t cont = hisT->GetBinContent(ibin);
630 if (cont<minEntries) continue;
633 TVectorD vecTime(entries);
634 TVectorD vecMean0(entries);
635 TVectorD vecRMS0(entries);
636 TVectorD vecMean1(entries);
637 TVectorD vecRMS1(entries);
639 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
640 Double_t cont = hisT->GetBinContent(ibin);
641 if (cont<minEntries) continue;
642 //hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
643 Int_t minBin = ibin-1;
644 Int_t maxBin = ibin+1;
645 if(minBin <= 0) minBin = 1;
646 if(maxBin >= hisN->GetAxis(itime)->GetNbins()) maxBin = hisN->GetAxis(itime)->GetNbins()-1;
647 hisN->GetAxis(itime)->SetRange(minBin,maxBin);
649 Double_t time = hisT->GetBinCenter(ibin);
650 TH1 * his = hisN->Projection(ival);
651 Double_t nentries0= his->GetBinContent(his->FindBin(0));
652 if (cont-nentries0<minEntries) continue;
654 his->SetBinContent(his->FindBin(0),0);
655 vecTime[entries]=time;
656 vecMean0[entries]=his->GetMean()+offset;
657 vecMean1[entries]=his->GetMeanError();
658 vecRMS0[entries] =his->GetRMS();
659 vecRMS1[entries] =his->GetRMSError();
665 TGraphErrors * graph = new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(), 0, vecMean1.GetMatrixArray());
676 void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
678 // Set default style for QA views
680 graph->GetXaxis()->SetTimeDisplay(kTRUE);
681 graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
682 graph->SetMaximum( 0.025);
683 graph->SetMinimum(-0.025);
684 graph->GetXaxis()->SetTitle("Time");
685 graph->GetYaxis()->SetTitle("v_{dcorr}");
687 graph->GetYaxis()->SetLabelSize(0.03);
688 graph->GetXaxis()->SetLabelSize(0.03);
690 graph->GetXaxis()->SetNdivisions(10,5,0);
691 graph->GetYaxis()->SetNdivisions(10,5,0);
693 graph->GetXaxis()->SetLabelOffset(0.02);
694 graph->GetYaxis()->SetLabelOffset(0.005);
696 graph->GetXaxis()->SetTitleOffset(1.3);
697 graph->GetYaxis()->SetTitleOffset(1.2);
699 graph->SetMarkerColor(color);
700 graph->SetLineColor(color);
701 graph->SetMarkerStyle(style);
704 void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
706 // Set default pad style for QA
709 pad->SetMargin(mx0,mx1,my0,my1);
713 void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray * /*picArray*/){
715 // 0. make a default QA plots
716 // 1. Store them in the array
719 Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
721 TGraphErrors* laserA =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
722 TGraphErrors* laserC =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
723 TGraphErrors* cosmic =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
724 TGraphErrors* cross =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
725 TGraphErrors* itstpcP =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
726 TGraphErrors* itstpcM =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
727 TGraphErrors* itstpcB =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
729 if (laserA) SetDefaultGraphDrift(laserA,2,25);
730 if (laserC) SetDefaultGraphDrift(laserC,4,26);
731 if (cosmic) SetDefaultGraphDrift(cosmic,3,27);
732 if (cross) SetDefaultGraphDrift(cross,4,28);
733 if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
734 if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
735 if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
743 pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
745 SetPadStyle(pad,mx0,mx1,my0,my1);
748 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
749 legend->AddEntry(laserA,"Laser A side");
750 legend->AddEntry(laserC,"Laser C side");
752 //picArray->AddLast(pad);
755 if (itstpcP&&itstpcM&&itstpcB){
756 pad = new TCanvas("ITSTPC","ITSTPC");
757 itstpcP->Draw("alp");
758 SetPadStyle(pad,mx0,mx1,my0,my1);
759 itstpcP->Draw("alp");
761 itstpcM->Draw("apl");
764 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
765 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
766 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
767 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
769 //picArray->AddLast(pad);
772 if (itstpcB&&laserA&&itstpcP&&itstpcM){
773 pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
774 SetPadStyle(pad,mx0,mx1,my0,my1);
779 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
780 legend->AddEntry(laserA,"TPC laser");
781 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
782 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
783 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
785 //picArray->AddLast(pad);
789 pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
790 SetPadStyle(pad,mx0,mx1,my0,my1);
791 itstpcP->Draw("alp");
796 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
798 legend->AddEntry(cross,"TPC cross tracks");
799 legend->AddEntry(itstpcB,"ITS-TPC smooth");
801 //picArray->AddLast(pad);
803 if (itstpcP&&cosmic){
804 pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
805 SetPadStyle(pad,mx0,mx1,my0,my1);
806 itstpcP->Draw("alp");
811 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
813 legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
814 legend->AddEntry(itstpcB,"ITS-TPC smooth");
816 //picArray->AddLast(pad);
823 void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber, TString pocdbStorage){
827 if (pocdbStorage.Length()==0) pocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
830 // 1. Read gain values
832 ReadGainGlobal(fileName);
835 // 2. Extract calibration values
837 AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
838 AnalyzeAttachment(startRunNumber,endRunNumber);
839 AnalyzePadRegionGain();
840 AnalyzeGainMultiplicity();
842 // 3. Make control plots
847 // 4. validate OCDB entries
849 if(fSwitchOnValidation==kTRUE && ValidateTimeGain()==kFALSE) {
850 Printf("TPC time gain OCDB parameters out of range!");
857 UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage.Data());
860 void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
862 // read calibration entries from file
864 TFile fcalib(fileName);
865 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
867 fGainMIP = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
868 fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
869 fGainMult = ( AliTPCcalibGainMult *)array->FindObject("calibGainMult");
871 fGainMIP = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
872 fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
873 fGainMult = ( AliTPCcalibGainMult *)fcalib.Get("calibGainMult");
876 TFile fcalibMult("TPCMultObjects.root");
877 fGainMult = ( AliTPCcalibGainMult *)fcalibMult.Get("calibGainMult");
880 Int_t firstBinA =0, lastBinA=0;
883 hisT= fGainCosmic->GetHistGainTime()->Projection(1);
884 firstBinA = hisT->FindFirstBinAbove(2);
885 lastBinA = hisT->FindLastBinAbove(2);
886 fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
891 hisT= fGainMIP->GetHistGainTime()->Projection(1);
892 firstBinA = hisT->FindFirstBinAbove(2);
893 lastBinA = hisT->FindLastBinAbove(2);
894 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
902 Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit, Float_t FPtoMIPratio){
904 // Analyze gain - produce the calibration graphs
907 // 1.) try to create MIP spline
910 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
911 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
912 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
914 fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
915 if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
916 if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
917 if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
918 fGainArray->AddAt(fFitMIP,0);
921 // 2.) try to create Cosmic spline
924 fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
925 fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
927 fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
928 if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
931 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
932 fGraphCosmic->GetY()[i] /= FPtoMIPratio;
933 fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
937 if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
938 if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
939 fGainArray->AddAt(fFitCosmic,1);
941 // with naming convention and backward compatibility
942 fGainArray->AddAt(fGraphMIP,2);
943 fGainArray->AddAt(fGraphCosmic,3);
944 cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
949 Bool_t AliTPCPreprocessorOffline::AnalyzeAttachment(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesFit) {
951 // determine slope as a function of mean driftlength
953 if(!fGainMIP) return kFALSE;
955 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
957 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
958 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
960 fGainMIP->GetHistGainTime()->GetAxis(3)->SetRangeUser(125,250);// only full tracking region (driftlength)
961 fGainMIP->GetHistGainTime()->GetAxis(0)->SetRangeUser(1.5,3.5);// only full tracking region (driftlength)
963 TH3D * hist = fGainMIP->GetHistGainTime()->Projection(1, 0, 3);
965 Double_t *xvec = new Double_t[hist->GetNbinsX()];
966 Double_t *yvec = new Double_t[hist->GetNbinsX()];
967 Double_t *xerr = new Double_t[hist->GetNbinsX()];
968 Double_t *yerr = new Double_t[hist->GetNbinsX()];
971 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
975 for (Int_t idelta=0; idelta<5; idelta++){
977 imin = TMath::Max(i-idelta,1);
978 imax = TMath::Min(i+idelta,hist->GetNbinsX());
979 nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
980 //if (nsum==0) break;
981 if (nsum>minEntriesFit) break;
983 if (nsum<minEntriesFit) continue;
985 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRangeUser(hist->GetXaxis()->GetBinCenter(imin-1),hist->GetXaxis()->GetBinCenter(imax+1)); // define time range
986 TH2D * histZdep = fGainMIP->GetHistGainTime()->Projection(0,3);
988 histZdep->FitSlicesY(0,0,-1,0,"QNR",&arr);
989 TH1D * driftDep = (TH1D*)arr.At(1);
991 //TGraphErrors * driftDep = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,3,100,1,0.,1);
992 /*if (driftDep->GetN() < 4) {
999 TF1 pol1("polynom1","pol1",125,240);
1000 //driftDep->Fit(&pol1,"QNRROB=0.8");
1001 driftDep->Fit(&pol1,"QNR");
1002 xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin-1)+hist->GetXaxis()->GetBinCenter(imax+1));
1003 yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
1004 xerr[counter] = hist->GetXaxis()->GetBinCenter(imax+1)-hist->GetXaxis()->GetBinCenter(imin-1);
1005 yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
1011 fGraphAttachmentMIP = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
1012 if (fGraphAttachmentMIP) fGraphAttachmentMIP->SetName("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");// set proper names according to naming convention
1013 fGainArray->AddLast(fGraphAttachmentMIP);
1021 if (counter < 1) return kFALSE;
1027 Bool_t AliTPCPreprocessorOffline::AnalyzePadRegionGain(){
1029 // Analyze gain for different pad regions - produce the calibration graphs 0,1,2
1033 TH2D * histQmax = (TH2D*) fGainMult->GetHistPadEqual()->Projection(0,2);
1034 TH2D * histQtot = (TH2D*) fGainMult->GetHistPadEqual()->Projection(1,2);
1037 histQmax->FitSlicesY(0,0,-1,0,"QNR",&arr);
1038 Double_t xMax[3] = {0,1,2};
1039 Double_t yMax[3] = {((TH1D*)arr.At(1))->GetBinContent(1),
1040 ((TH1D*)arr.At(1))->GetBinContent(2),
1041 ((TH1D*)arr.At(1))->GetBinContent(3)};
1042 Double_t yMaxErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1043 ((TH1D*)arr.At(1))->GetBinError(2),
1044 ((TH1D*)arr.At(1))->GetBinError(3)};
1045 TGraphErrors * fitPadRegionQmax = new TGraphErrors(3, xMax, yMax, 0, yMaxErr);
1047 histQtot->FitSlicesY(0,0,-1,0,"QNR",&arr);
1048 Double_t xTot[3] = {0,1,2};
1049 Double_t yTot[3] = {((TH1D*)arr.At(1))->GetBinContent(1),
1050 ((TH1D*)arr.At(1))->GetBinContent(2),
1051 ((TH1D*)arr.At(1))->GetBinContent(3)};
1052 Double_t yTotErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1053 ((TH1D*)arr.At(1))->GetBinError(2),
1054 ((TH1D*)arr.At(1))->GetBinError(3)};
1055 TGraphErrors * fitPadRegionQtot = new TGraphErrors(3, xTot, yTot, 0, yTotErr);
1057 fitPadRegionQtot->SetName("TGRAPHERRORS_MEANQTOT_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1058 fitPadRegionQmax->SetName("TGRAPHERRORS_MEANQMAX_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1060 fGainArray->AddLast(fitPadRegionQtot);
1061 fGainArray->AddLast(fitPadRegionQmax);
1069 Bool_t AliTPCPreprocessorOffline::AnalyzeGainMultiplicity() {
1071 // Analyze gain as a function of multiplicity and produce calibration graphs
1073 if (!fGainMult) return kFALSE;
1074 fGainMult->GetHistGainMult()->GetAxis(3)->SetRangeUser(3,3);
1075 TH2D * histMultMax = fGainMult->GetHistGainMult()->Projection(0,4);
1076 TH2D * histMultTot = fGainMult->GetHistGainMult()->Projection(1,4);
1077 histMultMax->RebinX(4);
1078 histMultTot->RebinX(4);
1082 histMultMax->FitSlicesY(0,0,-1,0,"QNR",&arrMax);
1083 histMultTot->FitSlicesY(0,0,-1,0,"QNR",&arrTot);
1085 TH1D * meanMax = (TH1D*)arrMax.At(1);
1086 TH1D * meanTot = (TH1D*)arrTot.At(1);
1087 Float_t meanMult = histMultMax->GetMean();
1088 if(meanMax->GetBinContent(meanMax->FindBin(meanMult))) {
1089 meanMax->Scale(1./meanMax->GetBinContent(meanMax->FindBin(meanMult)));
1094 if(meanTot->GetBinContent(meanTot->FindBin(meanMult))) {
1095 meanTot->Scale(1./meanTot->GetBinContent(meanTot->FindBin(meanMult)));
1100 Float_t xMultMax[50];
1101 Float_t yMultMax[50];
1102 Float_t yMultErrMax[50];
1103 Float_t xMultTot[50];
1104 Float_t yMultTot[50];
1105 Float_t yMultErrTot[50];
1107 Int_t nCountMax = 0;
1108 for(Int_t iBin = 1; iBin < meanMax->GetXaxis()->GetNbins(); iBin++) {
1109 Float_t yValMax = meanMax->GetBinContent(iBin);
1110 if (yValMax < 0.7) continue;
1111 if (yValMax > 1.3) continue;
1112 if (meanMax->GetBinError(iBin)/yValMax > 0.01) continue;
1113 xMultMax[nCountMax] = meanMax->GetXaxis()->GetBinCenter(iBin);
1114 yMultMax[nCountMax] = yValMax;
1115 yMultErrMax[nCountMax] = meanMax->GetBinError(iBin);
1119 if (nCountMax < 10) return kFALSE;
1120 TGraphErrors * fitMultMax = new TGraphErrors(nCountMax, xMultMax, yMultMax, 0, yMultErrMax);
1121 fitMultMax->SetName("TGRAPHERRORS_MEANQMAX_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1123 Int_t nCountTot = 0;
1124 for(Int_t iBin = 1; iBin < meanTot->GetXaxis()->GetNbins(); iBin++) {
1125 Float_t yValTot = meanTot->GetBinContent(iBin);
1126 if (yValTot < 0.7) continue;
1127 if (yValTot > 1.3) continue;
1128 if (meanTot->GetBinError(iBin)/yValTot > 0.1) continue;
1129 xMultTot[nCountTot] = meanTot->GetXaxis()->GetBinCenter(iBin);
1130 yMultTot[nCountTot] = yValTot;
1131 yMultErrTot[nCountTot] = meanTot->GetBinError(iBin);
1135 if (nCountTot < 10) return kFALSE;
1136 TGraphErrors * fitMultTot = new TGraphErrors(nCountTot, xMultTot, yMultTot, 0, yMultErrTot);
1137 fitMultTot->SetName("TGRAPHERRORS_MEANQTOT_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1139 fGainArray->AddLast(fitMultMax);
1140 fGainArray->AddLast(fitMultTot);
1147 void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
1149 // Update OCDB entry
1151 AliCDBMetaData *metaData= new AliCDBMetaData();
1152 metaData->SetObjectClassName("TObjArray");
1153 metaData->SetResponsible("Alexander Kalweit");
1154 metaData->SetBeamPeriod(1);
1155 metaData->SetAliRootVersion("05-24-00"); //root version
1156 metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
1157 AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
1158 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
1159 gStorage->Put(fGainArray, id1, metaData);
1162 void AliTPCPreprocessorOffline::MakeQAPlot(Float_t FPtoMIPratio) {
1164 // Make QA plot to visualize results
1169 TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
1171 TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
1172 gainHistoCosmic->SetDirectory(0);
1173 gainHistoCosmic->SetName("GainHistoCosmic");
1174 gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
1175 gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1176 gainHistoCosmic->Draw("colz");
1177 fGraphCosmic->SetMarkerStyle(25);
1178 fGraphCosmic->Draw("lp");
1179 fGraphCosmic->SetMarkerStyle(25);
1180 TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
1182 for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
1183 grfFitCosmic->GetY()[i] *= FPtoMIPratio;
1185 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
1186 fGraphCosmic->GetY()[i] *= FPtoMIPratio;
1189 fGraphCosmic->Draw("lp");
1191 grfFitCosmic->SetLineColor(2);
1192 grfFitCosmic->Draw("lu");
1194 fGainArray->AddLast(gainHistoCosmic);
1195 //fGainArray->AddLast(canvasCosmic->Clone());
1196 delete canvasCosmic;
1199 TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
1201 TH2D * gainHistoMIP = fGainMIP->GetHistGainTime()->Projection(0,1);
1202 gainHistoMIP->SetName("GainHistoCosmic");
1203 gainHistoMIP->SetDirectory(0);
1204 gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
1205 gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1206 gainHistoMIP->Draw("colz");
1207 fGraphMIP->SetMarkerStyle(25);
1208 fGraphMIP->Draw("lp");
1209 TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
1210 grfFitMIP->SetLineColor(2);
1211 grfFitMIP->Draw("lu");
1212 fGainArray->AddLast(gainHistoMIP);
1213 //fGainArray->AddLast(canvasMIP->Clone());
1218 void AliTPCPreprocessorOffline::MakeFitTime(){
1220 // mak aligment fit - store results in the file
1222 const Int_t kMinEntries=1000;
1224 MakePrimitivesTime();
1225 if (!fAlignTree) return;
1226 if (fAlignTree->GetEntries()<kMinEntries) return;
1227 fAlignTree->SetAlias("ptype","type");
1228 fAlignTree->SetAlias("hasITS","(1+0)");
1229 fAlignTree->SetAlias("dITS","1-2*(refX<40)");
1230 fAlignTree->SetAlias("isITS","refX>10");
1231 fAlignTree->SetAlias("isVertex","refX<10");
1233 Int_t npointsMax=30000000;
1234 TStatToolkit toolkit;
1240 TString fstringFast="";
1241 fstringFast+="FExBTwistX++";
1242 fstringFast+="FExBTwistY++";
1243 fstringFast+="FAlignRot0D++";
1244 fstringFast+="FAlignTrans0D++";
1245 fstringFast+="FAlignTrans1D++";
1247 fstringFast+="hasITS*FAlignTrans0++";
1248 fstringFast+="hasITS*FAlignTrans1++";
1249 fstringFast+="hasITS*FAlignRot0++";
1250 fstringFast+="hasITS*FAlignRot1++";
1251 fstringFast+="hasITS*FAlignRot2++";
1253 fstringFast+="dITS*FAlignTrans0++";
1254 fstringFast+="dITS*FAlignTrans1++";
1255 fstringFast+="dITS*FAlignRot0++";
1256 fstringFast+="dITS*FAlignRot1++";
1257 fstringFast+="dITS*FAlignRot2++";
1259 TCut cutFit="entries>10&&abs(mean)>0.00001";
1260 fAlignTree->SetAlias("err","rms");
1262 TString *strDeltaITS = TStatToolkit::FitPlaneConstrain(fAlignTree,"mean:err", fstringFast.Data(),cutFit, chi2,npoints,param,covar,-1,0, npointsMax, 1);
1263 strDeltaITS->Tokenize("++")->Print();
1264 fAlignTree->SetAlias("fitYFast",strDeltaITS->Data());
1266 TVectorD paramC= param;
1267 TMatrixD covarC= covar;
1268 TStatToolkit::Constrain1D(fstringFast,"Trans0D",paramC,covarC,0, 0.1);
1269 TStatToolkit::Constrain1D(fstringFast,"Trans1D",paramC,covarC,0, 0.1);
1270 TStatToolkit::Constrain1D(fstringFast,"TwistX",paramC,covarC,0, 0.1);
1271 TStatToolkit::Constrain1D(fstringFast,"TwistY",paramC,covarC,0, 0.1);
1272 TString strFitConst=TStatToolkit::MakeFitString(fstringFast, paramC,covar);
1273 fAlignTree->SetAlias("fitYFastC",strFitConst.Data());
1274 CreateAlignTime(fstringFast,paramC);
1280 void AliTPCPreprocessorOffline::MakeChainTime(){
1282 TFile f("CalibObjects.root");
1283 // const char *cdtype[7]={"ITS","TRD","Vertex","TOF","TPC","TPC0","TPC1"};
1284 //const char *cptype[5]={"dy","dz","dsnp","dtheta","d1pt"};
1285 const char * hname[5]={"dy","dz","dsnp","dtheta","d1pt"};
1287 AliTPCcalibTime *calibTime= (AliTPCcalibTime*) f.Get("calibTime");
1288 if (!calibTime) return;
1289 TTreeSRedirector *pcstream = new TTreeSRedirector("meanITSVertex.root");
1292 THnSparse *his = calibTime->GetResHistoTPCITS(ihis);
1294 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1295 his->GetAxis(2)->SetRange(0,1000000);
1296 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1297 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,5);
1300 his = calibTime->GetResHistoTPCITS(ihis);
1302 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1303 his->GetAxis(2)->SetRange(0,1000000);
1304 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1305 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,5);
1308 his = calibTime->GetResHistoTPCITS(ihis);
1310 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1311 his->GetAxis(2)->SetRange(0,1000000);
1312 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1313 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,5);
1316 his = calibTime->GetResHistoTPCvertex(ihis);
1318 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1319 his->GetAxis(2)->SetRange(0,1000000);
1320 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1321 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,5);
1324 his = calibTime->GetResHistoTPCvertex(ihis);
1326 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1327 his->GetAxis(2)->SetRange(0,1000000);
1328 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1329 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,5);
1333 his = calibTime->GetResHistoTPCvertex(ihis);
1335 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1336 his->GetAxis(2)->SetRange(0,1000000);
1337 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1338 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,5);
1342 his = calibTime->GetResHistoTPCTOF(ihis);
1344 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1345 his->GetAxis(2)->SetRange(0,1000000);
1346 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
1347 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("TOF%s",hname[ihis]),run,0.,ihis,10);
1351 his = calibTime->GetResHistoTPCTRD(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("TRD%s",hname[ihis]),run,0.,ihis,10);
1363 Double_t AliTPCPreprocessorOffline::EvalAt(Double_t phi, Double_t refX, Double_t theta, Int_t corr, Int_t ptype){
1367 Double_t sector = 9*phi/TMath::Pi();
1368 if (sector<0) sector+=18;
1369 Double_t y85=AliTPCCorrection::GetCorrSector(sector,85,theta,1,corr);
1370 Double_t y245=AliTPCCorrection::GetCorrSector(sector,245,theta,1,corr);
1371 if (ptype==0) return y85+(y245-y85)*(refX-85.)/(245.-85.);
1372 if (ptype==2) return (y245-y85)/(245.-85.);
1378 void AliTPCPreprocessorOffline::MakePrimitivesTime(){
1380 // Create primitive transformation to fit
1382 fAlignTree=new TChain("fit","fit");
1383 fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdy");
1384 fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdsnp");
1385 fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdy");
1386 fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdsnp");
1388 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1389 Double_t bzField=AliTrackerBase::GetBz();
1390 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
1391 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
1392 Double_t wtP = -10.0 * (bzField) * vdrift / ezField ;
1393 AliTPCExBTwist *fitExBTwistX= new AliTPCExBTwist;
1394 AliTPCExBTwist *fitExBTwistY= new AliTPCExBTwist;
1395 AliTPCCalibGlobalMisalignment *trans0 =new AliTPCCalibGlobalMisalignment;
1396 AliTPCCalibGlobalMisalignment *trans1 =new AliTPCCalibGlobalMisalignment;
1397 AliTPCCalibGlobalMisalignment *trans0D =new AliTPCCalibGlobalMisalignment;
1398 AliTPCCalibGlobalMisalignment *trans1D =new AliTPCCalibGlobalMisalignment;
1399 AliTPCCalibGlobalMisalignment *rot0 =new AliTPCCalibGlobalMisalignment;
1400 AliTPCCalibGlobalMisalignment *rot1 =new AliTPCCalibGlobalMisalignment;
1401 AliTPCCalibGlobalMisalignment *rot2 =new AliTPCCalibGlobalMisalignment;
1402 AliTPCCalibGlobalMisalignment *rot3 =new AliTPCCalibGlobalMisalignment;
1405 fitExBTwistX->SetXTwist(0.001);
1406 fitExBTwistX->SetOmegaTauT1T2(wtP,1,1);
1408 fitExBTwistY->SetYTwist(0.001);
1409 fitExBTwistY->SetOmegaTauT1T2(wtP,1,1);
1411 TGeoHMatrix *matrixRot = new TGeoHMatrix;
1412 TGeoHMatrix *matrixX = new TGeoHMatrix;
1413 TGeoHMatrix *matrixY = new TGeoHMatrix;
1414 matrixX->SetDx(0.1);
1415 matrixY->SetDy(0.1);
1416 Double_t rotAngles0[9]={0};
1417 Double_t rotAngles1[9]={0};
1418 Double_t rotAngles2[9]={0};
1420 Double_t rotAngles3[9]={0};
1422 rotAngles0[0]=1; rotAngles0[4]=1; rotAngles0[8]=1;
1423 rotAngles1[0]=1; rotAngles1[4]=1; rotAngles1[8]=1;
1424 rotAngles2[0]=1; rotAngles2[4]=1; rotAngles2[8]=1;
1425 rotAngles3[0]=1; rotAngles3[4]=1; rotAngles3[8]=1;
1427 rotAngles0[1]=-0.001;rotAngles0[3]=0.001;
1428 rotAngles1[5]=-0.001;rotAngles1[7]=0.001;
1429 rotAngles2[2]=0.001;rotAngles2[6]=-0.001;
1430 rotAngles3[1]=0.001;rotAngles3[3]=-0.001;
1431 matrixRot->SetRotation(rotAngles0);
1432 rot0->SetAlignGlobal(matrixRot);
1433 matrixRot->SetRotation(rotAngles1);
1434 rot1->SetAlignGlobal(matrixRot);
1435 matrixRot->SetRotation(rotAngles2);
1436 rot2->SetAlignGlobal(matrixRot);
1437 matrixRot->SetRotation(rotAngles3);
1438 rot3->SetAlignGlobalDelta(matrixRot);
1440 trans0->SetAlignGlobal(matrixX);
1441 trans1->SetAlignGlobal(matrixY);
1442 trans0D->SetAlignGlobalDelta(matrixX);
1443 trans1D->SetAlignGlobalDelta(matrixY);
1444 fitExBTwistX->Init();
1445 fitExBTwistY->Init();
1447 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistX->Clone()),100);
1448 fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistY->Clone()),101);
1450 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot0->Clone()),102);
1451 fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(rot1->Clone()),103);
1452 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot2->Clone()),104);
1453 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot3->Clone()),105);
1455 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0->Clone()),106);
1456 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1->Clone()),107);
1457 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0D->Clone()),108);
1458 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1D->Clone()),109);
1460 fAlignTree->SetAlias("FExBTwistX", "AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,100,ptype)+0");
1461 fAlignTree->SetAlias("FExBTwistY","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,101,ptype)+0");
1462 fAlignTree->SetAlias("FAlignRot0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,102,ptype)+0");
1463 fAlignTree->SetAlias("FAlignRot0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,105,ptype)+0");
1464 fAlignTree->SetAlias("FAlignRot1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,103,ptype)+0");
1465 fAlignTree->SetAlias("FAlignRot2","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,104,ptype)+0");
1466 fAlignTree->SetAlias("FAlignTrans0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,106,ptype)+0");
1467 fAlignTree->SetAlias("FAlignTrans1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,107,ptype)+0");
1468 fAlignTree->SetAlias("FAlignTrans0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,108,ptype)+0");
1469 fAlignTree->SetAlias("FAlignTrans1D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,109,ptype)+0");
1471 // test fast function
1473 // fAlignTree->Draw("FExBTwistX:ExBTwistX","isITS&&ptype==0&&abs(snp)<0.05","");
1474 // fAlignTree->Draw("FExBTwistY:ExBTwistY","isITS&&ptype==0&&abs(snp)<0.05","");
1475 // fAlignTree->Draw("FAlignRot0:alignRot0","isITS&&ptype==0&&abs(snp)<0.05","");
1476 // fAlignTree->Draw("FAlignRot1:alignRot1","isITS&&ptype==0&&abs(snp)<0.05","");
1477 // fAlignTree->Draw("FAlignRot2:alignRot2","isITS&&ptype==0&&abs(snp)<0.05","");
1479 // fAlignTree->Draw("FAlignTrans0:alignTrans0","isITS&&ptype==0&&abs(snp)<0.05","");
1480 // fAlignTree->Draw("FAlignTrans1:alignTrans1","isITS&&ptype==0&&abs(snp)<0.05","");
1485 void AliTPCPreprocessorOffline::CreateAlignTime(TString fstring, TVectorD paramC){
1490 TGeoHMatrix *matrixDelta = new TGeoHMatrix;
1491 TGeoHMatrix *matrixGlobal = new TGeoHMatrix;
1492 Double_t rAngles[9];
1495 index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans0D");
1496 if (index>=0) matrixDelta->SetDx(paramC[index+1]*0.1);
1497 index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans1D");
1498 if (index>=0) matrixDelta->SetDy(paramC[index+1]*0.1);
1499 rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1500 index=TStatToolkit::GetFitIndex(fstring,"FAlignRot0D");
1501 rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1502 rAngles[5]=0; rAngles[7] =0;
1503 rAngles[2]=0; rAngles[6] =0;
1504 matrixDelta->SetRotation(rAngles);
1508 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans0");
1509 if (index>=0) matrixGlobal->SetDx(paramC[index+1]*0.1);
1510 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans1");
1511 if (index>=0) matrixGlobal->SetDy(paramC[index+1]*0.1);
1512 rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1513 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot0");
1514 rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1515 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot1");
1516 rAngles[5]=-paramC[index+1]*0.001; rAngles[7]=paramC[index+1]*0.001;
1517 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot2");
1518 rAngles[2]=paramC[index+1]*0.001; rAngles[6] =-paramC[index+1]*0.001;
1519 matrixGlobal->SetRotation(rAngles);
1521 AliTPCCalibGlobalMisalignment *fitAlignTime =0;
1522 fitAlignTime =new AliTPCCalibGlobalMisalignment;
1523 fitAlignTime->SetName("FitAlignTime");
1524 fitAlignTime->SetTitle("FitAlignTime");
1525 fitAlignTime->SetAlignGlobalDelta(matrixDelta);
1526 fitAlignTime->SetAlignGlobal(matrixGlobal);
1528 AliTPCExBTwist * fitExBTwist= new AliTPCExBTwist;
1529 Int_t indexX=TStatToolkit::GetFitIndex(fstring,"ExBTwistX");
1530 Int_t indexY=TStatToolkit::GetFitIndex(fstring,"ExBTwistY");
1531 fitExBTwist->SetXTwist(0.001*paramC[indexX+1]); // 1 mrad twist in x
1532 fitExBTwist->SetYTwist(0.001*paramC[indexY+1]); // 1 mrad twist in x
1533 fitExBTwist->SetName("FitExBTwistTime");
1534 fitExBTwist->SetTitle("FitExBTwistTime");
1535 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1536 Double_t bzField=AliTrackerBase::GetBz();
1537 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
1539 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
1540 Double_t wt = -10.0 * (bzField) * vdrift / ezField ;
1542 fitExBTwist->SetOmegaTauT1T2(wt,1,1);
1543 fitExBTwist->Init();
1545 AliTPCComposedCorrection *corrTime = new AliTPCComposedCorrection;
1546 TObjArray *arr = new TObjArray;
1547 corrTime->SetCorrections(arr);
1549 corrTime->GetCorrections()->Add(fitExBTwist);
1550 corrTime->GetCorrections()->Add(fitAlignTime);
1551 corrTime->SetName("FitCorrectionTime");
1552 corrTime->SetTitle("FitCorrectionTime");
1554 fitExBTwist->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwist->Clone()),1001);
1555 fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(fitAlignTime->Clone()),1002);
1556 fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(corrTime->Clone()),1003);
1559 fAlignTree->SetAlias("ExBTwistTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1001,ptype)+0");
1560 fAlignTree->SetAlias("AlignTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1002,ptype)+0");
1561 fAlignTree->SetAlias("FitCorrectionTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1003,ptype)+0");
1564 TFile *f = new TFile("fitITSVertex.root","update");
1565 corrTime->Write("FitCorrectionTime");