]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCPreprocessorOffline.cxx
fix from previous commit
[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;
b322e06a 28 TString ocdbPath="local:////"
29 ocdbPath+=gSystem->GetFromPipe("pwd");
30
31 proces.CalibTimeGain("CalibObjects.root",run0,run1,ocdbPath);
32 proces.CalibTimeVdrift("CalibObjects.root",run0,run1,ocdbPath);
f12abd89 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
da268e11 37
68d461b2 38 e.g.
6d64657a 39 gSystem->Load("libANALYSIS");
40 gSystem->Load("libTPCcalib");
41 AliTPCPreprocessorOffline proces;
42 proces.CalibTimeGain("TPCMultObjects.root",114000,140040,0);
68d461b2 43 TFile oo("OCDB/TPC/Calib/TimeGain/Run114000_121040_v0_s0.root")
44 TObjArray * arr = AliCDBEntry->GetObject()
45 arr->At(4)->Draw("alp")
46
da268e11 47*/
f12abd89 48#include "Riostream.h"
49#include <fstream>
da268e11 50#include "TMap.h"
51#include "TGraphErrors.h"
52#include "AliExternalTrackParam.h"
53#include "TFile.h"
08669268 54#include "TDirectory.h"
da268e11 55#include "TGraph.h"
56#include "TMultiGraph.h"
57#include "TCanvas.h"
58#include "THnSparse.h"
59#include "TLegend.h"
60#include "TPad.h"
61#include "TH2D.h"
391ffdb2 62#include "TH3D.h"
da268e11 63#include "AliTPCROC.h"
64#include "AliTPCCalROC.h"
da268e11 65#include "AliESDfriend.h"
da268e11 66#include "AliTPCcalibTime.h"
67#include "AliSplineFit.h"
68#include "AliCDBMetaData.h"
69#include "AliCDBId.h"
70#include "AliCDBManager.h"
71#include "AliCDBStorage.h"
72#include "AliTPCcalibBase.h"
73#include "AliTPCcalibDB.h"
74#include "AliTPCcalibDButil.h"
75#include "AliRelAlignerKalman.h"
76#include "AliTPCParamSR.h"
77#include "AliTPCcalibTimeGain.h"
56b51ad9 78#include "AliTPCcalibGainMult.h"
3db6421b 79#include "AliTPCcalibAlign.h"
da268e11 80#include "AliSplineFit.h"
f3023796 81#include "AliTPCComposedCorrection.h"
82#include "AliTPCExBTwist.h"
83#include "AliTPCCalibGlobalMisalignment.h"
84#include "TStatToolkit.h"
85#include "TChain.h"
86#include "TCut.h"
87#include "AliTrackerBase.h"
3db6421b 88#include "AliTracker.h"
da268e11 89#include "AliTPCPreprocessorOffline.h"
3db6421b 90#include "AliTPCCorrectionFit.h"
da268e11 91
c82bb898 92using std::endl;
93using std::cout;
da268e11 94ClassImp(AliTPCPreprocessorOffline)
95
96AliTPCPreprocessorOffline::AliTPCPreprocessorOffline():
f12abd89 97 TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
40120adc 98 fMinEntries(500), // minimal number of entries for fit
a9f9c69b 99 fStartRun(0), // start Run - used to make fast selection in THnSparse
100 fEndRun(0), // end Run - used to make fast selection in THnSparse
101 fStartTime(0), // fStartTime - used to make fast selection in THnSparse
102 fEndTime(0), // fEndTime - used to make fast selection in THnSparse
103 fOCDBstorage(0), // OCDB storage
da268e11 104 fVdriftArray(new TObjArray),
105 fTimeDrift(0),
106 fGraphMIP(0), // graph time dependence of MIP
107 fGraphCosmic(0), // graph time dependence at Plateu
391ffdb2 108 fGraphAttachmentMIP(0),
da268e11 109 fFitMIP(0), // fit of dependence - MIP
110 fFitCosmic(0), // fit of dependence - Plateu
111 fGainArray(new TObjArray), // array to be stored in the OCDB
112 fGainMIP(0), // calibration component for MIP
68d461b2 113 fGainCosmic(0), // calibration component for cosmic
56b51ad9 114 fGainMult(0),
f3023796 115 fAlignTree(0), // alignment tree
88b09aaa 116 fSwitchOnValidation(kFALSE), // flag to switch on validation of OCDB parameters
117 fMinGain(2.0),
118 fMaxGain(3.0),
3f76dadb 119 fMaxVdriftCorr(0.03),
120 fNtracksVdrift(0),
121 fMinTracksVdrift(0),
122 fNeventsVdrift(0),
123 fMinEventsVdrift(0),
124 fCalibrationStatus(0)
da268e11 125{
126 //
f12abd89 127 // default constructor
128 //
da268e11 129}
130
131AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
132 //
133 // Destructor
134 //
135}
136
137
138
139
40120adc 140void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const timeDrift){
da268e11 141 //
142 // find the fist and last run
143 //
f12abd89 144 TObjArray *hisArray =timeDrift->GetHistoDrift();
da268e11 145 {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
146 THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
da268e11 147 if (!addHist) continue;
2a3fb9f6 148 if (addHist->GetEntries()<fMinEntries) continue;
da268e11 149 TH1D* histo =addHist->Projection(3);
150 TH1D* histoTime=addHist->Projection(0);
151 printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
152
a9f9c69b 153 if (fStartRun<=0){
154 fStartRun=histo->FindFirstBinAbove(0);
155 fEndRun =histo->FindLastBinAbove(0);
da268e11 156 }else{
a9f9c69b 157 fStartRun=TMath::Min(histo->FindFirstBinAbove(0),fStartRun);
158 fEndRun =TMath::Max(histo->FindLastBinAbove(0),fEndRun);
da268e11 159 }
a9f9c69b 160 if (fStartTime==0){
161 fStartTime=histoTime->FindFirstBinAbove(0);
162 fEndTime =histoTime->FindLastBinAbove(0);
da268e11 163 }else{
a9f9c69b 164 fStartTime=TMath::Min(histoTime->FindFirstBinAbove(0),fStartTime);
165 fEndTime =TMath::Max(histoTime->FindLastBinAbove(0),fEndTime);
da268e11 166 }
167 delete histo;
168 delete histoTime;
169 }}
a9f9c69b 170 if (fStartRun<0) fStartRun=0;
171 if (fEndRun<0) fEndRun=100000000;
172 printf("Run range :\t%d-%d\n", fStartRun, fEndRun);
173 printf("Time range :\t%d-%d\n", fStartTime, fEndTime);
da268e11 174
175}
176
177
178
a9f9c69b 179void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, AliCDBStorage* pocdbStorage){
da268e11 180 //
40120adc 181 // make calibration of the drift velocity
182 // Input parameters:
183 // file - the location of input file
184 // ustartRun, uendrun - run validity period
185 // pocdbStorage - path to hte OCDB storage
186 // - if empty - local storage 'pwd' uesed
a9f9c69b 187 if (pocdbStorage) fOCDBstorage=pocdbStorage;
188 else {
189 TString localStorage = "local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
190 fOCDBstorage=AliCDBManager::Instance()->GetStorage(localStorage.Data());
191 }
192
da268e11 193 //
194 // 1. Initialization and run range setting
195 TFile fcalib(file);
08669268 196 TObject* obj = dynamic_cast<TObject*>(fcalib.Get("TPCCalib"));
197 TObjArray* array = dynamic_cast<TObjArray*>(obj);
198 TDirectory* dir = dynamic_cast<TDirectory*>(obj);
199 if (dir) {
200 fTimeDrift = dynamic_cast<AliTPCcalibTime*>(dir->Get("calibTime"));
201 }
202 else if (array){
68d461b2 203 fTimeDrift = (AliTPCcalibTime *)array->FindObject("calibTime");
204 } else {
205 fTimeDrift = (AliTPCcalibTime*)fcalib.Get("calibTime");
206 }
207 if(!fTimeDrift) return;
208
3f76dadb 209 //extract statistics
210 fNtracksVdrift = TMath::Nint(fTimeDrift->GetResHistoTPCITS(0)->GetEntries());
48b76105 211 //if we have 0 ITS TPC matches it means we have no ITS tracks and we try to use TPC-TOF matching for calibration
212 if (fNtracksVdrift==0) fNtracksVdrift=TMath::Nint(fTimeDrift->GetResHistoTPCTOF(0)->GetEntries());
3f76dadb 213 fNeventsVdrift = TMath::Nint(fTimeDrift->GetTPCVertexHisto(0)->GetEntries());
214
a9f9c69b 215 fStartRun=ustartRun;
216 fEndRun=ustartRun;
da268e11 217 TObjArray *hisArray =fTimeDrift->GetHistoDrift();
218 GetRunRange(fTimeDrift);
219 for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
220 THnSparse* addHist=(THnSparse*)hisArray->At(i);
221 if (!addHist) continue;
a9f9c69b 222 if (fStartTime<fEndTime) addHist->GetAxis(0)->SetRange(fStartTime-1,fEndTime+1);
223 if (fStartRun<fEndRun) addHist->GetAxis(3)->SetRange(fStartRun-1,fEndRun+1);
da268e11 224 }
225 //
226 //
227 // 2. extraction of the information
228 //
229 fVdriftArray = new TObjArray();
230 AddAlignmentGraphs(fVdriftArray,fTimeDrift);
40120adc 231 AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
da268e11 232 AddLaserGraphs(fVdriftArray,fTimeDrift);
3f76dadb 233
da268e11 234 //
235 // 3. Append QA plots
236 //
237 MakeDefaultPlots(fVdriftArray,fVdriftArray);
68d461b2 238
da268e11 239 //
68d461b2 240 // 4. validate OCDB entries
da268e11 241 //
68d461b2 242 if(fSwitchOnValidation==kTRUE && ValidateTimeDrift()==kFALSE) {
243 Printf("TPC time drift OCDB parameters out of range!");
244 return;
245 }
f3023796 246 //
247 //4.b make alignment
248 //
249 MakeFitTime();
250 TFile * ftime= TFile::Open("fitITSVertex.root");
251 if (ftime){
252 TObject * alignmentTime=ftime->Get("FitCorrectionTime");
253 if (alignmentTime) fVdriftArray->AddLast(alignmentTime);
254 }
68d461b2 255 //
256 //
257 // 5. update of OCDB
da268e11 258 //
259 //
a9f9c69b 260 UpdateOCDBDrift(ustartRun,uendRun,fOCDBstorage);
da268e11 261}
262
a9f9c69b 263void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun, AliCDBStorage* storage ){
da268e11 264 //
265 // Update OCDB
266 //
267 AliCDBMetaData *metaData= new AliCDBMetaData();
268 metaData->SetObjectClassName("TObjArray");
269 metaData->SetResponsible("Marian Ivanov");
270 metaData->SetBeamPeriod(1);
271 metaData->SetAliRootVersion("05-25-01"); //root version
272 metaData->SetComment("Calibration of the time dependence of the drift velocity");
273 AliCDBId* id1=NULL;
274 id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
a9f9c69b 275 storage->Put(fVdriftArray, (*id1), metaData);
da268e11 276}
277
88b09aaa 278Bool_t AliTPCPreprocessorOffline::ValidateTimeGain()
68d461b2 279{
280 //
281 // Validate time gain corrections
282 //
283 Printf("ValidateTimeGain..." );
88b09aaa 284 Float_t minGain = fMinGain;
285 Float_t maxGain = fMaxGain;
68d461b2 286
287 TGraphErrors *gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
adeab288 288 if (!gr) {
289 gr = (TGraphErrors*)fGainArray->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
3f76dadb 290 if (!gr)
291 {
292 fCalibrationStatus |= kCalibFailedTimeGain;
293 return kFALSE;
294 }
adeab288 295 Printf("Assuming given run is a cosmic run. Using gain calibration from Fermi-plateau muons.");
296 }
3f76dadb 297 if(gr->GetN()<1)
298 {
299 fCalibrationStatus |= kCalibFailedTimeGain;
300 return kFALSE;
301 }
68d461b2 302
303 // check whether gain in the range
304 for(Int_t iPoint=0; iPoint<gr->GetN(); iPoint++)
305 {
306 if(gr->GetY()[iPoint] < minGain || gr->GetY()[iPoint] > maxGain)
3f76dadb 307 {
308 fCalibrationStatus |= kCalibFailedTimeGain;
68d461b2 309 return kFALSE;
3f76dadb 310 }
68d461b2 311 }
312
313return kTRUE;
314}
315
da268e11 316
88b09aaa 317Bool_t AliTPCPreprocessorOffline::ValidateTimeDrift()
68d461b2 318{
319 //
320 // Validate time drift velocity corrections
321 //
322 Printf("ValidateTimeDrift..." );
323
88b09aaa 324 Float_t maxVDriftCorr = fMaxVdriftCorr;
325
68d461b2 326 TGraphErrors* gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
3dcbe287 327 Printf("ALIGN_ITSB_TPC_DRIFTVD graph = %p",gr);
48b76105 328 if (!gr)
3f76dadb 329 {
48b76105 330 gr = (TGraphErrors*)fVdriftArray->FindObject("ALIGN_TOFB_TPC_DRIFTVD");
331 Printf("ALIGN_TOFB_TPC_DRIFTVD graph = %p",gr);
3f76dadb 332 }
333
334 if(!gr)
335 {
336 fCalibrationStatus|=kCalibFailedTimeDrift;
337 return kFALSE;
338 }
48b76105 339
340 // for now we validate even with low statistics
341 ////check if we have enough statistics
342 //if (fNtracksVdrift<fMinTracksVdrift)
343 //{
344 // fCalibrationStatus|=kCalibFailedTimeDrift;
345 // return kFALSE;
346 //}
347
3dcbe287 348 if(gr->GetN()<1) {
349 Printf("ALIGN_ITSB_TPC_DRIFTVD number of points = %d",gr->GetN());
3f76dadb 350 {
351 fCalibrationStatus|=kCalibFailedTimeDrift;
352 return kFALSE;
353 }
3dcbe287 354 }
68d461b2 355
356 // check whether drift velocity corrections in the range
357 for(Int_t iPoint = 0; iPoint<gr->GetN(); iPoint++)
358 {
a9f9c69b 359 //Printf("Y value from the graph: %f",TMath::Abs(gr->GetY()[iPoint]));
68d461b2 360 if(TMath::Abs(gr->GetY()[iPoint]) > maxVDriftCorr)
3f76dadb 361 {
362 fCalibrationStatus|=kCalibFailedTimeDrift;
68d461b2 363 return kFALSE;
3f76dadb 364 }
68d461b2 365 }
366
367return kTRUE;
368}
da268e11 369
40120adc 370void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
da268e11 371 //
372 // update the OCDB entry for the nominal time0
373 //
374 //
375 // AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
376 AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
377 TGraphErrors *grT = (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
378 Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
379 Double_t deltaT = deltaTcm/param->GetDriftV();
380 paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
381 paramNew->Update();
382
383 AliCDBMetaData *metaData= new AliCDBMetaData();
384 metaData->SetObjectClassName("TObjArray");
385 metaData->SetResponsible("Marian Ivanov");
386 metaData->SetBeamPeriod(1);
387 metaData->SetAliRootVersion("05-25-02"); //root version
388 metaData->SetComment("Updated calibration of nominal time 0");
389 AliCDBId* id1=NULL;
f12abd89 390 id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
a9f9c69b 391 fOCDBstorage->Put(param, (*id1), metaData);
da268e11 392
393}
394
395
396void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
397 //
40120adc 398 // Print the names of the entries in array
da268e11 399 //
400 Int_t entries = array->GetEntries();
401 for (Int_t i=0; i<entries; i++){
402 if (!array->At(i)) continue;
403 printf("%d\t %s\n", i, array->At(i)->GetName());
404 }
405}
406
407
408
409TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
410 // 2 filters:
411 // 1. filter graph - error cut errSigmaCut
412 // 2. filter graph - medianCutAbs around median
413 //
414 // errSigmaCut - cut on error
415 // medianCutAbs - cut on value around median
416 Double_t dummy=0; //
417 //
418 // 1. filter graph - error cut errSigmaCut
419 //
420 TGraphErrors *graphF;
421 graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
422 delete graph;
423 if (!graphF) return 0;
424 graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
425 delete graphF;
426 if (!graph) return 0;
427 //
428 // filter graph - kMedianCutAbs around median
429 //
430 graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
431 delete graph;
432 if (!graphF) return 0;
433 graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
434 delete graphF;
435 if (!graph) return 0;
436 return graph;
437}
438
439
440
441TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
442 //
443 // filter outlyer measurement
444 // Only points around median +- cut filtered
445 //
446 if (!graph) return 0;
447 Int_t kMinPoints=2;
448 Int_t npoints0 = graph->GetN();
449 Int_t npoints=0;
450 Float_t rmsY=0;
451 Double_t *outx=new Double_t[npoints0];
452 Double_t *outy=new Double_t[npoints0];
453 Double_t *errx=new Double_t[npoints0];
454 Double_t *erry=new Double_t[npoints0];
455 //
456 //
0bc13e06 457 if (npoints0<kMinPoints) {
458 delete []outx;
459 delete []outy;
460 delete []errx;
461 delete []erry;
462 return 0;
463 }
da268e11 464 for (Int_t iter=0; iter<3; iter++){
465 npoints=0;
466 for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
467 if (graph->GetY()[ipoint]==0) continue;
468 if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;
469 outx[npoints] = graph->GetX()[ipoint];
470 outy[npoints] = graph->GetY()[ipoint];
471 errx[npoints] = graph->GetErrorX(ipoint);
472 erry[npoints] = graph->GetErrorY(ipoint);
473 npoints++;
474 }
475 if (npoints<=1) break;
476 medianY =TMath::Median(npoints,outy);
477 rmsY =TMath::RMS(npoints,outy);
478 }
479 TGraphErrors *graphOut=0;
480 if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry);
831e26ce 481 delete []outx;
482 delete []outy;
483 delete []errx;
484 delete []erry;
da268e11 485 return graphOut;
486}
487
488
40120adc 489void AliTPCPreprocessorOffline::AddHistoGraphs( TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
da268e11 490 //
491 // Add graphs corresponding to the alignment
492 //
493 const Double_t kErrSigmaCut=5; // error sigma cut - for filtering
494 const Double_t kMedianCutAbs=0.03; // error sigma cut - for filtering
495 //
496 TObjArray * array=timeDrift->GetHistoDrift();
497 if (array){
498 THnSparse* hist=NULL;
499 // 2.a) cosmics with different triggers
500 for (Int_t i=0; i<array->GetEntriesFast();i++){
501 hist=(THnSparseF*)array->UncheckedAt(i);
502 if(!hist) continue;
503 if (hist->GetEntries()<minEntries) continue;
504 //hist->Print();
505 TString name=hist->GetName();
506 Int_t dim[4]={0,1,2,3};
507 THnSparse* newHist=hist->Projection(4,dim);
508 newHist->SetName(name);
509 TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
0bc13e06 510 if (!graph) {
511 printf("Graph =%s filtered out\n", name.Data());
512 continue;
513 }
da268e11 514 printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
515 Int_t pos=name.Index("_");
516 name=name(pos,name.Capacity()-pos);
517 TString graphName=graph->ClassName();
518 graphName+=name;
519 graphName.ToUpper();
520 //
521 graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
56b51ad9 522 //
831e26ce 523 if (graph){
56b51ad9 524 graph->SetMarkerStyle(i%8+20);
525 graph->SetMarkerColor(i%7);
526 graph->GetXaxis()->SetTitle("Time");
527 graph->GetYaxis()->SetTitle("v_{dcor}");
528 graph->SetName(graphName);
529 graph->SetTitle(graphName);
530 printf("Graph %d\t=\t%s\n", i, graphName.Data());
531 vdriftArray->Add(graph);
da268e11 532 }
da268e11 533 }
534 }
535}
536
537
538
539
40120adc 540void AliTPCPreprocessorOffline::AddAlignmentGraphs( TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
da268e11 541 //
542 // Add graphs corresponding to alignment to the object array
543 //
544 TObjArray *arrayITS=0;
545 TObjArray *arrayTOF=0;
546 TObjArray *arrayTRD=0;
547 TMatrixD *mstatITS=0;
548 TMatrixD *mstatTOF=0;
549 TMatrixD *mstatTRD=0;
550 //
551 arrayITS=timeDrift->GetAlignITSTPC();
552 arrayTRD=timeDrift->GetAlignTRDTPC();
553 arrayTOF=timeDrift->GetAlignTOFTPC();
554
938ce618 555 if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.7,50,fMaxVdriftCorr);
556 if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.7,1000,fMaxVdriftCorr);
557 if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.7,50,fMaxVdriftCorr);
da268e11 558 //
559 TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
560 TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
561 TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
562 TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
563 TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
564 TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
565
566 TObjArray * arrayTRDP= 0x0;
567 TObjArray * arrayTRDM= 0x0;
568 TObjArray * arrayTRDB= 0x0;
569 arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
570 arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
571 arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
572 //
573 //
574 Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
575 TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
576 arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
577 arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
578 TString grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
579 "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
580 "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
581 TString grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
582 "DELTAX", "DELTAY", "DELTAZ",
583 "DRIFTVD", "T0", "VDGY"};
584
585
586 TVectorD vX(entries);
587 TVectorD vY(entries);
588 TVectorD vEx(entries);
589 TVectorD vEy(entries);
590 TObjArray *arr=0;
591 for (Int_t iarray=0; iarray<12; iarray++){
592 arr = arrays[iarray];
593 if (arr==0) continue;
594 for (Int_t ipar=0; ipar<9; ipar++){
595 Int_t counter=0;
596 for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
597 AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
598 if (!kalman) continue;
599 vX[counter]=kalman->GetTimeStamp();
600 vY[counter]=(*(kalman->GetState()))[ipar];
601 if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
602 vEx[counter]=0;
603 vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
604 counter++;
605 }
606
607 TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
608 vY.GetMatrixArray(),
609 vEx.GetMatrixArray(),
610 vEy.GetMatrixArray());
611 TString grName=grnames[iarray];
612 grName+="_TPC_";
613 grName+=grpar[ipar];
614 graph->SetName(grName.Data());
615 vdriftArray->AddLast(graph);
616 }
617 }
618}
619
620
621
622
623void AliTPCPreprocessorOffline::AddLaserGraphs( TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
624 //
625 // add graphs for laser
626 //
627 const Double_t delayL0L1 = 0.071; //this is hack for 1/2 weeks
0bc13e06 628 //THnSparse *hisN=0;
da268e11 629 TGraphErrors *grLaser[6]={0,0,0,0,0,0};
0bc13e06 630 //hisN = timeDrift->GetHistVdriftLaserA(0);
da268e11 631 if (timeDrift->GetHistVdriftLaserA(0)){
632 grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
633 grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
634 vdriftArray->AddLast(grLaser[0]);
635 }
636 if (timeDrift->GetHistVdriftLaserA(1)){
637 grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
638 grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
639 vdriftArray->AddLast(grLaser[1]);
640 }
641 if (timeDrift->GetHistVdriftLaserA(2)){
642 grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
643 grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
644 vdriftArray->AddLast(grLaser[2]);
645 }
646 if (timeDrift->GetHistVdriftLaserC(0)){
647 grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
648 grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
649 vdriftArray->AddLast(grLaser[3]);
650 }
651 if (timeDrift->GetHistVdriftLaserC(1)){
652 grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
653 grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
654 vdriftArray->AddLast(grLaser[4]);
655 }
656 if (timeDrift->GetHistVdriftLaserC(2)){
657 grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
658 grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
659 vdriftArray->AddLast(grLaser[5]);
660 }
661 for (Int_t i=0; i<6;i++){
662 if (grLaser[i]) {
663 SetDefaultGraphDrift(grLaser[i], 1,(i+20));
664 grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
665 }
666 }
667}
668
669
670TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
671 //
672 // Make graph with mean values and rms
673 //
674 hisN->GetAxis(itime)->SetRange(0,100000000);
675 hisN->GetAxis(ival)->SetRange(0,100000000);
676 TH1 * hisT = hisN->Projection(itime);
677 TH1 * hisV = hisN->Projection(ival);
678 //
679 Int_t firstBinA = hisT->FindFirstBinAbove(2);
680 Int_t lastBinA = hisT->FindLastBinAbove(2);
681 Int_t firstBinV = hisV->FindFirstBinAbove(0);
682 Int_t lastBinV = hisV->FindLastBinAbove(0);
683 hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
684 hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
685 Int_t entries=0;
68d461b2 686 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
da268e11 687 Double_t cont = hisT->GetBinContent(ibin);
688 if (cont<minEntries) continue;
689 entries++;
690 }
691 TVectorD vecTime(entries);
692 TVectorD vecMean0(entries);
693 TVectorD vecRMS0(entries);
694 TVectorD vecMean1(entries);
695 TVectorD vecRMS1(entries);
696 entries=0;
68d461b2 697 for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
da268e11 698 Double_t cont = hisT->GetBinContent(ibin);
699 if (cont<minEntries) continue;
68d461b2 700 //hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
701 Int_t minBin = ibin-1;
702 Int_t maxBin = ibin+1;
703 if(minBin <= 0) minBin = 1;
704 if(maxBin >= hisN->GetAxis(itime)->GetNbins()) maxBin = hisN->GetAxis(itime)->GetNbins()-1;
705 hisN->GetAxis(itime)->SetRange(minBin,maxBin);
706
da268e11 707 Double_t time = hisT->GetBinCenter(ibin);
708 TH1 * his = hisN->Projection(ival);
709 Double_t nentries0= his->GetBinContent(his->FindBin(0));
710 if (cont-nentries0<minEntries) continue;
711 //
712 his->SetBinContent(his->FindBin(0),0);
713 vecTime[entries]=time;
714 vecMean0[entries]=his->GetMean()+offset;
715 vecMean1[entries]=his->GetMeanError();
716 vecRMS0[entries] =his->GetRMS();
717 vecRMS1[entries] =his->GetRMSError();
718 delete his;
719 entries++;
68d461b2 720 }
da268e11 721 delete hisT;
722 delete hisV;
723 TGraphErrors * graph = new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(), 0, vecMean1.GetMatrixArray());
724 return graph;
725}
726
727
728
729
730
731
732
733
734void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
735 //
40120adc 736 // Set default style for QA views
da268e11 737 //
738 graph->GetXaxis()->SetTimeDisplay(kTRUE);
739 graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
740 graph->SetMaximum( 0.025);
741 graph->SetMinimum(-0.025);
742 graph->GetXaxis()->SetTitle("Time");
743 graph->GetYaxis()->SetTitle("v_{dcorr}");
744 //
745 graph->GetYaxis()->SetLabelSize(0.03);
746 graph->GetXaxis()->SetLabelSize(0.03);
747 //
748 graph->GetXaxis()->SetNdivisions(10,5,0);
749 graph->GetYaxis()->SetNdivisions(10,5,0);
750 //
751 graph->GetXaxis()->SetLabelOffset(0.02);
752 graph->GetYaxis()->SetLabelOffset(0.005);
753 //
754 graph->GetXaxis()->SetTitleOffset(1.3);
755 graph->GetYaxis()->SetTitleOffset(1.2);
756 //
757 graph->SetMarkerColor(color);
758 graph->SetLineColor(color);
759 graph->SetMarkerStyle(style);
760}
761
762void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
40120adc 763 //
764 // Set default pad style for QA
765 //
da268e11 766 pad->SetTicks(1,1);
767 pad->SetMargin(mx0,mx1,my0,my1);
768}
769
770
68d461b2 771void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray * /*picArray*/){
da268e11 772 //
40120adc 773 // 0. make a default QA plots
774 // 1. Store them in the array
da268e11 775 //
776 //
da268e11 777 Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
778 //
779 TGraphErrors* laserA =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
780 TGraphErrors* laserC =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
781 TGraphErrors* cosmic =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
782 TGraphErrors* cross =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
783 TGraphErrors* itstpcP =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
784 TGraphErrors* itstpcM =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
785 TGraphErrors* itstpcB =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
786 //
787 if (laserA) SetDefaultGraphDrift(laserA,2,25);
788 if (laserC) SetDefaultGraphDrift(laserC,4,26);
789 if (cosmic) SetDefaultGraphDrift(cosmic,3,27);
790 if (cross) SetDefaultGraphDrift(cross,4,28);
791 if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
792 if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
793 if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
794 //
795 //
796 TPad *pad=0;
797 //
798 // Laser-Laser
799 //
800 if (laserA&&laserC){
801 pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
802 laserA->Draw("alp");
803 SetPadStyle(pad,mx0,mx1,my0,my1);
804 laserA->Draw("apl");
805 laserC->Draw("p");
806 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
807 legend->AddEntry(laserA,"Laser A side");
808 legend->AddEntry(laserC,"Laser C side");
809 legend->Draw();
68d461b2 810 //picArray->AddLast(pad);
da268e11 811 }
812
831e26ce 813 if (itstpcP&&itstpcM&&itstpcB){
da268e11 814 pad = new TCanvas("ITSTPC","ITSTPC");
815 itstpcP->Draw("alp");
816 SetPadStyle(pad,mx0,mx1,my0,my1);
817 itstpcP->Draw("alp");
818 gPad->Clear();
819 itstpcM->Draw("apl");
820 itstpcP->Draw("p");
821 itstpcB->Draw("p");
822 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
823 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
824 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
825 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
826 legend->Draw();
68d461b2 827 //picArray->AddLast(pad);
da268e11 828 }
829
831e26ce 830 if (itstpcB&&laserA&&itstpcP&&itstpcM){
da268e11 831 pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
832 SetPadStyle(pad,mx0,mx1,my0,my1);
833 laserA->Draw("alp");
834 itstpcP->Draw("p");
835 itstpcM->Draw("p");
836 itstpcB->Draw("p");
837 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
838 legend->AddEntry(laserA,"TPC laser");
839 legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
840 legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
841 legend->AddEntry(itstpcB,"ITS-TPC smooth ");
842 legend->Draw();
68d461b2 843 //picArray->AddLast(pad);
da268e11 844 }
845
846 if (itstpcP&&cross){
847 pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
848 SetPadStyle(pad,mx0,mx1,my0,my1);
849 itstpcP->Draw("alp");
850 pad->Clear();
851 cross->Draw("ap");
852 itstpcP->Draw("p");
853 //
854 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
855
856 legend->AddEntry(cross,"TPC cross tracks");
857 legend->AddEntry(itstpcB,"ITS-TPC smooth");
858 legend->Draw();
68d461b2 859 //picArray->AddLast(pad);
da268e11 860 }
861 if (itstpcP&&cosmic){
862 pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
863 SetPadStyle(pad,mx0,mx1,my0,my1);
864 itstpcP->Draw("alp");
865 pad->Clear();
866 cosmic->Draw("ap");
867 itstpcP->Draw("p");
868 //
869 TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
870
871 legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
872 legend->AddEntry(itstpcB,"ITS-TPC smooth");
873 legend->Draw();
68d461b2 874 //picArray->AddLast(pad);
da268e11 875 }
876}
877
878
879
880
a9f9c69b 881void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage* pocdbStorage){
da268e11 882 //
883 // Update OCDB gain
884 //
a9f9c69b 885 if (pocdbStorage==0) {
886 TString localStorage = "local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
887 pocdbStorage = AliCDBManager::Instance()->GetStorage(localStorage.Data());
888 }
68d461b2 889
890 //
891 // 1. Read gain values
892 //
da268e11 893 ReadGainGlobal(fileName);
68d461b2 894
895 //
896 // 2. Extract calibration values
897 //
da268e11 898 AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
391ffdb2 899 AnalyzeAttachment(startRunNumber,endRunNumber);
56b51ad9 900 AnalyzePadRegionGain();
6d64657a 901 AnalyzeGainMultiplicity();
7e3e1a9c 902 AnalyzeGainChamberByChamber();
68d461b2 903 //
904 // 3. Make control plots
905 //
906 MakeQAPlot(1.43);
da268e11 907
68d461b2 908 //
909 // 4. validate OCDB entries
910 //
911 if(fSwitchOnValidation==kTRUE && ValidateTimeGain()==kFALSE) {
912 Printf("TPC time gain OCDB parameters out of range!");
913 return;
914 }
da268e11 915
68d461b2 916 //
917 // 5. Update OCDB
918 //
a9f9c69b 919 UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage);
68d461b2 920}
da268e11 921
40120adc 922void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
da268e11 923 //
924 // read calibration entries from file
925 //
926 TFile fcalib(fileName);
08669268 927 TObject* obj = dynamic_cast<TObject*>(fcalib.Get("TPCCalib"));
928 TObjArray * array = dynamic_cast<TObjArray*>(obj);
929 TDirectory * dir = dynamic_cast<TDirectory*>(obj);
930 if (dir) {
931 fGainMIP = dynamic_cast<AliTPCcalibTimeGain *>(dir->Get("calibTimeGain"));
932 fGainCosmic = dynamic_cast<AliTPCcalibTimeGain *>(dir->Get("calibTimeGainCosmic"));
933 fGainMult = dynamic_cast<AliTPCcalibGainMult *>(dir->Get("calibGainMult"));
934 }
935 else if (array){
da268e11 936 fGainMIP = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
937 fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
56b51ad9 938 fGainMult = ( AliTPCcalibGainMult *)array->FindObject("calibGainMult");
da268e11 939 }else{
940 fGainMIP = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
941 fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
56b51ad9 942 fGainMult = ( AliTPCcalibGainMult *)fcalib.Get("calibGainMult");
da268e11 943 }
f8d8b429 944 if (!fGainMult){
a9f9c69b 945 TFile calibMultFile("TPCMultObjects.root");
946 fGainMult = ( AliTPCcalibGainMult *)calibMultFile.Get("calibGainMult");
f8d8b429 947 }
da268e11 948 TH1 * hisT=0;
949 Int_t firstBinA =0, lastBinA=0;
950
951 if (fGainCosmic){
952 hisT= fGainCosmic->GetHistGainTime()->Projection(1);
953 firstBinA = hisT->FindFirstBinAbove(2);
954 lastBinA = hisT->FindLastBinAbove(2);
955 fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
956 delete hisT;
957 }
958
959 if (fGainMIP){
960 hisT= fGainMIP->GetHistGainTime()->Projection(1);
961 firstBinA = hisT->FindFirstBinAbove(2);
962 lastBinA = hisT->FindLastBinAbove(2);
963 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
964 delete hisT;
965 }
966
967}
968
da268e11 969Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit, Float_t FPtoMIPratio){
970 //
40120adc 971 // Analyze gain - produce the calibration graphs
da268e11 972 //
68d461b2 973
da268e11 974 // 1.) try to create MIP spline
68d461b2 975 if (fGainMIP)
976 {
977 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
978 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
979 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
980 //
981 fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
982 if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
983 if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
984 if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
985 fGainArray->AddAt(fFitMIP,0);
986 }
da268e11 987
988 // 2.) try to create Cosmic spline
68d461b2 989 if (fGainCosmic)
990 {
6c8b3862 991 fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
992 fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
993 //
994 fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
995 if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
996 //
997 if (fGraphCosmic) {
998 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
999 fGraphCosmic->GetY()[i] /= FPtoMIPratio;
1000 fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
1001 }
da268e11 1002 }
6c8b3862 1003 //
1004 if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
1005 if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
1006 fGainArray->AddAt(fFitCosmic,1);
da268e11 1007 }
da268e11 1008 // with naming convention and backward compatibility
1009 fGainArray->AddAt(fGraphMIP,2);
1010 fGainArray->AddAt(fGraphCosmic,3);
1011 cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
1012 return kTRUE;
1013
1014}
1015
391ffdb2 1016Bool_t AliTPCPreprocessorOffline::AnalyzeAttachment(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesFit) {
1017 //
1018 // determine slope as a function of mean driftlength
1019 //
68d461b2 1020 if(!fGainMIP) return kFALSE;
1021
391ffdb2 1022 fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
1023 //
1024 fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
1025 fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
1026 //
68d461b2 1027 fGainMIP->GetHistGainTime()->GetAxis(3)->SetRangeUser(125,250);// only full tracking region (driftlength)
1028 fGainMIP->GetHistGainTime()->GetAxis(0)->SetRangeUser(1.5,3.5);// only full tracking region (driftlength)
1029 //
391ffdb2 1030 TH3D * hist = fGainMIP->GetHistGainTime()->Projection(1, 0, 3);
1031 //
1032 Double_t *xvec = new Double_t[hist->GetNbinsX()];
1033 Double_t *yvec = new Double_t[hist->GetNbinsX()];
1034 Double_t *xerr = new Double_t[hist->GetNbinsX()];
1035 Double_t *yerr = new Double_t[hist->GetNbinsX()];
1036 Int_t counter = 0;
1037 //
1038 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
1039 Int_t nsum=0;
1040 Int_t imin = i;
1041 Int_t imax = i;
68d461b2 1042 for (Int_t idelta=0; idelta<5; idelta++){
391ffdb2 1043 //
1044 imin = TMath::Max(i-idelta,1);
1045 imax = TMath::Min(i+idelta,hist->GetNbinsX());
1046 nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
1047 //if (nsum==0) break;
1048 if (nsum>minEntriesFit) break;
1049 }
1050 if (nsum<minEntriesFit) continue;
1051 //
68d461b2 1052 fGainMIP->GetHistGainTime()->GetAxis(1)->SetRangeUser(hist->GetXaxis()->GetBinCenter(imin-1),hist->GetXaxis()->GetBinCenter(imax+1)); // define time range
1053 TH2D * histZdep = fGainMIP->GetHistGainTime()->Projection(0,3);
1054 TObjArray arr;
1055 histZdep->FitSlicesY(0,0,-1,0,"QNR",&arr);
1056 TH1D * driftDep = (TH1D*)arr.At(1);
1057 delete histZdep;
1058 //TGraphErrors * driftDep = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,3,100,1,0.,1);
1059 /*if (driftDep->GetN() < 4) {
391ffdb2 1060 delete driftDep;
68d461b2 1061 continue;
1062 }*/
391ffdb2 1063 //
68d461b2 1064 //TObjArray arr;
391ffdb2 1065 //
68d461b2 1066 TF1 pol1("polynom1","pol1",125,240);
391ffdb2 1067 //driftDep->Fit(&pol1,"QNRROB=0.8");
1068 driftDep->Fit(&pol1,"QNR");
68d461b2 1069 xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin-1)+hist->GetXaxis()->GetBinCenter(imax+1));
391ffdb2 1070 yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
68d461b2 1071 xerr[counter] = hist->GetXaxis()->GetBinCenter(imax+1)-hist->GetXaxis()->GetBinCenter(imin-1);
391ffdb2 1072 yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
1073 counter++;
1074 //
68d461b2 1075 //delete driftDep;
391ffdb2 1076 }
1077 //
1078 fGraphAttachmentMIP = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
1079 if (fGraphAttachmentMIP) fGraphAttachmentMIP->SetName("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");// set proper names according to naming convention
1080 fGainArray->AddLast(fGraphAttachmentMIP);
1081 //
1082 delete [] xvec;
1083 delete [] yvec;
1084 delete [] xerr;
1085 delete [] yerr;
1086 delete hist;
1087 //
1088 if (counter < 1) return kFALSE;
1089 return kTRUE;
391ffdb2 1090
68d461b2 1091}
391ffdb2 1092
da268e11 1093
56b51ad9 1094Bool_t AliTPCPreprocessorOffline::AnalyzePadRegionGain(){
1095 //
1096 // Analyze gain for different pad regions - produce the calibration graphs 0,1,2
1097 //
1098 if (fGainMult)
1099 {
1100 TH2D * histQmax = (TH2D*) fGainMult->GetHistPadEqual()->Projection(0,2);
1101 TH2D * histQtot = (TH2D*) fGainMult->GetHistPadEqual()->Projection(1,2);
1102 //
1103 TObjArray arr;
1104 histQmax->FitSlicesY(0,0,-1,0,"QNR",&arr);
1105 Double_t xMax[3] = {0,1,2};
1106 Double_t yMax[3] = {((TH1D*)arr.At(1))->GetBinContent(1),
1107 ((TH1D*)arr.At(1))->GetBinContent(2),
1108 ((TH1D*)arr.At(1))->GetBinContent(3)};
1109 Double_t yMaxErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1110 ((TH1D*)arr.At(1))->GetBinError(2),
1111 ((TH1D*)arr.At(1))->GetBinError(3)};
1112 TGraphErrors * fitPadRegionQmax = new TGraphErrors(3, xMax, yMax, 0, yMaxErr);
1113 //
1114 histQtot->FitSlicesY(0,0,-1,0,"QNR",&arr);
1115 Double_t xTot[3] = {0,1,2};
1116 Double_t yTot[3] = {((TH1D*)arr.At(1))->GetBinContent(1),
1117 ((TH1D*)arr.At(1))->GetBinContent(2),
1118 ((TH1D*)arr.At(1))->GetBinContent(3)};
1119 Double_t yTotErr[3] = {((TH1D*)arr.At(1))->GetBinError(1),
1120 ((TH1D*)arr.At(1))->GetBinError(2),
1121 ((TH1D*)arr.At(1))->GetBinError(3)};
1122 TGraphErrors * fitPadRegionQtot = new TGraphErrors(3, xTot, yTot, 0, yTotErr);
1123 //
1124 fitPadRegionQtot->SetName("TGRAPHERRORS_MEANQTOT_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1125 fitPadRegionQmax->SetName("TGRAPHERRORS_MEANQMAX_PADREGIONGAIN_BEAM_ALL");// set proper names according to naming convention
1126 //
1127 fGainArray->AddLast(fitPadRegionQtot);
1128 fGainArray->AddLast(fitPadRegionQmax);
1129 return kTRUE;
1130 }
1131 return kFALSE;
1132
1133}
1134
1135
6d64657a 1136Bool_t AliTPCPreprocessorOffline::AnalyzeGainMultiplicity() {
1137 //
1138 // Analyze gain as a function of multiplicity and produce calibration graphs
1139 //
1140 if (!fGainMult) return kFALSE;
1141 fGainMult->GetHistGainMult()->GetAxis(3)->SetRangeUser(3,3);
1142 TH2D * histMultMax = fGainMult->GetHistGainMult()->Projection(0,4);
1143 TH2D * histMultTot = fGainMult->GetHistGainMult()->Projection(1,4);
1144 histMultMax->RebinX(4);
1145 histMultTot->RebinX(4);
1146 //
1147 TObjArray arrMax;
1148 TObjArray arrTot;
1149 histMultMax->FitSlicesY(0,0,-1,0,"QNR",&arrMax);
1150 histMultTot->FitSlicesY(0,0,-1,0,"QNR",&arrTot);
1151 //
1152 TH1D * meanMax = (TH1D*)arrMax.At(1);
1153 TH1D * meanTot = (TH1D*)arrTot.At(1);
1154 Float_t meanMult = histMultMax->GetMean();
f8d8b429 1155 if(meanMax->GetBinContent(meanMax->FindBin(meanMult))) {
1156 meanMax->Scale(1./meanMax->GetBinContent(meanMax->FindBin(meanMult)));
1157 }
1158 else {
1159 return kFALSE;
1160 }
1161 if(meanTot->GetBinContent(meanTot->FindBin(meanMult))) {
1162 meanTot->Scale(1./meanTot->GetBinContent(meanTot->FindBin(meanMult)));
1163 }
1164 else {
1165 return kFALSE;
1166 }
6d64657a 1167 Float_t xMultMax[50];
1168 Float_t yMultMax[50];
1169 Float_t yMultErrMax[50];
1170 Float_t xMultTot[50];
1171 Float_t yMultTot[50];
1172 Float_t yMultErrTot[50];
1173 //
1174 Int_t nCountMax = 0;
1175 for(Int_t iBin = 1; iBin < meanMax->GetXaxis()->GetNbins(); iBin++) {
1176 Float_t yValMax = meanMax->GetBinContent(iBin);
b5547fe4 1177 if (yValMax < 0.7) continue;
1178 if (yValMax > 1.3) continue;
6d64657a 1179 if (meanMax->GetBinError(iBin)/yValMax > 0.01) continue;
1180 xMultMax[nCountMax] = meanMax->GetXaxis()->GetBinCenter(iBin);
1181 yMultMax[nCountMax] = yValMax;
1182 yMultErrMax[nCountMax] = meanMax->GetBinError(iBin);
1183 nCountMax++;
1184 }
1185 //
1186 if (nCountMax < 10) return kFALSE;
1187 TGraphErrors * fitMultMax = new TGraphErrors(nCountMax, xMultMax, yMultMax, 0, yMultErrMax);
1188 fitMultMax->SetName("TGRAPHERRORS_MEANQMAX_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1189 //
1190 Int_t nCountTot = 0;
1191 for(Int_t iBin = 1; iBin < meanTot->GetXaxis()->GetNbins(); iBin++) {
1192 Float_t yValTot = meanTot->GetBinContent(iBin);
3dcbe287 1193 if (yValTot < 0.7) continue;
1194 if (yValTot > 1.3) continue;
6d64657a 1195 if (meanTot->GetBinError(iBin)/yValTot > 0.1) continue;
1196 xMultTot[nCountTot] = meanTot->GetXaxis()->GetBinCenter(iBin);
1197 yMultTot[nCountTot] = yValTot;
1198 yMultErrTot[nCountTot] = meanTot->GetBinError(iBin);
1199 nCountTot++;
1200 }
1201 //
1202 if (nCountTot < 10) return kFALSE;
1203 TGraphErrors * fitMultTot = new TGraphErrors(nCountTot, xMultTot, yMultTot, 0, yMultErrTot);
1204 fitMultTot->SetName("TGRAPHERRORS_MEANQTOT_MULTIPLICITYDEPENDENCE_BEAM_ALL");
1205 //
1206 fGainArray->AddLast(fitMultMax);
1207 fGainArray->AddLast(fitMultTot);
1208 //
1209 return kTRUE;
1210
1211}
1212
7e3e1a9c 1213Bool_t AliTPCPreprocessorOffline::AnalyzeGainChamberByChamber(){
1214 //
1215 // get chamber by chamber gain
1216 //
a9f9c69b 1217 if (!fGainMult) return kFALSE;
7e3e1a9c 1218 TGraphErrors *grShort = fGainMult->GetGainPerChamber(0);
1219 TGraphErrors *grMedium = fGainMult->GetGainPerChamber(1);
1220 TGraphErrors *grLong = fGainMult->GetGainPerChamber(2);
1221 if (grShort==0x0 || grMedium==0x0 || grLong==0x0) {
1222 delete grShort;
1223 delete grMedium;
1224 delete grLong;
1225 return kFALSE;
1226 }
1227
1228 fGainArray->AddLast(grShort);
1229 fGainArray->AddLast(grMedium);
1230 fGainArray->AddLast(grLong);
1231
1232 return kTRUE;
1233}
6d64657a 1234
a9f9c69b 1235void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage *storage){
da268e11 1236 //
1237 // Update OCDB entry
1238 //
1239 AliCDBMetaData *metaData= new AliCDBMetaData();
1240 metaData->SetObjectClassName("TObjArray");
1241 metaData->SetResponsible("Alexander Kalweit");
1242 metaData->SetBeamPeriod(1);
1243 metaData->SetAliRootVersion("05-24-00"); //root version
1244 metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
1245 AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
a9f9c69b 1246 storage->Put(fGainArray, id1, metaData);
da268e11 1247}
1248
1249void AliTPCPreprocessorOffline::MakeQAPlot(Float_t FPtoMIPratio) {
1250 //
1251 // Make QA plot to visualize results
1252 //
1253 //
1254 //
1255 if (fGraphCosmic) {
1256 TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
1257 canvasCosmic->cd();
1258 TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
1259 gainHistoCosmic->SetDirectory(0);
1260 gainHistoCosmic->SetName("GainHistoCosmic");
1261 gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
1262 gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1263 gainHistoCosmic->Draw("colz");
1264 fGraphCosmic->SetMarkerStyle(25);
1265 fGraphCosmic->Draw("lp");
1266 fGraphCosmic->SetMarkerStyle(25);
1267 TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
1268 if (grfFitCosmic) {
1269 for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
1270 grfFitCosmic->GetY()[i] *= FPtoMIPratio;
1271 }
1272 for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
1273 fGraphCosmic->GetY()[i] *= FPtoMIPratio;
1274 }
1275 }
0bc13e06 1276 fGraphCosmic->Draw("lp");
1277 if (grfFitCosmic) {
1278 grfFitCosmic->SetLineColor(2);
1279 grfFitCosmic->Draw("lu");
1280 }
56b51ad9 1281 fGainArray->AddLast(gainHistoCosmic);
1282 //fGainArray->AddLast(canvasCosmic->Clone());
1283 delete canvasCosmic;
da268e11 1284 }
1285 if (fFitMIP) {
1286 TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
1287 canvasMIP->cd();
1288 TH2D * gainHistoMIP = fGainMIP->GetHistGainTime()->Projection(0,1);
1289 gainHistoMIP->SetName("GainHistoCosmic");
1290 gainHistoMIP->SetDirectory(0);
1291 gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
1292 gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
1293 gainHistoMIP->Draw("colz");
1294 fGraphMIP->SetMarkerStyle(25);
1295 fGraphMIP->Draw("lp");
1296 TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
1297 grfFitMIP->SetLineColor(2);
1298 grfFitMIP->Draw("lu");
1299 fGainArray->AddLast(gainHistoMIP);
68d461b2 1300 //fGainArray->AddLast(canvasMIP->Clone());
da268e11 1301 delete canvasMIP;
1302 }
1303}
1304
f3023796 1305void AliTPCPreprocessorOffline::MakeFitTime(){
1306 //
3db6421b 1307 // make aligment fit - store results in the file
f3023796 1308 //
1309 const Int_t kMinEntries=1000;
1310 MakeChainTime();
1311 MakePrimitivesTime();
1312 if (!fAlignTree) return;
1313 if (fAlignTree->GetEntries()<kMinEntries) return;
1314 fAlignTree->SetAlias("ptype","type");
1315 fAlignTree->SetAlias("hasITS","(1+0)");
1316 fAlignTree->SetAlias("dITS","1-2*(refX<40)");
1317 fAlignTree->SetAlias("isITS","refX>10");
1318 fAlignTree->SetAlias("isVertex","refX<10");
1319 //
1320 Int_t npointsMax=30000000;
1321 TStatToolkit toolkit;
1322 Double_t chi2=0;
1323 Int_t npoints=0;
1324 TVectorD param;
1325 TMatrixD covar;
1326
1327 TString fstringFast="";
1328 fstringFast+="FExBTwistX++";
1329 fstringFast+="FExBTwistY++";
1330 fstringFast+="FAlignRot0D++";
1331 fstringFast+="FAlignTrans0D++";
1332 fstringFast+="FAlignTrans1D++";
1333 //
1334 fstringFast+="hasITS*FAlignTrans0++";
1335 fstringFast+="hasITS*FAlignTrans1++";
1336 fstringFast+="hasITS*FAlignRot0++";
1337 fstringFast+="hasITS*FAlignRot1++";
1338 fstringFast+="hasITS*FAlignRot2++";
1339 //
1340 fstringFast+="dITS*FAlignTrans0++";
1341 fstringFast+="dITS*FAlignTrans1++";
1342 fstringFast+="dITS*FAlignRot0++";
1343 fstringFast+="dITS*FAlignRot1++";
1344 fstringFast+="dITS*FAlignRot2++";
1345
39c42ea2 1346 TCut cutFit="entries>10&&abs(mean)>0.00001&&rms>0";
f3023796 1347 fAlignTree->SetAlias("err","rms");
1348
1349 TString *strDeltaITS = TStatToolkit::FitPlaneConstrain(fAlignTree,"mean:err", fstringFast.Data(),cutFit, chi2,npoints,param,covar,-1,0, npointsMax, 1);
1350 strDeltaITS->Tokenize("++")->Print();
1351 fAlignTree->SetAlias("fitYFast",strDeltaITS->Data());
1352 //
1353 TVectorD paramC= param;
1354 TMatrixD covarC= covar;
1355 TStatToolkit::Constrain1D(fstringFast,"Trans0D",paramC,covarC,0, 0.1);
1356 TStatToolkit::Constrain1D(fstringFast,"Trans1D",paramC,covarC,0, 0.1);
1357 TStatToolkit::Constrain1D(fstringFast,"TwistX",paramC,covarC,0, 0.1);
1358 TStatToolkit::Constrain1D(fstringFast,"TwistY",paramC,covarC,0, 0.1);
1359 TString strFitConst=TStatToolkit::MakeFitString(fstringFast, paramC,covar);
1360 fAlignTree->SetAlias("fitYFastC",strFitConst.Data());
1361 CreateAlignTime(fstringFast,paramC);
1362
1363
1364}
1365
1366
1367void AliTPCPreprocessorOffline::MakeChainTime(){
3db6421b 1368 //
1369 //
f3023796 1370 //
3dcbe287 1371 TFile f("CalibObjects.root");
3db6421b 1372
f3023796 1373 // const char *cdtype[7]={"ITS","TRD","Vertex","TOF","TPC","TPC0","TPC1"};
1374 //const char *cptype[5]={"dy","dz","dsnp","dtheta","d1pt"};
1375 const char * hname[5]={"dy","dz","dsnp","dtheta","d1pt"};
1376 Int_t run=0;
27c79aef 1377 AliTPCcalibTime *calibTime = 0;
08669268 1378 TObject* obj = dynamic_cast<TObject*>(f.Get("TPCCalib"));
1379 TObjArray * array = dynamic_cast<TObjArray*>(obj);
1380 TDirectory * dir = dynamic_cast<TDirectory*>(obj);
1381 if (dir) {
1382 calibTime = dynamic_cast<AliTPCcalibTime*>(dir->Get("calibTime"));
1383 }
1384 else if (array){
27c79aef 1385 calibTime = (AliTPCcalibTime *)array->FindObject("calibTime");
1386 } else {
1387 calibTime = (AliTPCcalibTime*)f.Get("calibTime");
1388 }
f3023796 1389 if (!calibTime) return;
3db6421b 1390 AliTPCCorrectionFit::CreateAlignMaps(AliTracker::GetBz(), run);
f3023796 1391 TTreeSRedirector *pcstream = new TTreeSRedirector("meanITSVertex.root");
0bc13e06 1392 //
f3023796 1393 Int_t ihis=0;
1394 THnSparse *his = calibTime->GetResHistoTPCITS(ihis);
1395 if (his){
1396 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1397 his->GetAxis(2)->SetRange(0,1000000);
1398 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1399 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
f3023796 1400 }
0bc13e06 1401 ihis=1;
1402 his = calibTime->GetResHistoTPCITS(ihis);
1403 if (his){
1404 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1405 his->GetAxis(2)->SetRange(0,1000000);
1406 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1407 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
0bc13e06 1408 }
f3023796 1409 ihis=2;
1410 his = calibTime->GetResHistoTPCITS(ihis);
1411 if (his){
1412 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1413 his->GetAxis(2)->SetRange(0,1000000);
1414 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1415 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("ITS%s",hname[ihis]),run,85.,ihis,3);
f3023796 1416 }
1417 ihis=0;
1418 his = calibTime->GetResHistoTPCvertex(ihis);
1419 if (his){
1420 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1421 his->GetAxis(2)->SetRange(0,1000000);
1422 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1423 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
f3023796 1424 }
1425 ihis=2;
1426 his = calibTime->GetResHistoTPCvertex(ihis);
1427 if (his){
1428 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1429 his->GetAxis(2)->SetRange(0,1000000);
1430 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1431 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
f3023796 1432
0bc13e06 1433 }
1434 ihis=1;
1435 his = calibTime->GetResHistoTPCvertex(ihis);
1436 if (his){
1437 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1438 his->GetAxis(2)->SetRange(0,1000000);
1439 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1440 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("Vertex%s",hname[ihis]),run,0.,ihis,3);
0bc13e06 1441
1442 }
1443 ihis=0;
1444 his = calibTime->GetResHistoTPCTOF(ihis);
1445 if (his){
1446 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1447 his->GetAxis(2)->SetRange(0,1000000);
1448 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1449 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("TOF%s",hname[ihis]),run,0.,ihis,3);
0bc13e06 1450
1451 }
1452 ihis=0;
1453 his = calibTime->GetResHistoTPCTRD(ihis);
1454 if (his){
1455 his->GetAxis(1)->SetRangeUser(-1.1,1.1);
1456 his->GetAxis(2)->SetRange(0,1000000);
1457 his->GetAxis(3)->SetRangeUser(-0.35,0.35);
39c42ea2 1458 AliTPCCorrection::MakeDistortionMap(his,pcstream, Form("TRD%s",hname[ihis]),run,0.,ihis,3);
0bc13e06 1459
f3023796 1460 }
1461 delete pcstream;
1462}
1463
1464
1465Double_t AliTPCPreprocessorOffline::EvalAt(Double_t phi, Double_t refX, Double_t theta, Int_t corr, Int_t ptype){
1466 //
1467 //
1468 //
1469 Double_t sector = 9*phi/TMath::Pi();
1470 if (sector<0) sector+=18;
1471 Double_t y85=AliTPCCorrection::GetCorrSector(sector,85,theta,1,corr);
1472 Double_t y245=AliTPCCorrection::GetCorrSector(sector,245,theta,1,corr);
1473 if (ptype==0) return y85+(y245-y85)*(refX-85.)/(245.-85.);
1474 if (ptype==2) return (y245-y85)/(245.-85.);
1475 return 0;
1476}
1477
1478
39c42ea2 1479Double_t AliTPCPreprocessorOffline::EvalAtPar(Double_t phi0, Double_t snp, Double_t refX, Double_t theta, Int_t corr, Int_t ptype, Int_t nsteps){
1480 //
1481 // Fit the distortion along the line with the parabolic model
1482 // Parameters:
1483 // phi0 - phi at the entrance of the TPC
1484 // snp - local inclination angle at the entrance of the TPC
1485 // refX - ref X where the distortion is evanluated
1486 // theta
1487 //
1488 static TLinearFitter fitter(3,"pol2");
1489 fitter.ClearPoints();
1490 if (nsteps<3) nsteps=3;
1491 Double_t deltaX=(245-85)/(nsteps);
1492 for (Int_t istep=0; istep<(nsteps+1); istep++){
1493 //
1494 Double_t localX =85.+deltaX*istep;
1495 Double_t localPhi=phi0+deltaX*snp*istep;
1496 Double_t sector = 9*localPhi/TMath::Pi();
1497 if (sector<0) sector+=18;
1498 Double_t y=AliTPCCorrection::GetCorrSector(sector,localX,theta,1,corr);
1499 Double_t dlocalX=AliTPCCorrection::GetCorrSector(sector,localX,theta,0,corr);
1500 Double_t x[1]={localX-dlocalX};
1501 fitter.AddPoint(x,y);
1502 }
1503 fitter.Eval();
1504 Double_t par[3];
1505 par[0]=fitter.GetParameter(0);
1506 par[1]=fitter.GetParameter(1);
1507 par[2]=fitter.GetParameter(2);
1508
1509 if (ptype==0) return par[0]+par[1]*refX+par[2]*refX*refX;
1510 if (ptype==2) return par[1]+2*par[2]*refX;
1511 if (ptype==4) return par[2];
1512 return 0;
1513}
1514
1515
1516
1517
1518
f3023796 1519
1520void AliTPCPreprocessorOffline::MakePrimitivesTime(){
1521 //
1522 // Create primitive transformation to fit
1523 //
1524 fAlignTree=new TChain("fit","fit");
1525 fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdy");
1526 fAlignTree->AddFile("meanITSVertex.root",10000000,"ITSdsnp");
1527 fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdy");
1528 fAlignTree->AddFile("meanITSVertex.root",10000000,"Vertexdsnp");
1529 //
1530 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1531 Double_t bzField=AliTrackerBase::GetBz();
1532 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
1533 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
1534 Double_t wtP = -10.0 * (bzField) * vdrift / ezField ;
1535 AliTPCExBTwist *fitExBTwistX= new AliTPCExBTwist;
1536 AliTPCExBTwist *fitExBTwistY= new AliTPCExBTwist;
1537 AliTPCCalibGlobalMisalignment *trans0 =new AliTPCCalibGlobalMisalignment;
1538 AliTPCCalibGlobalMisalignment *trans1 =new AliTPCCalibGlobalMisalignment;
1539 AliTPCCalibGlobalMisalignment *trans0D =new AliTPCCalibGlobalMisalignment;
1540 AliTPCCalibGlobalMisalignment *trans1D =new AliTPCCalibGlobalMisalignment;
1541 AliTPCCalibGlobalMisalignment *rot0 =new AliTPCCalibGlobalMisalignment;
1542 AliTPCCalibGlobalMisalignment *rot1 =new AliTPCCalibGlobalMisalignment;
1543 AliTPCCalibGlobalMisalignment *rot2 =new AliTPCCalibGlobalMisalignment;
1544 AliTPCCalibGlobalMisalignment *rot3 =new AliTPCCalibGlobalMisalignment;
1545 //
1546 //
1547 fitExBTwistX->SetXTwist(0.001);
1548 fitExBTwistX->SetOmegaTauT1T2(wtP,1,1);
1549 //
1550 fitExBTwistY->SetYTwist(0.001);
1551 fitExBTwistY->SetOmegaTauT1T2(wtP,1,1);
1552 //
1553 TGeoHMatrix *matrixRot = new TGeoHMatrix;
1554 TGeoHMatrix *matrixX = new TGeoHMatrix;
1555 TGeoHMatrix *matrixY = new TGeoHMatrix;
1556 matrixX->SetDx(0.1);
1557 matrixY->SetDy(0.1);
1558 Double_t rotAngles0[9]={0};
1559 Double_t rotAngles1[9]={0};
1560 Double_t rotAngles2[9]={0};
1561 //
1562 Double_t rotAngles3[9]={0};
1563
1564 rotAngles0[0]=1; rotAngles0[4]=1; rotAngles0[8]=1;
1565 rotAngles1[0]=1; rotAngles1[4]=1; rotAngles1[8]=1;
1566 rotAngles2[0]=1; rotAngles2[4]=1; rotAngles2[8]=1;
1567 rotAngles3[0]=1; rotAngles3[4]=1; rotAngles3[8]=1;
1568
1569 rotAngles0[1]=-0.001;rotAngles0[3]=0.001;
1570 rotAngles1[5]=-0.001;rotAngles1[7]=0.001;
1571 rotAngles2[2]=0.001;rotAngles2[6]=-0.001;
1572 rotAngles3[1]=0.001;rotAngles3[3]=-0.001;
1573 matrixRot->SetRotation(rotAngles0);
1574 rot0->SetAlignGlobal(matrixRot);
1575 matrixRot->SetRotation(rotAngles1);
1576 rot1->SetAlignGlobal(matrixRot);
1577 matrixRot->SetRotation(rotAngles2);
1578 rot2->SetAlignGlobal(matrixRot);
1579 matrixRot->SetRotation(rotAngles3);
1580 rot3->SetAlignGlobalDelta(matrixRot);
1581 //
1582 trans0->SetAlignGlobal(matrixX);
1583 trans1->SetAlignGlobal(matrixY);
1584 trans0D->SetAlignGlobalDelta(matrixX);
1585 trans1D->SetAlignGlobalDelta(matrixY);
1586 fitExBTwistX->Init();
1587 fitExBTwistY->Init();
1588 //
1589 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistX->Clone()),100);
1590 fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwistY->Clone()),101);
1591 //
1592 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot0->Clone()),102);
1593 fitExBTwistY->AddVisualCorrection((AliTPCExBTwist*)(rot1->Clone()),103);
1594 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot2->Clone()),104);
1595 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(rot3->Clone()),105);
1596
1597 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0->Clone()),106);
1598 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1->Clone()),107);
1599 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans0D->Clone()),108);
1600 fitExBTwistX->AddVisualCorrection((AliTPCExBTwist*)(trans1D->Clone()),109);
1601 //
1602 fAlignTree->SetAlias("FExBTwistX", "AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,100,ptype)+0");
1603 fAlignTree->SetAlias("FExBTwistY","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,101,ptype)+0");
1604 fAlignTree->SetAlias("FAlignRot0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,102,ptype)+0");
1605 fAlignTree->SetAlias("FAlignRot0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,105,ptype)+0");
1606 fAlignTree->SetAlias("FAlignRot1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,103,ptype)+0");
1607 fAlignTree->SetAlias("FAlignRot2","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,104,ptype)+0");
1608 fAlignTree->SetAlias("FAlignTrans0","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,106,ptype)+0");
1609 fAlignTree->SetAlias("FAlignTrans1","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,107,ptype)+0");
1610 fAlignTree->SetAlias("FAlignTrans0D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,108,ptype)+0");
1611 fAlignTree->SetAlias("FAlignTrans1D","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,109,ptype)+0");
1612 //
1613 // test fast function
1614 //
1615// fAlignTree->Draw("FExBTwistX:ExBTwistX","isITS&&ptype==0&&abs(snp)<0.05","");
1616// fAlignTree->Draw("FExBTwistY:ExBTwistY","isITS&&ptype==0&&abs(snp)<0.05","");
1617// fAlignTree->Draw("FAlignRot0:alignRot0","isITS&&ptype==0&&abs(snp)<0.05","");
1618// fAlignTree->Draw("FAlignRot1:alignRot1","isITS&&ptype==0&&abs(snp)<0.05","");
1619// fAlignTree->Draw("FAlignRot2:alignRot2","isITS&&ptype==0&&abs(snp)<0.05","");
1620// //
1621// fAlignTree->Draw("FAlignTrans0:alignTrans0","isITS&&ptype==0&&abs(snp)<0.05","");
1622// fAlignTree->Draw("FAlignTrans1:alignTrans1","isITS&&ptype==0&&abs(snp)<0.05","");
1623
1624}
1625
1626
1627void AliTPCPreprocessorOffline::CreateAlignTime(TString fstring, TVectorD paramC){
1628 //
1629 //
1630 //
1631 //
1632 TGeoHMatrix *matrixDelta = new TGeoHMatrix;
1633 TGeoHMatrix *matrixGlobal = new TGeoHMatrix;
1634 Double_t rAngles[9];
1635 Int_t index=0;
1636 //
1637 index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans0D");
1638 if (index>=0) matrixDelta->SetDx(paramC[index+1]*0.1);
1639 index=TStatToolkit::GetFitIndex(fstring,"FAlignTrans1D");
1640 if (index>=0) matrixDelta->SetDy(paramC[index+1]*0.1);
1641 rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1642 index=TStatToolkit::GetFitIndex(fstring,"FAlignRot0D");
1643 rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1644 rAngles[5]=0; rAngles[7] =0;
1645 rAngles[2]=0; rAngles[6] =0;
1646 matrixDelta->SetRotation(rAngles);
1647 //
1648 //
1649 //
1650 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans0");
1651 if (index>=0) matrixGlobal->SetDx(paramC[index+1]*0.1);
1652 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignTrans1");
1653 if (index>=0) matrixGlobal->SetDy(paramC[index+1]*0.1);
1654 rAngles[0]=1; rAngles[4]=1; rAngles[8]=1;
1655 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot0");
1656 rAngles[1]=-paramC[index+1]*0.001; rAngles[3]=paramC[index+1]*0.001;
1657 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot1");
1658 rAngles[5]=-paramC[index+1]*0.001; rAngles[7]=paramC[index+1]*0.001;
1659 index=TStatToolkit::GetFitIndex(fstring,"hasITS*FAlignRot2");
1660 rAngles[2]=paramC[index+1]*0.001; rAngles[6] =-paramC[index+1]*0.001;
1661 matrixGlobal->SetRotation(rAngles);
1662 //
0bc13e06 1663 AliTPCCalibGlobalMisalignment *fitAlignTime =0;
f3023796 1664 fitAlignTime =new AliTPCCalibGlobalMisalignment;
1665 fitAlignTime->SetName("FitAlignTime");
1666 fitAlignTime->SetTitle("FitAlignTime");
1667 fitAlignTime->SetAlignGlobalDelta(matrixDelta);
1668 fitAlignTime->SetAlignGlobal(matrixGlobal);
1669 //
1670 AliTPCExBTwist * fitExBTwist= new AliTPCExBTwist;
1671 Int_t indexX=TStatToolkit::GetFitIndex(fstring,"ExBTwistX");
1672 Int_t indexY=TStatToolkit::GetFitIndex(fstring,"ExBTwistY");
1673 fitExBTwist->SetXTwist(0.001*paramC[indexX+1]); // 1 mrad twist in x
1674 fitExBTwist->SetYTwist(0.001*paramC[indexY+1]); // 1 mrad twist in x
1675 fitExBTwist->SetName("FitExBTwistTime");
1676 fitExBTwist->SetTitle("FitExBTwistTime");
1677 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
1678 Double_t bzField=AliTrackerBase::GetBz();
1679 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
1680
1681 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
1682 Double_t wt = -10.0 * (bzField) * vdrift / ezField ;
1683 //
1684 fitExBTwist->SetOmegaTauT1T2(wt,1,1);
1685 fitExBTwist->Init();
1686
1687 AliTPCComposedCorrection *corrTime = new AliTPCComposedCorrection;
1688 TObjArray *arr = new TObjArray;
1689 corrTime->SetCorrections(arr);
1690
1691 corrTime->GetCorrections()->Add(fitExBTwist);
1692 corrTime->GetCorrections()->Add(fitAlignTime);
1693 corrTime->SetName("FitCorrectionTime");
1694 corrTime->SetTitle("FitCorrectionTime");
1695
1696 fitExBTwist->AddVisualCorrection((AliTPCExBTwist*)(fitExBTwist->Clone()),1001);
1697 fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(fitAlignTime->Clone()),1002);
1698 fitAlignTime->AddVisualCorrection((AliTPCExBTwist*)(corrTime->Clone()),1003);
1699
1700
1701 fAlignTree->SetAlias("ExBTwistTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1001,ptype)+0");
1702 fAlignTree->SetAlias("AlignTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1002,ptype)+0");
1703 fAlignTree->SetAlias("FitCorrectionTime","AliTPCPreprocessorOffline::EvalAt(phi,refX,theta,1003,ptype)+0");
1704
1705
1706 TFile *f = new TFile("fitITSVertex.root","update");
1707 corrTime->Write("FitCorrectionTime");
1708 f->Close();
1709}
1710
3f76dadb 1711//_____________________________________________________________________________
1712Int_t AliTPCPreprocessorOffline::GetStatus()
1713{
1714 //get the calibration status
1715 // 0 means OK
1716 // positive numbers invalidate for unknown reasons.
1717 // negative numbers invalidate with a known reason (e.g. low statistics).
1718 // the returned integer has one bit set for every component that failed.
1719
1720 Bool_t enoughStatistics = (fNtracksVdrift>fMinTracksVdrift && fNeventsVdrift>fMinEventsVdrift);
1721
1722 if (!enoughStatistics)
1723 {
1724 fCalibrationStatus=-TMath::Abs(fCalibrationStatus);
1725 }
1726
1727 return fCalibrationStatus;
1728}