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