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 AliTPCPreprocessorOffline process;
40 process.CalibTimeGain("CalibObjects.root",114000,121040,0);
41 TFile oo("OCDB/TPC/Calib/TimeGain/Run114000_121040_v0_s0.root")
42 TObjArray * arr = AliCDBEntry->GetObject()
43 arr->At(4)->Draw("alp")
46 #include "Riostream.h"
49 #include "TGraphErrors.h"
50 #include "AliExternalTrackParam.h"
53 #include "TMultiGraph.h"
55 #include "THnSparse.h"
60 #include "AliTPCROC.h"
61 #include "AliTPCCalROC.h"
62 #include "AliESDfriend.h"
63 #include "AliTPCcalibTime.h"
64 #include "AliSplineFit.h"
65 #include "AliCDBMetaData.h"
67 #include "AliCDBManager.h"
68 #include "AliCDBStorage.h"
69 #include "AliTPCcalibBase.h"
70 #include "AliTPCcalibDB.h"
71 #include "AliTPCcalibDButil.h"
72 #include "AliRelAlignerKalman.h"
73 #include "AliTPCParamSR.h"
74 #include "AliTPCcalibTimeGain.h"
75 #include "AliSplineFit.h"
76 #include "AliTPCPreprocessorOffline.h"
79 ClassImp(AliTPCPreprocessorOffline)
81 AliTPCPreprocessorOffline::AliTPCPreprocessorOffline():
82 TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
83 fMinEntries(500), // minimal number of entries for fit
84 startRun(0), // start Run - used to make fast selection in THnSparse
85 endRun(0), // end Run - used to make fast selection in THnSparse
86 startTime(0), // startTime - used to make fast selection in THnSparse
87 endTime(0), // endTime - used to make fast selection in THnSparse
88 ocdbStorage(""), // path to the OCDB storage
89 fVdriftArray(new TObjArray),
91 fGraphMIP(0), // graph time dependence of MIP
92 fGraphCosmic(0), // graph time dependence at Plateu
93 fGraphAttachmentMIP(0),
94 fFitMIP(0), // fit of dependence - MIP
95 fFitCosmic(0), // fit of dependence - Plateu
96 fGainArray(new TObjArray), // array to be stored in the OCDB
97 fGainMIP(0), // calibration component for MIP
98 fGainCosmic(0), // calibration component for cosmic
99 fSwitchOnValidation(kFALSE) // flag to switch on validation of OCDB parameters
102 // default constructor
106 AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
115 void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const timeDrift){
117 // find the fist and last run
119 TObjArray *hisArray =timeDrift->GetHistoDrift();
120 {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
121 THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
122 if (addHist->GetEntries()<fMinEntries) continue;
123 if (!addHist) continue;
124 TH1D* histo =addHist->Projection(3);
125 TH1D* histoTime=addHist->Projection(0);
126 printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
129 startRun=histo->FindFirstBinAbove(0);
130 endRun =histo->FindLastBinAbove(0);
132 startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun);
133 endRun =TMath::Max(histo->FindLastBinAbove(0),endRun);
136 startTime=histoTime->FindFirstBinAbove(0);
137 endTime =histoTime->FindLastBinAbove(0);
139 startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime);
140 endTime =TMath::Max(histoTime->FindLastBinAbove(0),endTime);
145 if (startRun<0) startRun=0;
146 if (endRun<0) endRun=100000000;
147 printf("Run range :\t%d-%d\n", startRun, endRun);
148 printf("Time range :\t%d-%d\n", startTime, endTime);
154 void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, TString pocdbStorage){
156 // make calibration of the drift velocity
158 // file - the location of input file
159 // ustartRun, uendrun - run validity period
160 // pocdbStorage - path to hte OCDB storage
161 // - if empty - local storage 'pwd' uesed
162 if (pocdbStorage.Length()>0) ocdbStorage=pocdbStorage;
164 ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
166 // 1. Initialization and run range setting
168 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
170 fTimeDrift = (AliTPCcalibTime *)array->FindObject("calibTime");
172 fTimeDrift = (AliTPCcalibTime*)fcalib.Get("calibTime");
174 if(!fTimeDrift) return;
178 TObjArray *hisArray =fTimeDrift->GetHistoDrift();
179 GetRunRange(fTimeDrift);
180 for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
181 THnSparse* addHist=(THnSparse*)hisArray->At(i);
182 if (!addHist) continue;
183 if (startTime<endTime) addHist->GetAxis(0)->SetRange(startTime-1,endTime+1);
184 if (startRun<endRun) addHist->GetAxis(3)->SetRange(startRun-1,endRun+1);
188 // 2. extraction of the information
190 fVdriftArray = new TObjArray();
191 AddAlignmentGraphs(fVdriftArray,fTimeDrift);
192 AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
193 AddLaserGraphs(fVdriftArray,fTimeDrift);
195 // 3. Append QA plots
197 MakeDefaultPlots(fVdriftArray,fVdriftArray);
200 // 4. validate OCDB entries
202 if(fSwitchOnValidation==kTRUE && ValidateTimeDrift()==kFALSE) {
203 Printf("TPC time drift OCDB parameters out of range!");
212 UpdateOCDBDrift(ustartRun,uendRun,ocdbStorage);
215 void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun, const char* storagePath ){
219 AliCDBMetaData *metaData= new AliCDBMetaData();
220 metaData->SetObjectClassName("TObjArray");
221 metaData->SetResponsible("Marian Ivanov");
222 metaData->SetBeamPeriod(1);
223 metaData->SetAliRootVersion("05-25-01"); //root version
224 metaData->SetComment("Calibration of the time dependence of the drift velocity");
226 id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
227 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
228 gStorage->Put(fVdriftArray, (*id1), metaData);
231 Bool_t AliTPCPreprocessorOffline::ValidateTimeGain(Double_t minGain, Double_t maxGain)
234 // Validate time gain corrections
236 Printf("ValidateTimeGain..." );
238 TGraphErrors *gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
239 if(!gr) return kFALSE;
240 if(gr->GetN()<1) return kFALSE;
242 // check whether gain in the range
243 for(Int_t iPoint=0; iPoint<gr->GetN(); iPoint++)
245 if(gr->GetY()[iPoint] < minGain || gr->GetY()[iPoint] > maxGain)
253 Bool_t AliTPCPreprocessorOffline::ValidateTimeDrift(Double_t maxVDriftCorr)
256 // Validate time drift velocity corrections
258 Printf("ValidateTimeDrift..." );
260 TGraphErrors* gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
261 if(!gr) return kFALSE;
262 if(gr->GetN()<1) return kFALSE;
264 // check whether drift velocity corrections in the range
265 for(Int_t iPoint = 0; iPoint<gr->GetN(); iPoint++)
267 if(TMath::Abs(gr->GetY()[iPoint]) > maxVDriftCorr)
274 void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
276 // update the OCDB entry for the nominal time0
279 // AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
280 AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
281 TGraphErrors *grT = (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
282 Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
283 Double_t deltaT = deltaTcm/param->GetDriftV();
284 paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
287 AliCDBMetaData *metaData= new AliCDBMetaData();
288 metaData->SetObjectClassName("TObjArray");
289 metaData->SetResponsible("Marian Ivanov");
290 metaData->SetBeamPeriod(1);
291 metaData->SetAliRootVersion("05-25-02"); //root version
292 metaData->SetComment("Updated calibration of nominal time 0");
294 id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
295 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
296 gStorage->Put(param, (*id1), metaData);
301 void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
303 // Print the names of the entries in array
305 Int_t entries = array->GetEntries();
306 for (Int_t i=0; i<entries; i++){
307 if (!array->At(i)) continue;
308 printf("%d\t %s\n", i, array->At(i)->GetName());
314 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
316 // 1. filter graph - error cut errSigmaCut
317 // 2. filter graph - medianCutAbs around median
319 // errSigmaCut - cut on error
320 // medianCutAbs - cut on value around median
323 // 1. filter graph - error cut errSigmaCut
325 TGraphErrors *graphF;
326 graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
328 if (!graphF) return 0;
329 graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
331 if (!graph) return 0;
333 // filter graph - kMedianCutAbs around median
335 graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
337 if (!graphF) return 0;
338 graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
340 if (!graph) return 0;
346 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
348 // filter outlyer measurement
349 // Only points around median +- cut filtered
351 if (!graph) return 0;
353 Int_t npoints0 = graph->GetN();
356 Double_t *outx=new Double_t[npoints0];
357 Double_t *outy=new Double_t[npoints0];
358 Double_t *errx=new Double_t[npoints0];
359 Double_t *erry=new Double_t[npoints0];
362 if (npoints0<kMinPoints) return 0;
363 for (Int_t iter=0; iter<3; iter++){
365 for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
366 if (graph->GetY()[ipoint]==0) continue;
367 if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;
368 outx[npoints] = graph->GetX()[ipoint];
369 outy[npoints] = graph->GetY()[ipoint];
370 errx[npoints] = graph->GetErrorX(ipoint);
371 erry[npoints] = graph->GetErrorY(ipoint);
374 if (npoints<=1) break;
375 medianY =TMath::Median(npoints,outy);
376 rmsY =TMath::RMS(npoints,outy);
378 TGraphErrors *graphOut=0;
379 if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry);
384 void AliTPCPreprocessorOffline::AddHistoGraphs( TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
386 // Add graphs corresponding to the alignment
388 const Double_t kErrSigmaCut=5; // error sigma cut - for filtering
389 const Double_t kMedianCutAbs=0.03; // error sigma cut - for filtering
391 TObjArray * array=timeDrift->GetHistoDrift();
393 THnSparse* hist=NULL;
394 // 2.a) cosmics with different triggers
395 for (Int_t i=0; i<array->GetEntriesFast();i++){
396 hist=(THnSparseF*)array->UncheckedAt(i);
398 if (hist->GetEntries()<minEntries) continue;
400 TString name=hist->GetName();
401 Int_t dim[4]={0,1,2,3};
402 THnSparse* newHist=hist->Projection(4,dim);
403 newHist->SetName(name);
404 TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
405 printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
406 Int_t pos=name.Index("_");
407 name=name(pos,name.Capacity()-pos);
408 TString graphName=graph->ClassName();
412 graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
414 printf("Graph =%s filtered out\n", name.Data());
418 graph->SetMarkerStyle(i%8+20);
419 graph->SetMarkerColor(i%7);
420 graph->GetXaxis()->SetTitle("Time");
421 graph->GetYaxis()->SetTitle("v_{dcor}");
422 graph->SetName(graphName);
423 graph->SetTitle(graphName);
424 printf("Graph %d\t=\t%s\n", i, graphName.Data());
425 vdriftArray->Add(graph);
433 void AliTPCPreprocessorOffline::AddAlignmentGraphs( TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
435 // Add graphs corresponding to alignment to the object array
437 TObjArray *arrayITS=0;
438 TObjArray *arrayTOF=0;
439 TObjArray *arrayTRD=0;
440 TMatrixD *mstatITS=0;
441 TMatrixD *mstatTOF=0;
442 TMatrixD *mstatTRD=0;
444 arrayITS=timeDrift->GetAlignITSTPC();
445 arrayTRD=timeDrift->GetAlignTRDTPC();
446 arrayTOF=timeDrift->GetAlignTOFTPC();
448 if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.9,50,0.025);
449 if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.9,1000,0.025);
450 if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.9,50,0.025);
452 TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
453 TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
454 TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
455 TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
456 TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
457 TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
459 TObjArray * arrayTRDP= 0x0;
460 TObjArray * arrayTRDM= 0x0;
461 TObjArray * arrayTRDB= 0x0;
462 arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
463 arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
464 arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
467 Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
468 TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
469 arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
470 arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
471 TString grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
472 "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
473 "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
474 TString grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
475 "DELTAX", "DELTAY", "DELTAZ",
476 "DRIFTVD", "T0", "VDGY"};
479 TVectorD vX(entries);
480 TVectorD vY(entries);
481 TVectorD vEx(entries);
482 TVectorD vEy(entries);
484 for (Int_t iarray=0; iarray<12; iarray++){
485 arr = arrays[iarray];
486 if (arr==0) continue;
487 for (Int_t ipar=0; ipar<9; ipar++){
489 for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
490 AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
491 if (!kalman) continue;
492 vX[counter]=kalman->GetTimeStamp();
493 vY[counter]=(*(kalman->GetState()))[ipar];
494 if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
496 vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
500 TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
502 vEx.GetMatrixArray(),
503 vEy.GetMatrixArray());
504 TString grName=grnames[iarray];
507 graph->SetName(grName.Data());
508 vdriftArray->AddLast(graph);
516 void AliTPCPreprocessorOffline::AddLaserGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
518 // add graphs for laser
520 const Double_t delayL0L1 = 0.071; //this is hack for 1/2 weeks
522 TGraphErrors *grLaser[6]={0,0,0,0,0,0};
523 hisN = timeDrift->GetHistVdriftLaserA(0);
524 if (timeDrift->GetHistVdriftLaserA(0)){
525 grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
526 grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
527 vdriftArray->AddLast(grLaser[0]);
529 if (timeDrift->GetHistVdriftLaserA(1)){
530 grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
531 grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
532 vdriftArray->AddLast(grLaser[1]);
534 if (timeDrift->GetHistVdriftLaserA(2)){
535 grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
536 grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
537 vdriftArray->AddLast(grLaser[2]);
539 if (timeDrift->GetHistVdriftLaserC(0)){
540 grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
541 grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
542 vdriftArray->AddLast(grLaser[3]);
544 if (timeDrift->GetHistVdriftLaserC(1)){
545 grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
546 grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
547 vdriftArray->AddLast(grLaser[4]);
549 if (timeDrift->GetHistVdriftLaserC(2)){
550 grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
551 grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
552 vdriftArray->AddLast(grLaser[5]);
554 for (Int_t i=0; i<6;i++){
556 SetDefaultGraphDrift(grLaser[i], 1,(i+20));
557 grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
563 TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
565 // Make graph with mean values and rms
567 hisN->GetAxis(itime)->SetRange(0,100000000);
568 hisN->GetAxis(ival)->SetRange(0,100000000);
569 TH1 * hisT = hisN->Projection(itime);
570 TH1 * hisV = hisN->Projection(ival);
572 Int_t firstBinA = hisT->FindFirstBinAbove(2);
573 Int_t lastBinA = hisT->FindLastBinAbove(2);
574 Int_t firstBinV = hisV->FindFirstBinAbove(0);
575 Int_t lastBinV = hisV->FindLastBinAbove(0);
576 hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
577 hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
579 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
580 Double_t cont = hisT->GetBinContent(ibin);
581 if (cont<minEntries) continue;
584 TVectorD vecTime(entries);
585 TVectorD vecMean0(entries);
586 TVectorD vecRMS0(entries);
587 TVectorD vecMean1(entries);
588 TVectorD vecRMS1(entries);
590 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
591 Double_t cont = hisT->GetBinContent(ibin);
592 if (cont<minEntries) continue;
593 //hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
594 Int_t minBin = ibin-1;
595 Int_t maxBin = ibin+1;
596 if(minBin <= 0) minBin = 1;
597 if(maxBin >= hisN->GetAxis(itime)->GetNbins()) maxBin = hisN->GetAxis(itime)->GetNbins()-1;
598 hisN->GetAxis(itime)->SetRange(minBin,maxBin);
600 Double_t time = hisT->GetBinCenter(ibin);
601 TH1 * his = hisN->Projection(ival);
602 Double_t nentries0= his->GetBinContent(his->FindBin(0));
603 if (cont-nentries0<minEntries) continue;
605 his->SetBinContent(his->FindBin(0),0);
606 vecTime[entries]=time;
607 vecMean0[entries]=his->GetMean()+offset;
608 vecMean1[entries]=his->GetMeanError();
609 vecRMS0[entries] =his->GetRMS();
610 vecRMS1[entries] =his->GetRMSError();
616 TGraphErrors * graph = new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(), 0, vecMean1.GetMatrixArray());
627 void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
629 // Set default style for QA views
631 graph->GetXaxis()->SetTimeDisplay(kTRUE);
632 graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
633 graph->SetMaximum( 0.025);
634 graph->SetMinimum(-0.025);
635 graph->GetXaxis()->SetTitle("Time");
636 graph->GetYaxis()->SetTitle("v_{dcorr}");
638 graph->GetYaxis()->SetLabelSize(0.03);
639 graph->GetXaxis()->SetLabelSize(0.03);
641 graph->GetXaxis()->SetNdivisions(10,5,0);
642 graph->GetYaxis()->SetNdivisions(10,5,0);
644 graph->GetXaxis()->SetLabelOffset(0.02);
645 graph->GetYaxis()->SetLabelOffset(0.005);
647 graph->GetXaxis()->SetTitleOffset(1.3);
648 graph->GetYaxis()->SetTitleOffset(1.2);
650 graph->SetMarkerColor(color);
651 graph->SetLineColor(color);
652 graph->SetMarkerStyle(style);
655 void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
657 // Set default pad style for QA
660 pad->SetMargin(mx0,mx1,my0,my1);
664 void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray * /*picArray*/){
666 // 0. make a default QA plots
667 // 1. Store them in the array
670 Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
672 TGraphErrors* laserA =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
673 TGraphErrors* laserC =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
674 TGraphErrors* cosmic =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
675 TGraphErrors* cross =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
676 TGraphErrors* itstpcP =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
677 TGraphErrors* itstpcM =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
678 TGraphErrors* itstpcB =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
680 if (laserA) SetDefaultGraphDrift(laserA,2,25);
681 if (laserC) SetDefaultGraphDrift(laserC,4,26);
682 if (cosmic) SetDefaultGraphDrift(cosmic,3,27);
683 if (cross) SetDefaultGraphDrift(cross,4,28);
684 if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
685 if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
686 if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
694 pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
696 SetPadStyle(pad,mx0,mx1,my0,my1);
699 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
700 legend->AddEntry(laserA,"Laser A side");
701 legend->AddEntry(laserC,"Laser C side");
703 //picArray->AddLast(pad);
706 if (itstpcP&&itstpcM){
707 pad = new TCanvas("ITSTPC","ITSTPC");
708 itstpcP->Draw("alp");
709 SetPadStyle(pad,mx0,mx1,my0,my1);
710 itstpcP->Draw("alp");
712 itstpcM->Draw("apl");
715 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
716 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
717 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
718 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
720 //picArray->AddLast(pad);
723 if (itstpcB&&laserA){
724 pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
725 SetPadStyle(pad,mx0,mx1,my0,my1);
730 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
731 legend->AddEntry(laserA,"TPC laser");
732 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
733 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
734 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
736 //picArray->AddLast(pad);
740 pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
741 SetPadStyle(pad,mx0,mx1,my0,my1);
742 itstpcP->Draw("alp");
747 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
749 legend->AddEntry(cross,"TPC cross tracks");
750 legend->AddEntry(itstpcB,"ITS-TPC smooth");
752 //picArray->AddLast(pad);
754 if (itstpcP&&cosmic){
755 pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
756 SetPadStyle(pad,mx0,mx1,my0,my1);
757 itstpcP->Draw("alp");
762 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
764 legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
765 legend->AddEntry(itstpcB,"ITS-TPC smooth");
767 //picArray->AddLast(pad);
774 void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber, TString pocdbStorage){
778 if (pocdbStorage.Length()==0) pocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
781 // 1. Read gain values
783 ReadGainGlobal(fileName);
786 // 2. Extract calibration values
788 AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
789 AnalyzeAttachment(startRunNumber,endRunNumber);
792 // 3. Make control plots
797 // 4. validate OCDB entries
799 if(fSwitchOnValidation==kTRUE && ValidateTimeGain()==kFALSE) {
800 Printf("TPC time gain OCDB parameters out of range!");
807 UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage.Data());
810 void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
812 // read calibration entries from file
814 TFile fcalib(fileName);
815 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
817 fGainMIP = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
818 fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
820 fGainMIP = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
821 fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
824 Int_t firstBinA =0, lastBinA=0;
827 hisT= fGainCosmic->GetHistGainTime()->Projection(1);
828 firstBinA = hisT->FindFirstBinAbove(2);
829 lastBinA = hisT->FindLastBinAbove(2);
830 fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
835 hisT= fGainMIP->GetHistGainTime()->Projection(1);
836 firstBinA = hisT->FindFirstBinAbove(2);
837 lastBinA = hisT->FindLastBinAbove(2);
838 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
846 Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit, Float_t FPtoMIPratio){
848 // Analyze gain - produce the calibration graphs
851 // 1.) try to create MIP spline
854 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
855 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
856 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
858 fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
859 if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
860 if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
861 if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
862 fGainArray->AddAt(fFitMIP,0);
865 // 2.) try to create Cosmic spline
868 fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
869 fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
871 fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
872 if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
875 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
876 fGraphCosmic->GetY()[i] /= FPtoMIPratio;
877 fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
881 if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
882 if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
883 fGainArray->AddAt(fFitCosmic,1);
885 // with naming convention and backward compatibility
886 fGainArray->AddAt(fGraphMIP,2);
887 fGainArray->AddAt(fGraphCosmic,3);
888 cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
893 Bool_t AliTPCPreprocessorOffline::AnalyzeAttachment(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesFit) {
895 // determine slope as a function of mean driftlength
897 if(!fGainMIP) return kFALSE;
899 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
901 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
902 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
904 fGainMIP->GetHistGainTime()->GetAxis(3)->SetRangeUser(125,250);// only full tracking region (driftlength)
905 fGainMIP->GetHistGainTime()->GetAxis(0)->SetRangeUser(1.5,3.5);// only full tracking region (driftlength)
907 TH3D * hist = fGainMIP->GetHistGainTime()->Projection(1, 0, 3);
909 Double_t *xvec = new Double_t[hist->GetNbinsX()];
910 Double_t *yvec = new Double_t[hist->GetNbinsX()];
911 Double_t *xerr = new Double_t[hist->GetNbinsX()];
912 Double_t *yerr = new Double_t[hist->GetNbinsX()];
915 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
919 for (Int_t idelta=0; idelta<5; idelta++){
921 imin = TMath::Max(i-idelta,1);
922 imax = TMath::Min(i+idelta,hist->GetNbinsX());
923 nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
924 //if (nsum==0) break;
925 if (nsum>minEntriesFit) break;
927 if (nsum<minEntriesFit) continue;
929 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRangeUser(hist->GetXaxis()->GetBinCenter(imin-1),hist->GetXaxis()->GetBinCenter(imax+1)); // define time range
930 TH2D * histZdep = fGainMIP->GetHistGainTime()->Projection(0,3);
932 histZdep->FitSlicesY(0,0,-1,0,"QNR",&arr);
933 TH1D * driftDep = (TH1D*)arr.At(1);
935 //TGraphErrors * driftDep = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,3,100,1,0.,1);
936 /*if (driftDep->GetN() < 4) {
943 TF1 pol1("polynom1","pol1",125,240);
944 //driftDep->Fit(&pol1,"QNRROB=0.8");
945 driftDep->Fit(&pol1,"QNR");
946 xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin-1)+hist->GetXaxis()->GetBinCenter(imax+1));
947 yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
948 xerr[counter] = hist->GetXaxis()->GetBinCenter(imax+1)-hist->GetXaxis()->GetBinCenter(imin-1);
949 yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
955 fGraphAttachmentMIP = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
956 if (fGraphAttachmentMIP) fGraphAttachmentMIP->SetName("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");// set proper names according to naming convention
957 fGainArray->AddLast(fGraphAttachmentMIP);
965 if (counter < 1) return kFALSE;
971 void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
975 AliCDBMetaData *metaData= new AliCDBMetaData();
976 metaData->SetObjectClassName("TObjArray");
977 metaData->SetResponsible("Alexander Kalweit");
978 metaData->SetBeamPeriod(1);
979 metaData->SetAliRootVersion("05-24-00"); //root version
980 metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
981 AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
982 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
983 gStorage->Put(fGainArray, id1, metaData);
986 void AliTPCPreprocessorOffline::MakeQAPlot(Float_t FPtoMIPratio) {
988 // Make QA plot to visualize results
993 TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
995 TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
996 gainHistoCosmic->SetDirectory(0);
997 gainHistoCosmic->SetName("GainHistoCosmic");
998 gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
999 gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1000 gainHistoCosmic->Draw("colz");
1001 fGraphCosmic->SetMarkerStyle(25);
1002 fGraphCosmic->Draw("lp");
1003 fGraphCosmic->SetMarkerStyle(25);
1004 TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
1006 for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
1007 grfFitCosmic->GetY()[i] *= FPtoMIPratio;
1009 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
1010 fGraphCosmic->GetY()[i] *= FPtoMIPratio;
1013 fGraphCosmic->Draw("lp");
1014 grfFitCosmic->SetLineColor(2);
1015 grfFitCosmic->Draw("lu");
1016 fGainArray->AddLast(gainHistoCosmic);
1017 //fGainArray->AddLast(canvasCosmic->Clone());
1018 delete canvasCosmic;
1021 TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
1023 TH2D * gainHistoMIP = fGainMIP->GetHistGainTime()->Projection(0,1);
1024 gainHistoMIP->SetName("GainHistoCosmic");
1025 gainHistoMIP->SetDirectory(0);
1026 gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
1027 gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1028 gainHistoMIP->Draw("colz");
1029 fGraphMIP->SetMarkerStyle(25);
1030 fGraphMIP->Draw("lp");
1031 TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
1032 grfFitMIP->SetLineColor(2);
1033 grfFitMIP->Draw("lu");
1034 fGainArray->AddLast(gainHistoMIP);
1035 //fGainArray->AddLast(canvasMIP->Clone());