Make slewing correction independent of the TDC threshold. Checked by comparing 900GeV...
[u/mrichter/AliRoot.git] / TPC / AliTPCPreprocessorOffline.cxx
CommitLineData
da268e11 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16
17
18/*
f12abd89 19 Responsible: marian.ivanov@cern.ch
20 Code to analyze the TPC calibration and to produce OCDB entries
21
22
da268e11 23 .x ~/rootlogon.C
24 gSystem->Load("libANALYSIS");
25 gSystem->Load("libTPCcalib");
da268e11 26
27 AliTPCPreprocessorOffline proces;
f12abd89 28 process.CalibTimeGain("CalibObjects.root",run0,run1,ocdbPath);
29 // take the raw calibration data from the file CalibObjects.root
30 // and make a OCDB entry with run validity run0-run1
31 // results are stored at the ocdbPath - local or alien ...
32 // default storage ""- data stored at current working directory
da268e11 33
34*/
f12abd89 35#include "Riostream.h"
36#include <fstream>
da268e11 37#include "TMap.h"
38#include "TGraphErrors.h"
39#include "AliExternalTrackParam.h"
40#include "TFile.h"
41#include "TGraph.h"
42#include "TMultiGraph.h"
43#include "TCanvas.h"
44#include "THnSparse.h"
45#include "TLegend.h"
46#include "TPad.h"
47#include "TH2D.h"
48#include "AliTPCROC.h"
49#include "AliTPCCalROC.h"
da268e11 50#include "AliESDfriend.h"
da268e11 51#include "AliTPCcalibTime.h"
52#include "AliSplineFit.h"
53#include "AliCDBMetaData.h"
54#include "AliCDBId.h"
55#include "AliCDBManager.h"
56#include "AliCDBStorage.h"
57#include "AliTPCcalibBase.h"
58#include "AliTPCcalibDB.h"
59#include "AliTPCcalibDButil.h"
60#include "AliRelAlignerKalman.h"
61#include "AliTPCParamSR.h"
62#include "AliTPCcalibTimeGain.h"
63#include "AliSplineFit.h"
64#include "AliTPCPreprocessorOffline.h"
da268e11 65
66
67ClassImp(AliTPCPreprocessorOffline)
68
69AliTPCPreprocessorOffline::AliTPCPreprocessorOffline():
f12abd89 70 TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
40120adc 71 fMinEntries(500), // minimal number of entries for fit
da268e11 72 startRun(0), // start Run - used to make fast selection in THnSparse
73 endRun(0), // end Run - used to make fast selection in THnSparse
74 startTime(0), // startTime - used to make fast selection in THnSparse
75 endTime(0), // endTime - used to make fast selection in THnSparse
76 ocdbStorage(""), // path to the OCDB storage
77 fVdriftArray(new TObjArray),
78 fTimeDrift(0),
79 fGraphMIP(0), // graph time dependence of MIP
80 fGraphCosmic(0), // graph time dependence at Plateu
81 fFitMIP(0), // fit of dependence - MIP
82 fFitCosmic(0), // fit of dependence - Plateu
83 fGainArray(new TObjArray), // array to be stored in the OCDB
84 fGainMIP(0), // calibration component for MIP
85 fGainCosmic(0) // calibration component for cosmic
86{
87 //
f12abd89 88 // default constructor
89 //
da268e11 90}
91
92AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
93 //
94 // Destructor
95 //
96}
97
98
99
100
40120adc 101void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const timeDrift){
da268e11 102 //
103 // find the fist and last run
104 //
f12abd89 105 TObjArray *hisArray =timeDrift->GetHistoDrift();
da268e11 106 {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
107 THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
40120adc 108 if (addHist->GetEntries()<fMinEntries) continue;
da268e11 109 if (!addHist) continue;
110 TH1D* histo =addHist->Projection(3);
111 TH1D* histoTime=addHist->Projection(0);
112 printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
113
114 if (startRun<=0){
115 startRun=histo->FindFirstBinAbove(0);
116 endRun =histo->FindLastBinAbove(0);
117 }else{
118 startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun);
119 endRun =TMath::Max(histo->FindLastBinAbove(0),endRun);
120 }
121 if (startTime==0){
122 startTime=histoTime->FindFirstBinAbove(0);
123 endTime =histoTime->FindLastBinAbove(0);
124 }else{
125 startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime);
126 endTime =TMath::Max(histoTime->FindLastBinAbove(0),endTime);
127 }
128 delete histo;
129 delete histoTime;
130 }}
131 if (startRun<0) startRun=0;
132 if (endRun<0) endRun=100000000;
133 printf("Run range :\t%d-%d\n", startRun, endRun);
134 printf("Time range :\t%d-%d\n", startTime, endTime);
135
136}
137
138
139
40120adc 140void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, TString pocdbStorage){
da268e11 141 //
40120adc 142 // make calibration of the drift velocity
143 // Input parameters:
144 // file - the location of input file
145 // ustartRun, uendrun - run validity period
146 // pocdbStorage - path to hte OCDB storage
147 // - if empty - local storage 'pwd' uesed
da268e11 148 if (pocdbStorage.Length()>0) ocdbStorage=pocdbStorage;
149 else
150 ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
151 //
152 // 1. Initialization and run range setting
153 TFile fcalib(file);
154 fTimeDrift=(AliTPCcalibTime*)fcalib.Get("calibTime");
155 startRun=ustartRun;
156 endRun=ustartRun;
157 TObjArray *hisArray =fTimeDrift->GetHistoDrift();
158 GetRunRange(fTimeDrift);
159 for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
160 THnSparse* addHist=(THnSparse*)hisArray->At(i);
161 if (!addHist) continue;
162 if (startTime<endTime) addHist->GetAxis(0)->SetRange(startTime-1,endTime+1);
163 if (startRun<endRun) addHist->GetAxis(3)->SetRange(startRun-1,endRun+1);
164 }
165 //
166 //
167 // 2. extraction of the information
168 //
169 fVdriftArray = new TObjArray();
170 AddAlignmentGraphs(fVdriftArray,fTimeDrift);
40120adc 171 AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
da268e11 172 AddLaserGraphs(fVdriftArray,fTimeDrift);
173 //
174 // 3. Append QA plots
175 //
176 MakeDefaultPlots(fVdriftArray,fVdriftArray);
177 //
178 //
179 // 4. update of OCDB
180 //
181 //
182
183 UpdateOCDBDrift(ustartRun,uendRun,ocdbStorage);
184}
185
186void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun, const char* storagePath ){
187 //
188 // Update OCDB
189 //
190 AliCDBMetaData *metaData= new AliCDBMetaData();
191 metaData->SetObjectClassName("TObjArray");
192 metaData->SetResponsible("Marian Ivanov");
193 metaData->SetBeamPeriod(1);
194 metaData->SetAliRootVersion("05-25-01"); //root version
195 metaData->SetComment("Calibration of the time dependence of the drift velocity");
196 AliCDBId* id1=NULL;
197 id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
198 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
199 gStorage->Put(fVdriftArray, (*id1), metaData);
200}
201
202
203
40120adc 204void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
da268e11 205 //
206 // update the OCDB entry for the nominal time0
207 //
208 //
209 // AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
210 AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
211 TGraphErrors *grT = (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
212 Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
213 Double_t deltaT = deltaTcm/param->GetDriftV();
214 paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
215 paramNew->Update();
216
217 AliCDBMetaData *metaData= new AliCDBMetaData();
218 metaData->SetObjectClassName("TObjArray");
219 metaData->SetResponsible("Marian Ivanov");
220 metaData->SetBeamPeriod(1);
221 metaData->SetAliRootVersion("05-25-02"); //root version
222 metaData->SetComment("Updated calibration of nominal time 0");
223 AliCDBId* id1=NULL;
f12abd89 224 id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
da268e11 225 AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
226 gStorage->Put(param, (*id1), metaData);
227
228}
229
230
231void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
232 //
40120adc 233 // Print the names of the entries in array
da268e11 234 //
235 Int_t entries = array->GetEntries();
236 for (Int_t i=0; i<entries; i++){
237 if (!array->At(i)) continue;
238 printf("%d\t %s\n", i, array->At(i)->GetName());
239 }
240}
241
242
243
244TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
245 // 2 filters:
246 // 1. filter graph - error cut errSigmaCut
247 // 2. filter graph - medianCutAbs around median
248 //
249 // errSigmaCut - cut on error
250 // medianCutAbs - cut on value around median
251 Double_t dummy=0; //
252 //
253 // 1. filter graph - error cut errSigmaCut
254 //
255 TGraphErrors *graphF;
256 graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
257 delete graph;
258 if (!graphF) return 0;
259 graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
260 delete graphF;
261 if (!graph) return 0;
262 //
263 // filter graph - kMedianCutAbs around median
264 //
265 graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
266 delete graph;
267 if (!graphF) return 0;
268 graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
269 delete graphF;
270 if (!graph) return 0;
271 return graph;
272}
273
274
275
276TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
277 //
278 // filter outlyer measurement
279 // Only points around median +- cut filtered
280 //
281 if (!graph) return 0;
282 Int_t kMinPoints=2;
283 Int_t npoints0 = graph->GetN();
284 Int_t npoints=0;
285 Float_t rmsY=0;
286 Double_t *outx=new Double_t[npoints0];
287 Double_t *outy=new Double_t[npoints0];
288 Double_t *errx=new Double_t[npoints0];
289 Double_t *erry=new Double_t[npoints0];
290 //
291 //
292 if (npoints0<kMinPoints) return 0;
293 for (Int_t iter=0; iter<3; iter++){
294 npoints=0;
295 for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
296 if (graph->GetY()[ipoint]==0) continue;
297 if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;
298 outx[npoints] = graph->GetX()[ipoint];
299 outy[npoints] = graph->GetY()[ipoint];
300 errx[npoints] = graph->GetErrorX(ipoint);
301 erry[npoints] = graph->GetErrorY(ipoint);
302 npoints++;
303 }
304 if (npoints<=1) break;
305 medianY =TMath::Median(npoints,outy);
306 rmsY =TMath::RMS(npoints,outy);
307 }
308 TGraphErrors *graphOut=0;
309 if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry);
310 return graphOut;
311}
312
313
40120adc 314void AliTPCPreprocessorOffline::AddHistoGraphs( TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
da268e11 315 //
316 // Add graphs corresponding to the alignment
317 //
318 const Double_t kErrSigmaCut=5; // error sigma cut - for filtering
319 const Double_t kMedianCutAbs=0.03; // error sigma cut - for filtering
320 //
321 TObjArray * array=timeDrift->GetHistoDrift();
322 if (array){
323 THnSparse* hist=NULL;
324 // 2.a) cosmics with different triggers
325 for (Int_t i=0; i<array->GetEntriesFast();i++){
326 hist=(THnSparseF*)array->UncheckedAt(i);
327 if(!hist) continue;
328 if (hist->GetEntries()<minEntries) continue;
329 //hist->Print();
330 TString name=hist->GetName();
331 Int_t dim[4]={0,1,2,3};
332 THnSparse* newHist=hist->Projection(4,dim);
333 newHist->SetName(name);
334 TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
335 printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
336 Int_t pos=name.Index("_");
337 name=name(pos,name.Capacity()-pos);
338 TString graphName=graph->ClassName();
339 graphName+=name;
340 graphName.ToUpper();
341 //
342 graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
343 if (!graph) {
344 printf("Graph =%s filtered out\n", name.Data());
345 continue;
346 }
347 //
348 graph->SetMarkerStyle(i%8+20);
349 graph->SetMarkerColor(i%7);
350 graph->GetXaxis()->SetTitle("Time");
351 graph->GetYaxis()->SetTitle("v_{dcor}");
352 graph->SetName(graphName);
353 graph->SetTitle(graphName);
354 printf("Graph %d\t=\t%s\n", i, graphName.Data());
355 vdriftArray->Add(graph);
356 }
357 }
358}
359
360
361
362
40120adc 363void AliTPCPreprocessorOffline::AddAlignmentGraphs( TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
da268e11 364 //
365 // Add graphs corresponding to alignment to the object array
366 //
367 TObjArray *arrayITS=0;
368 TObjArray *arrayTOF=0;
369 TObjArray *arrayTRD=0;
370 TMatrixD *mstatITS=0;
371 TMatrixD *mstatTOF=0;
372 TMatrixD *mstatTRD=0;
373 //
374 arrayITS=timeDrift->GetAlignITSTPC();
375 arrayTRD=timeDrift->GetAlignTRDTPC();
376 arrayTOF=timeDrift->GetAlignTOFTPC();
377
378 if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.9,50,0.025);
379 if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.9,1000,0.025);
380 if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.9,50,0.025);
381 //
382 TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
383 TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
384 TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
385 TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
386 TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
387 TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
388
389 TObjArray * arrayTRDP= 0x0;
390 TObjArray * arrayTRDM= 0x0;
391 TObjArray * arrayTRDB= 0x0;
392 arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
393 arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
394 arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
395 //
396 //
397 Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
398 TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
399 arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
400 arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
401 TString grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
402 "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
403 "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
404 TString grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
405 "DELTAX", "DELTAY", "DELTAZ",
406 "DRIFTVD", "T0", "VDGY"};
407
408
409 TVectorD vX(entries);
410 TVectorD vY(entries);
411 TVectorD vEx(entries);
412 TVectorD vEy(entries);
413 TObjArray *arr=0;
414 for (Int_t iarray=0; iarray<12; iarray++){
415 arr = arrays[iarray];
416 if (arr==0) continue;
417 for (Int_t ipar=0; ipar<9; ipar++){
418 Int_t counter=0;
419 for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
420 AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
421 if (!kalman) continue;
422 vX[counter]=kalman->GetTimeStamp();
423 vY[counter]=(*(kalman->GetState()))[ipar];
424 if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
425 vEx[counter]=0;
426 vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
427 counter++;
428 }
429
430 TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
431 vY.GetMatrixArray(),
432 vEx.GetMatrixArray(),
433 vEy.GetMatrixArray());
434 TString grName=grnames[iarray];
435 grName+="_TPC_";
436 grName+=grpar[ipar];
437 graph->SetName(grName.Data());
438 vdriftArray->AddLast(graph);
439 }
440 }
441}
442
443
444
445
446void AliTPCPreprocessorOffline::AddLaserGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
447 //
448 // add graphs for laser
449 //
450 const Double_t delayL0L1 = 0.071; //this is hack for 1/2 weeks
451 THnSparse *hisN=0;
452 TGraphErrors *grLaser[6]={0,0,0,0,0,0};
453 hisN = timeDrift->GetHistVdriftLaserA(0);
454 if (timeDrift->GetHistVdriftLaserA(0)){
455 grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
456 grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
457 vdriftArray->AddLast(grLaser[0]);
458 }
459 if (timeDrift->GetHistVdriftLaserA(1)){
460 grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
461 grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
462 vdriftArray->AddLast(grLaser[1]);
463 }
464 if (timeDrift->GetHistVdriftLaserA(2)){
465 grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
466 grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
467 vdriftArray->AddLast(grLaser[2]);
468 }
469 if (timeDrift->GetHistVdriftLaserC(0)){
470 grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
471 grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
472 vdriftArray->AddLast(grLaser[3]);
473 }
474 if (timeDrift->GetHistVdriftLaserC(1)){
475 grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
476 grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
477 vdriftArray->AddLast(grLaser[4]);
478 }
479 if (timeDrift->GetHistVdriftLaserC(2)){
480 grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
481 grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
482 vdriftArray->AddLast(grLaser[5]);
483 }
484 for (Int_t i=0; i<6;i++){
485 if (grLaser[i]) {
486 SetDefaultGraphDrift(grLaser[i], 1,(i+20));
487 grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
488 }
489 }
490}
491
492
493TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
494 //
495 // Make graph with mean values and rms
496 //
497 hisN->GetAxis(itime)->SetRange(0,100000000);
498 hisN->GetAxis(ival)->SetRange(0,100000000);
499 TH1 * hisT = hisN->Projection(itime);
500 TH1 * hisV = hisN->Projection(ival);
501 //
502 Int_t firstBinA = hisT->FindFirstBinAbove(2);
503 Int_t lastBinA = hisT->FindLastBinAbove(2);
504 Int_t firstBinV = hisV->FindFirstBinAbove(0);
505 Int_t lastBinV = hisV->FindLastBinAbove(0);
506 hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
507 hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
508 Int_t entries=0;
509 for (Int_t ibin=firstBinA; ibin<lastBinA; ibin++){
510 Double_t cont = hisT->GetBinContent(ibin);
511 if (cont<minEntries) continue;
512 entries++;
513 }
514 TVectorD vecTime(entries);
515 TVectorD vecMean0(entries);
516 TVectorD vecRMS0(entries);
517 TVectorD vecMean1(entries);
518 TVectorD vecRMS1(entries);
519 entries=0;
520 {for (Int_t ibin=firstBinA; ibin<lastBinA; ibin++){
521 Double_t cont = hisT->GetBinContent(ibin);
522 if (cont<minEntries) continue;
523 hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
524 Double_t time = hisT->GetBinCenter(ibin);
525 TH1 * his = hisN->Projection(ival);
526 Double_t nentries0= his->GetBinContent(his->FindBin(0));
527 if (cont-nentries0<minEntries) continue;
528 //
529 his->SetBinContent(his->FindBin(0),0);
530 vecTime[entries]=time;
531 vecMean0[entries]=his->GetMean()+offset;
532 vecMean1[entries]=his->GetMeanError();
533 vecRMS0[entries] =his->GetRMS();
534 vecRMS1[entries] =his->GetRMSError();
535 delete his;
536 entries++;
537 }}
538 delete hisT;
539 delete hisV;
540 TGraphErrors * graph = new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(), 0, vecMean1.GetMatrixArray());
541 return graph;
542}
543
544
545
546
547
548
549
550
551void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
552 //
40120adc 553 // Set default style for QA views
da268e11 554 //
555 graph->GetXaxis()->SetTimeDisplay(kTRUE);
556 graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
557 graph->SetMaximum( 0.025);
558 graph->SetMinimum(-0.025);
559 graph->GetXaxis()->SetTitle("Time");
560 graph->GetYaxis()->SetTitle("v_{dcorr}");
561 //
562 graph->GetYaxis()->SetLabelSize(0.03);
563 graph->GetXaxis()->SetLabelSize(0.03);
564 //
565 graph->GetXaxis()->SetNdivisions(10,5,0);
566 graph->GetYaxis()->SetNdivisions(10,5,0);
567 //
568 graph->GetXaxis()->SetLabelOffset(0.02);
569 graph->GetYaxis()->SetLabelOffset(0.005);
570 //
571 graph->GetXaxis()->SetTitleOffset(1.3);
572 graph->GetYaxis()->SetTitleOffset(1.2);
573 //
574 graph->SetMarkerColor(color);
575 graph->SetLineColor(color);
576 graph->SetMarkerStyle(style);
577}
578
579void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
40120adc 580 //
581 // Set default pad style for QA
582 //
da268e11 583 pad->SetTicks(1,1);
584 pad->SetMargin(mx0,mx1,my0,my1);
585}
586
587
40120adc 588void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray *picArray){
da268e11 589 //
40120adc 590 // 0. make a default QA plots
591 // 1. Store them in the array
da268e11 592 //
593 //
da268e11 594 Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
595 //
596 TGraphErrors* laserA =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
597 TGraphErrors* laserC =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
598 TGraphErrors* cosmic =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
599 TGraphErrors* cross =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
600 TGraphErrors* itstpcP =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
601 TGraphErrors* itstpcM =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
602 TGraphErrors* itstpcB =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
603 //
604 if (laserA) SetDefaultGraphDrift(laserA,2,25);
605 if (laserC) SetDefaultGraphDrift(laserC,4,26);
606 if (cosmic) SetDefaultGraphDrift(cosmic,3,27);
607 if (cross) SetDefaultGraphDrift(cross,4,28);
608 if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
609 if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
610 if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
611 //
612 //
613 TPad *pad=0;
614 //
615 // Laser-Laser
616 //
617 if (laserA&&laserC){
618 pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
619 laserA->Draw("alp");
620 SetPadStyle(pad,mx0,mx1,my0,my1);
621 laserA->Draw("apl");
622 laserC->Draw("p");
623 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
624 legend->AddEntry(laserA,"Laser A side");
625 legend->AddEntry(laserC,"Laser C side");
626 legend->Draw();
627 picArray->AddLast(pad);
628 }
629
630 if (itstpcP&&itstpcM){
631 pad = new TCanvas("ITSTPC","ITSTPC");
632 itstpcP->Draw("alp");
633 SetPadStyle(pad,mx0,mx1,my0,my1);
634 itstpcP->Draw("alp");
635 gPad->Clear();
636 itstpcM->Draw("apl");
637 itstpcP->Draw("p");
638 itstpcB->Draw("p");
639 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
640 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
641 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
642 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
643 legend->Draw();
644 picArray->AddLast(pad);
645 }
646
647 if (itstpcB&&laserA){
648 pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
649 SetPadStyle(pad,mx0,mx1,my0,my1);
650 laserA->Draw("alp");
651 itstpcP->Draw("p");
652 itstpcM->Draw("p");
653 itstpcB->Draw("p");
654 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
655 legend->AddEntry(laserA,"TPC laser");
656 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
657 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
658 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
659 legend->Draw();
660 picArray->AddLast(pad);
661 }
662
663 if (itstpcP&&cross){
664 pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
665 SetPadStyle(pad,mx0,mx1,my0,my1);
666 itstpcP->Draw("alp");
667 pad->Clear();
668 cross->Draw("ap");
669 itstpcP->Draw("p");
670 //
671 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
672
673 legend->AddEntry(cross,"TPC cross tracks");
674 legend->AddEntry(itstpcB,"ITS-TPC smooth");
675 legend->Draw();
676 picArray->AddLast(pad);
677 }
678 if (itstpcP&&cosmic){
679 pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
680 SetPadStyle(pad,mx0,mx1,my0,my1);
681 itstpcP->Draw("alp");
682 pad->Clear();
683 cosmic->Draw("ap");
684 itstpcP->Draw("p");
685 //
686 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
687
688 legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
689 legend->AddEntry(itstpcB,"ITS-TPC smooth");
690 legend->Draw();
691 picArray->AddLast(pad);
692 }
693}
694
695
696
697
40120adc 698void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber, TString pocdbStorage){
da268e11 699 //
700 // Update OCDB gain
701 //
702 ReadGainGlobal(fileName);
703 AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
704 MakeQAPlot(1.43);
f12abd89 705 if (pocdbStorage.Length()==0) pocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
706 UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage.Data());
da268e11 707}
708
709
710
711
40120adc 712void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
da268e11 713 //
714 // read calibration entries from file
715 //
716 TFile fcalib(fileName);
717 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
718 if (array){
719 fGainMIP = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
720 fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
721 }else{
722 fGainMIP = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
723 fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
724 }
725 TH1 * hisT=0;
726 Int_t firstBinA =0, lastBinA=0;
727
728 if (fGainCosmic){
729 hisT= fGainCosmic->GetHistGainTime()->Projection(1);
730 firstBinA = hisT->FindFirstBinAbove(2);
731 lastBinA = hisT->FindLastBinAbove(2);
732 fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
733 delete hisT;
734 }
735
736 if (fGainMIP){
737 hisT= fGainMIP->GetHistGainTime()->Projection(1);
738 firstBinA = hisT->FindFirstBinAbove(2);
739 lastBinA = hisT->FindLastBinAbove(2);
740 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
741 delete hisT;
742 }
743
744}
745
746
747
748Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit, Float_t FPtoMIPratio){
749 //
40120adc 750 // Analyze gain - produce the calibration graphs
da268e11 751 //
752 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
753 // 1.) try to create MIP spline
754 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
755 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
756 //
757 fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
758 if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
759 if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
760 if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
761 fGainArray->AddAt(fFitMIP,0);
762
763
764 // 2.) try to create Cosmic spline
765 fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
766 fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
767 //
768 fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
769 if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
770 //
771 if (fGraphCosmic) {
772 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
773 fGraphCosmic->GetY()[i] /= FPtoMIPratio;
774 fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
775 }
776 }
777 //
778 if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
779 if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
780 fGainArray->AddAt(fFitCosmic,1);
781 // with naming convention and backward compatibility
782 fGainArray->AddAt(fGraphMIP,2);
783 fGainArray->AddAt(fGraphCosmic,3);
784 cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
785 return kTRUE;
786
787}
788
789
790
791void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
792 //
793 // Update OCDB entry
794 //
795 AliCDBMetaData *metaData= new AliCDBMetaData();
796 metaData->SetObjectClassName("TObjArray");
797 metaData->SetResponsible("Alexander Kalweit");
798 metaData->SetBeamPeriod(1);
799 metaData->SetAliRootVersion("05-24-00"); //root version
800 metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
801 AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
802 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
803 gStorage->Put(fGainArray, id1, metaData);
804}
805
806void AliTPCPreprocessorOffline::MakeQAPlot(Float_t FPtoMIPratio) {
807 //
808 // Make QA plot to visualize results
809 //
810 //
811 //
812 if (fGraphCosmic) {
813 TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
814 canvasCosmic->cd();
815 TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
816 gainHistoCosmic->SetDirectory(0);
817 gainHistoCosmic->SetName("GainHistoCosmic");
818 gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
819 gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
820 gainHistoCosmic->Draw("colz");
821 fGraphCosmic->SetMarkerStyle(25);
822 fGraphCosmic->Draw("lp");
823 fGraphCosmic->SetMarkerStyle(25);
824 TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
825 if (grfFitCosmic) {
826 for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
827 grfFitCosmic->GetY()[i] *= FPtoMIPratio;
828 }
829 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
830 fGraphCosmic->GetY()[i] *= FPtoMIPratio;
831 }
832 }
833 fGraphCosmic->Draw("lp");
834 grfFitCosmic->SetLineColor(2);
835 grfFitCosmic->Draw("lu");
836 fGainArray->AddLast(gainHistoCosmic);
837 fGainArray->AddLast(canvasCosmic->Clone());
838 delete canvasCosmic;
839 }
840 if (fFitMIP) {
841 TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
842 canvasMIP->cd();
843 TH2D * gainHistoMIP = fGainMIP->GetHistGainTime()->Projection(0,1);
844 gainHistoMIP->SetName("GainHistoCosmic");
845 gainHistoMIP->SetDirectory(0);
846 gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
847 gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
848 gainHistoMIP->Draw("colz");
849 fGraphMIP->SetMarkerStyle(25);
850 fGraphMIP->Draw("lp");
851 TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
852 grfFitMIP->SetLineColor(2);
853 grfFitMIP->Draw("lu");
854 fGainArray->AddLast(gainHistoMIP);
855 fGainArray->AddLast(canvasMIP->Clone());
856 delete canvasMIP;
857 }
858}
859
860
861