]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCcalibTime.cxx
Coverity fix:
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibTime.cxx
CommitLineData
c74c5f6c 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/*
74235403 17Comments to be written here:
c74c5f6c 18
74235403 191. What do we calibrate.
20
21 Time dependence of gain and drift velocity in order to account for changes in: temperature, pressure, gas composition.
c74c5f6c 22
23 AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime",0, 1213.9e+06, 1213.96e+06, 0.04e+04, 0.04e+04);
24
c74c5f6c 25*/
26
c74c5f6c 27#include "Riostream.h"
2811495d 28#include "TDatabasePDG.h"
29#include "TGraphErrors.h"
c74c5f6c 30#include "TH1F.h"
c74c5f6c 31#include "THnSparse.h"
32#include "TList.h"
33#include "TMath.h"
2811495d 34#include "TTimeStamp.h"
35#include "TTree.h"
c74c5f6c 36#include "TVectorD.h"
2811495d 37//#include "TChain.h"
38//#include "TFile.h"
39
40#include "AliDCSSensor.h"
41#include "AliDCSSensorArray.h"
c74c5f6c 42#include "AliESDEvent.h"
c74c5f6c 43#include "AliESDInputHandler.h"
2811495d 44#include "AliESDVertex.h"
45#include "AliESDfriend.h"
c74c5f6c 46#include "AliLog.h"
817766d5 47#include "AliRelAlignerKalman.h"
2811495d 48#include "AliTPCCalROC.h"
49#include "AliTPCParam.h"
c74c5f6c 50#include "AliTPCTracklet.h"
da6c0bc9 51#include "AliTPCcalibDB.h"
52#include "AliTPCcalibLaser.h"
2811495d 53#include "AliTPCcalibTime.h"
54#include "AliTPCclusterMI.h"
55#include "AliTPCseed.h"
817766d5 56#include "AliTrackPointArray.h"
2811495d 57#include "AliTracker.h"
70036ef3 58#include "AliKFVertex.h"
660f0543 59#include <AliLog.h>
60
c74c5f6c 61ClassImp(AliTPCcalibTime)
62
63
64AliTPCcalibTime::AliTPCcalibTime()
660f0543 65 :AliTPCcalibBase(),
66 fMemoryMode(1), // 0 -do not fill THnSparse with residuals 1- fill only important QA THn 2 - Fill all THnsparse for calibration
da6c0bc9 67 fLaser(0), // pointer to laser calibration
68 fDz(0), // current delta z
d3ce44cb 69 fCutMaxD(3), // maximal distance in rfi ditection
70 fCutMaxDz(25), // maximal distance in rfi ditection
c74c5f6c 71 fCutTheta(0.03), // maximal distan theta
2be25cc9 72 fCutMinDir(-0.99), // direction vector products
b9ab8e40 73 fCutTracks(2500),
97d3fea0 74 fArrayLaserA(0), //laser fit parameters C
75 fArrayLaserC(0), //laser fit parameters A
dde68d36 76 fArrayDz(0), //NEW! Tmap of V drifts for different triggers
817766d5 77 fAlignITSTPC(0), //alignemnt array ITS TPC match
78 fAlignTRDTPC(0), //alignemnt array TRD TPC match
79 fAlignTOFTPC(0), //alignemnt array TOF TPC match
97d3fea0 80 fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
2be25cc9 81 fTimeBins(0),
82 fTimeStart(0),
83 fTimeEnd(0),
84 fPtBins(0),
85 fPtStart(0),
86 fPtEnd(0),
87 fVdriftBins(0),
88 fVdriftStart(0),
89 fVdriftEnd(0),
90 fRunBins(0),
91 fRunStart(0),
92 fRunEnd(0)
c74c5f6c 93{
6ba68367 94 //
95 // default constructor
96 //
c74c5f6c 97 AliInfo("Default Constructor");
2be25cc9 98 for (Int_t i=0;i<3;i++) {
99 fHistVdriftLaserA[i]=0;
100 fHistVdriftLaserC[i]=0;
101 }
d3ce44cb 102 for (Int_t i=0;i<10;i++) {
74235403 103 fCosmiMatchingHisto[i]=0;
104 }
391ffdb2 105 //
106 for (Int_t i=0;i<5;i++) {
b9908d0b 107 fResHistoTPCCE[i]=0;
391ffdb2 108 fResHistoTPCITS[i]=0;
109 fResHistoTPCTRD[i]=0;
0b736a46 110 fResHistoTPCTOF[i]=0;
391ffdb2 111 fResHistoTPCvertex[i]=0;
70036ef3 112 fTPCVertex[i]=0;
113 }
114 for (Int_t i=0;i<12;i++) {
115 fTPCVertex[i]=0;
391ffdb2 116 }
b1fea14f 117 for (Int_t i=0;i<5;i++) {
118 fTPCVertexCorrelation[i]=0;
119 }
660f0543 120 static Int_t counter=0;
121 if (1) {
122 TTimeStamp s;
123 Int_t time=s;
124 AliInfo(Form("Counter Constructor\t%d\t%d",counter,time));
125 counter++;
126 }
127
c74c5f6c 128}
129
660f0543 130AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeVdrift, Int_t memoryMode)
c74c5f6c 131 :AliTPCcalibBase(),
660f0543 132 fMemoryMode(memoryMode), // 0 -do not fill THnSparse with residuals 1- fill only important QA THn 2 - Fill all THnsparse for calibration
d3ce44cb 133 fLaser(0), // pointer to laser calibration
134 fDz(0), // current delta z
135 fCutMaxD(5*0.5356), // maximal distance in rfi ditection
dde68d36 136 fCutMaxDz(40), // maximal distance in rfi ditection
d3ce44cb 137 fCutTheta(5*0.004644),// maximal distan theta
138 fCutMinDir(-0.99), // direction vector products
b9ab8e40 139 fCutTracks(2500),
97d3fea0 140 fArrayLaserA(new TObjArray(1000)), //laser fit parameters C
141 fArrayLaserC(new TObjArray(1000)), //laser fit parameters A
dde68d36 142 fArrayDz(0), //Tmap of V drifts for different triggers
817766d5 143 fAlignITSTPC(0), //alignemnt array ITS TPC match
144 fAlignTRDTPC(0), //alignemnt array TRD TPC match
145 fAlignTOFTPC(0), //alignemnt array TOF TPC match
97d3fea0 146 fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
2be25cc9 147 fTimeBins(0),
148 fTimeStart(0),
149 fTimeEnd(0),
150 fPtBins(0),
151 fPtStart(0),
152 fPtEnd(0),
153 fVdriftBins(0),
154 fVdriftStart(0),
155 fVdriftEnd(0),
156 fRunBins(0),
157 fRunStart(0),
158 fRunEnd(0)
c74c5f6c 159{
6ba68367 160 //
161 // Non deafaul constructor - to be used in the Calibration setups
162 //
163
c74c5f6c 164 SetName(name);
165 SetTitle(title);
2be25cc9 166 for (Int_t i=0;i<3;i++) {
167 fHistVdriftLaserA[i]=0;
168 fHistVdriftLaserC[i]=0;
169 }
c74c5f6c 170
391ffdb2 171 for (Int_t i=0;i<5;i++) {
b9908d0b 172 fResHistoTPCCE[i]=0;
391ffdb2 173 fResHistoTPCITS[i]=0;
174 fResHistoTPCTRD[i]=0;
0b736a46 175 fResHistoTPCTOF[i]=0;
391ffdb2 176 fResHistoTPCvertex[i]=0;
177 }
178
179
c74c5f6c 180 AliInfo("Non Default Constructor");
74235403 181 fTimeBins =(EndTime-StartTime)/deltaIntegrationTimeVdrift;
2be25cc9 182 fTimeStart =StartTime; //(((TObjString*)(mapGRP->GetValue("fAliceStartTime")))->GetString()).Atoi();
183 fTimeEnd =EndTime; //(((TObjString*)(mapGRP->GetValue("fAliceStopTime")))->GetString()).Atoi();
dde68d36 184 fPtBins = 400;
2be25cc9 185 fPtStart = -0.04;
186 fPtEnd = 0.04;
dde68d36 187 fVdriftBins = 500;
188 fVdriftStart= -0.1;
189 fVdriftEnd = 0.1;
32100b1d 190 fRunBins = 1000001;
dde68d36 191 fRunStart = -1.5;
32100b1d 192 fRunEnd = 999999.5;
2be25cc9 193
194 Int_t binsVdriftLaser[4] = {fTimeBins , fPtBins , fVdriftBins*20, fRunBins };
195 Double_t xminVdriftLaser[4] = {fTimeStart, fPtStart, fVdriftStart , fRunStart};
196 Double_t xmaxVdriftLaser[4] = {fTimeEnd , fPtEnd , fVdriftEnd , fRunEnd };
dde68d36 197 TString axisTitle[4]={
198 "T",
199 "#delta_{P/T}",
200 "value",
201 "run"
202 };
203 TString histoName[3]={
204 "Loffset",
205 "Lcorr",
206 "Lgy"
207 };
2be25cc9 208
dde68d36 209
2be25cc9 210 for (Int_t i=0;i<3;i++) {
211 fHistVdriftLaserA[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
212 fHistVdriftLaserC[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
dde68d36 213 fHistVdriftLaserA[i]->SetName(histoName[i]);
214 fHistVdriftLaserC[i]->SetName(histoName[i]);
215 for (Int_t iaxis=0; iaxis<4;iaxis++){
216 fHistVdriftLaserA[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
217 fHistVdriftLaserC[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
218 }
2be25cc9 219 }
220 fBinsVdrift[0] = fTimeBins;
221 fBinsVdrift[1] = fPtBins;
222 fBinsVdrift[2] = fVdriftBins;
223 fBinsVdrift[3] = fRunBins;
224 fXminVdrift[0] = fTimeStart;
225 fXminVdrift[1] = fPtStart;
226 fXminVdrift[2] = fVdriftStart;
227 fXminVdrift[3] = fRunStart;
228 fXmaxVdrift[0] = fTimeEnd;
229 fXmaxVdrift[1] = fPtEnd;
230 fXmaxVdrift[2] = fVdriftEnd;
231 fXmaxVdrift[3] = fRunEnd;
232
dde68d36 233 fArrayDz=new TObjArray();
817766d5 234 fAlignITSTPC = new TObjArray; //alignemnt array ITS TPC match
235 fAlignTRDTPC = new TObjArray; //alignemnt array ITS TPC match
236 fAlignTOFTPC = new TObjArray; //alignemnt array ITS TPC match
3e55050f 237 fAlignITSTPC->SetOwner(kTRUE);
238 fAlignTRDTPC->SetOwner(kTRUE);
239 fAlignTOFTPC->SetOwner(kTRUE);
240
c74c5f6c 241
d3ce44cb 242 fCosmiMatchingHisto[0]=new TH1F("Cosmics matching","p0-all" ,100,-10*0.5356 ,10*0.5356 );
243 fCosmiMatchingHisto[1]=new TH1F("Cosmics matching","p1-all" ,100,-10*4.541 ,10*4.541 );
244 fCosmiMatchingHisto[2]=new TH1F("Cosmics matching","p2-all" ,100,-10*0.01134 ,10*0.01134 );
245 fCosmiMatchingHisto[3]=new TH1F("Cosmics matching","p3-all" ,100,-10*0.004644,10*0.004644);
246 fCosmiMatchingHisto[4]=new TH1F("Cosmics matching","p4-all" ,100,-10*0.03773 ,10*0.03773 );
247 fCosmiMatchingHisto[5]=new TH1F("Cosmics matching","p0-isPair",100,-10*0.5356 ,10*0.5356 );
248 fCosmiMatchingHisto[6]=new TH1F("Cosmics matching","p1-isPair",100,-10*4.541 ,10*4.541 );
249 fCosmiMatchingHisto[7]=new TH1F("Cosmics matching","p2-isPair",100,-10*0.01134 ,10*0.01134 );
250 fCosmiMatchingHisto[8]=new TH1F("Cosmics matching","p3-isPair",100,-10*0.004644,10*0.004644);
251 fCosmiMatchingHisto[9]=new TH1F("Cosmics matching","p4-isPair",100,-10*0.03773 ,10*0.03773 );
70036ef3 252 for (Int_t i=0;i<12;i++) {
253 fTPCVertex[i]=0;
254 }
b1fea14f 255 for (Int_t i=0;i<5;i++) {
256 fTPCVertexCorrelation[i]=0;
257 }
391ffdb2 258 BookDistortionMaps();
70036ef3 259
74235403 260}
c74c5f6c 261
262AliTPCcalibTime::~AliTPCcalibTime(){
263 //
6ba68367 264 // Virtual Destructor
c74c5f6c 265 //
660f0543 266 static Int_t counter=0;
267 if (1) {
268 TTimeStamp s;
269 Int_t time=s;
270 AliInfo(Form("Counter Destructor\t%s\t%d\t%d",GetName(),counter,time));
271 counter++;
272 }
d3ce44cb 273 for(Int_t i=0;i<3;i++){
274 if(fHistVdriftLaserA[i]){
275 delete fHistVdriftLaserA[i];
276 fHistVdriftLaserA[i]=NULL;
277 }
278 if(fHistVdriftLaserC[i]){
279 delete fHistVdriftLaserC[i];
280 fHistVdriftLaserC[i]=NULL;
281 }
282 }
dde68d36 283 if(fArrayDz){
284 fArrayDz->SetOwner();
285 fArrayDz->Delete();
286 delete fArrayDz;
287 fArrayDz=NULL;
2be25cc9 288 }
d3ce44cb 289 for(Int_t i=0;i<5;i++){
290 if(fCosmiMatchingHisto[i]){
291 delete fCosmiMatchingHisto[i];
292 fCosmiMatchingHisto[i]=NULL;
293 }
74235403 294 }
391ffdb2 295
296 for (Int_t i=0;i<5;i++) {
b9908d0b 297 delete fResHistoTPCCE[i];
391ffdb2 298 delete fResHistoTPCITS[i];
299 delete fResHistoTPCTRD[i];
0b736a46 300 delete fResHistoTPCTOF[i];
391ffdb2 301 delete fResHistoTPCvertex[i];
b9908d0b 302 fResHistoTPCCE[i]=0;
391ffdb2 303 fResHistoTPCITS[i]=0;
304 fResHistoTPCTRD[i]=0;
0b736a46 305 fResHistoTPCTOF[i]=0;
391ffdb2 306 fResHistoTPCvertex[i]=0;
307 }
308
472f0066 309 if (fTPCVertex[0]) {
70036ef3 310 for (Int_t i=0;i<12;i++) delete fTPCVertex[i];
311 }
472f0066 312 if (fTPCVertexCorrelation[0]) {
b1fea14f 313 for (Int_t i=0;i<5;i++) delete fTPCVertexCorrelation[i];
314 }
70036ef3 315
3e55050f 316 fAlignITSTPC->SetOwner(kTRUE);
317 fAlignTRDTPC->SetOwner(kTRUE);
318 fAlignTOFTPC->SetOwner(kTRUE);
319
817766d5 320 fAlignITSTPC->Delete();
321 fAlignTRDTPC->Delete();
322 fAlignTOFTPC->Delete();
3e55050f 323 delete fAlignITSTPC;
324 delete fAlignTRDTPC;
325 delete fAlignTOFTPC;
c74c5f6c 326}
2be25cc9 327
9170f0f0 328// Bool_t AliTPCcalibTime::IsLaser(const AliESDEvent *const /*event*/) const{
329// //
330// // Indicator is laser event not yet implemented - to be done using trigger info or event specie
331// //
332// return kTRUE; //More accurate creteria to be added
333// }
334// Bool_t AliTPCcalibTime::IsCosmics(const AliESDEvent *const /*event*/){
335// //
336// // Indicator is cosmic event not yet implemented - to be done using trigger info or event specie
337// //
338
339// return kTRUE; //More accurate creteria to be added
340// }
341// Bool_t AliTPCcalibTime::IsBeam(const AliESDEvent *const /*event*/) const{
342// //
343// // Indicator is physic event not yet implemented - to be done using trigger info or event specie
344// //
345
346// return kTRUE; //More accurate creteria to be added
347// }
348void AliTPCcalibTime::ResetCurrent(){
6ba68367 349 //
9170f0f0 350 //ResetCurrent
6ba68367 351 //
dde68d36 352 fDz=0; //Reset current dz
2be25cc9 353}
6ba68367 354
355
356
2be25cc9 357void AliTPCcalibTime::Process(AliESDEvent *event){
6ba68367 358 //
359 // main function to make calibration
360 //
2be25cc9 361 if(!event) return;
362 if (event->GetNumberOfTracks()<2) return;
363 ResetCurrent();
9170f0f0 364 //if(IsLaser (event))
365 ProcessLaser (event);
366 //if(IsCosmics(event))
367 ProcessCosmic(event);
368 //if(IsBeam (event))
369 ProcessBeam (event);
2be25cc9 370}
371
372void AliTPCcalibTime::ProcessLaser(AliESDEvent *event){
c74c5f6c 373 //
2be25cc9 374 // Fit drift velocity using laser
375 //
376 // 0. cuts
dde68d36 377 const Int_t kMinTracks = 40; // minimal number of laser tracks
378 const Int_t kMinTracksSide = 20; // minimal number of tracks per side
379 const Float_t kMaxDeltaZ = 30.; // maximal trigger delay
380 const Float_t kMaxDeltaV = 0.05; // maximal deltaV
2be25cc9 381 const Float_t kMaxRMS = 0.1; // maximal RMS of tracks
dde68d36 382 //
2be25cc9 383 /*
384 TCut cutRMS("sqrt(laserA.fElements[4])<0.1&&sqrt(laserC.fElements[4])<0.1");
385 TCut cutZ("abs(laserA.fElements[0]-laserC.fElements[0])<3");
386 TCut cutV("abs(laserA.fElements[1]-laserC.fElements[1])<0.01");
387 TCut cutY("abs(laserA.fElements[2]-laserC.fElements[2])<2");
388 TCut cutAll = cutRMS+cutZ+cutV+cutY;
389 */
390 if (event->GetNumberOfTracks()<kMinTracks) return;
c74c5f6c 391 //
2be25cc9 392 if(!fLaser) fLaser = new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
393 fLaser->Process(event);
394 if (fLaser->GetNtracks()<kMinTracks) return; // small amount of tracks cut
dde68d36 395 if (fLaser->fFitAside->GetNrows()==0 && fLaser->fFitCside->GetNrows()==0) return; // no fit neither a or C side
c74c5f6c 396 //
2be25cc9 397 // debug streamer - activate stream level
398 // Use it for tuning of the cuts
da6c0bc9 399 //
dde68d36 400 // cuts to be applied
401 //
402 Int_t isReject[2]={0,0};
403 //
404 // not enough tracks
405 if (TMath::Abs((*fLaser->fFitAside)[3]) < kMinTracksSide) isReject[0]|=1;
406 if (TMath::Abs((*fLaser->fFitCside)[3]) < kMinTracksSide) isReject[1]|=1;
407 // unreasonable z offset
408 if (TMath::Abs((*fLaser->fFitAside)[0])>kMaxDeltaZ) isReject[0]|=2;
409 if (TMath::Abs((*fLaser->fFitCside)[0])>kMaxDeltaZ) isReject[1]|=2;
410 // unreasonable drift velocity
411 if (TMath::Abs((*fLaser->fFitAside)[1]-1)>kMaxDeltaV) isReject[0]|=4;
412 if (TMath::Abs((*fLaser->fFitCside)[1]-1)>kMaxDeltaV) isReject[1]|=4;
413 // big chi2
414 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitAside)[4]))>kMaxRMS ) isReject[0]|=8;
415 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitCside)[4]))>kMaxRMS ) isReject[1]|=8;
416
417
418
419
2be25cc9 420 if (fStreamLevel>0){
421 printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
422
da6c0bc9 423 TTreeSRedirector *cstream = GetDebugStreamer();
424 if (cstream){
425 TTimeStamp tstamp(fTime);
2be25cc9 426 (*cstream)<<"laserInfo"<<
da6c0bc9 427 "run="<<fRun<< // run number
428 "event="<<fEvent<< // event number
429 "time="<<fTime<< // time stamp of event
430 "trigger="<<fTrigger<< // trigger
431 "mag="<<fMagF<< // magnetic field
da6c0bc9 432 //laser
dde68d36 433 "rejectA="<<isReject[0]<<
434 "rejectC="<<isReject[1]<<
2be25cc9 435 "laserA.="<<fLaser->fFitAside<<
436 "laserC.="<<fLaser->fFitCside<<
437 "laserAC.="<<fLaser->fFitACside<<
438 "trigger="<<event->GetFiredTriggerClasses()<<
da6c0bc9 439 "\n";
440 }
441 }
2be25cc9 442 //
2be25cc9 443 // fill histos
444 //
012c4694 445 TVectorD vdriftA(5), vdriftC(5),vdriftAC(6);
2be25cc9 446 vdriftA=*(fLaser->fFitAside);
447 vdriftC=*(fLaser->fFitCside);
448 vdriftAC=*(fLaser->fFitACside);
449 Int_t npointsA=0, npointsC=0;
450 Float_t chi2A=0, chi2C=0;
451 npointsA= TMath::Nint(vdriftA[3]);
452 chi2A= vdriftA[4];
453 npointsC= TMath::Nint(vdriftC[3]);
454 chi2C= vdriftC[4];
2be25cc9 455
97d3fea0 456 if (npointsA>kMinTracksSide || npointsC>kMinTracksSide){
457 TVectorD *fitA = new TVectorD(6);
458 TVectorD *fitC = new TVectorD(6);
459 for (Int_t ipar=0; ipar<5; ipar++){
460 (*fitA)[ipar]=vdriftA[ipar];
461 (*fitC)[ipar]=vdriftC[ipar];
462 }
463 (*fitA)[5]=fTime;
464 (*fitC)[5]=fTime;
465 fArrayLaserA->AddLast(fitA);
466 fArrayLaserC->AddLast(fitC);
467 }
468 //
469
2be25cc9 470 TTimeStamp tstamp(fTime);
471 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
472 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
dde68d36 473 Double_t driftA=0, driftC=0;
474 if (vdriftA[1]>1.-kMaxDeltaV) driftA = 1./vdriftA[1]-1.;
475 if (vdriftC[1]>1.-kMaxDeltaV) driftC = 1./vdriftC[1]-1.;
476 //
477 Double_t vecDriftLaserA[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftA,event->GetRunNumber()};
478 Double_t vecDriftLaserC[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftC,event->GetRunNumber()};
479 // Double_t vecDrift[4] ={fTime,(ptrelative0+ptrelative1)/2.0,1./((*(fLaser->fFitACside))[1])-1,event->GetRunNumber()};
480
481 for (Int_t icalib=0;icalib<3;icalib++){
482 if (icalib==0){ //z0 shift
483 vecDriftLaserA[2]=vdriftA[0]/250.;
484 vecDriftLaserC[2]=vdriftC[0]/250.;
2be25cc9 485 }
dde68d36 486 if (icalib==1){ //vdrel shift
487 vecDriftLaserA[2]=driftA;
488 vecDriftLaserC[2]=driftC;
2be25cc9 489 }
dde68d36 490 if (icalib==2){ //gy shift - full gy - full drift
491 vecDriftLaserA[2]=vdriftA[2]/250.;
492 vecDriftLaserC[2]=vdriftC[2]/250.;
2be25cc9 493 }
391ffdb2 494 //if (isReject[0]==0) fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
495 //if (isReject[1]==0) fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
496 fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
497 fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
2be25cc9 498 }
2be25cc9 499}
c74c5f6c 500
2811495d 501void AliTPCcalibTime::ProcessCosmic(const AliESDEvent *const event){
6ba68367 502 //
503 // process Cosmic event - track matching A side C side
504 //
c74c5f6c 505 if (!event) {
506 Printf("ERROR: ESD not available");
507 return;
508 }
509 if (event->GetTimeStamp() == 0 ) {
510 Printf("no time stamp!");
511 return;
512 }
74235403 513
2be25cc9 514 //fd
c74c5f6c 515 // Find cosmic pairs
516 //
517 // Track0 is choosen in upper TPC part
518 // Track1 is choosen in lower TPC part
519 //
cc65e4f5 520 const Int_t kMinClustersCross =30;
521 const Int_t kMinClusters =80;
c74c5f6c 522 Int_t ntracks=event->GetNumberOfTracks();
523 if (ntracks==0) return;
74235403 524 if (ntracks > fCutTracks) return;
525
a390f11f 526 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
3e0d9569 527 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
c74c5f6c 528
529 TObjArray tpcSeeds(ntracks);
530 Double_t vtxx[3]={0,0,0};
531 Double_t svtxx[3]={0.000001,0.000001,100.};
532 AliESDVertex vtx(vtxx,svtxx);
533 //
534 // track loop
535 //
cc65e4f5 536 TArrayI clusterSideA(ntracks);
537 TArrayI clusterSideC(ntracks);
c74c5f6c 538 for (Int_t i=0;i<ntracks;++i) {
cc65e4f5 539 clusterSideA[i]=0;
540 clusterSideC[i]=0;
74235403 541 AliESDtrack *track = event->GetTrack(i);
542
543 const AliExternalTrackParam * trackIn = track->GetInnerParam();
544 const AliExternalTrackParam * trackOut = track->GetOuterParam();
545 if (!trackIn) continue;
546 if (!trackOut) continue;
547
6ba68367 548 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
b1c27e4f 549 if (!friendTrack) continue;
3e55050f 550 if (friendTrack) ProcessSame(track,friendTrack,event);
6ba68367 551 if (friendTrack) ProcessAlignITS(track,friendTrack,event,esdFriend);
0dac7b3a 552 if (friendTrack) ProcessAlignTRD(track,friendTrack);
817766d5 553 if (friendTrack) ProcessAlignTOF(track,friendTrack);
74235403 554 TObject *calibObject;
555 AliTPCseed *seed = 0;
556 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
cc65e4f5 557 if (seed) {
558 tpcSeeds.AddAt(seed,i);
559 Int_t nA=0, nC=0;
560 for (Int_t irow=159;irow>0;irow--) {
561 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
562 if (!cl) continue;
563 if ((cl->GetDetector()%36)<18) nA++;
564 if ((cl->GetDetector()%36)>=18) nC++;
565 }
566 clusterSideA[i]=nA;
567 clusterSideC[i]=nC;
568 }
c74c5f6c 569 }
c74c5f6c 570 if (ntracks<2) return;
571 //
572 // Find pairs
573 //
d3ce44cb 574
c74c5f6c 575 for (Int_t i=0;i<ntracks;++i) {
74235403 576 AliESDtrack *track0 = event->GetTrack(i);
c74c5f6c 577 // track0 - choosen upper part
578 if (!track0) continue;
579 if (!track0->GetOuterParam()) continue;
580 if (track0->GetOuterParam()->GetAlpha()<0) continue;
581 Double_t d1[3];
582 track0->GetDirection(d1);
583 for (Int_t j=0;j<ntracks;++j) {
584 if (i==j) continue;
74235403 585 AliESDtrack *track1 = event->GetTrack(j);
586 //track 1 lower part
587 if (!track1) continue;
588 if (!track1->GetOuterParam()) continue;
cc65e4f5 589 if (track0->GetTPCNcls()+ track1->GetTPCNcls()< kMinClusters) continue;
590 Int_t nAC = TMath::Max( TMath::Min(clusterSideA[i], clusterSideC[j]),
591 TMath::Min(clusterSideC[i], clusterSideA[j]));
592 if (nAC<kMinClustersCross) continue;
593 Int_t nA0=clusterSideA[i];
594 Int_t nC0=clusterSideC[i];
595 Int_t nA1=clusterSideA[j];
596 Int_t nC1=clusterSideC[j];
3e55050f 597 // if (track1->GetOuterParam()->GetAlpha()>0) continue;
74235403 598 //
599 Double_t d2[3];
600 track1->GetDirection(d2);
601
602 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
603 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
604 if (! seed0) continue;
605 if (! seed1) continue;
606 Float_t dir = (d1[0]*d2[0] + d1[1]*d2[1] + d1[2]*d2[2]);
607 Float_t dist0 = track0->GetLinearD(0,0);
608 Float_t dist1 = track1->GetLinearD(0,0);
609 //
610 // conservative cuts - convergence to be guarantied
611 // applying before track propagation
3e55050f 612 if (TMath::Abs(TMath::Abs(dist0)-TMath::Abs(dist1))>fCutMaxD) continue; // distance to the 0,0
613 if (TMath::Abs(dir)<TMath::Abs(fCutMinDir)) continue; // direction vector product
74235403 614 Float_t bz = AliTracker::GetBz();
615 Float_t dvertex0[2]; //distance to 0,0
616 Float_t dvertex1[2]; //distance to 0,0
617 track0->GetDZ(0,0,0,bz,dvertex0);
618 track1->GetDZ(0,0,0,bz,dvertex1);
619 if (TMath::Abs(dvertex0[1])>250) continue;
620 if (TMath::Abs(dvertex1[1])>250) continue;
621 //
622 //
623 //
624 Float_t dmax = TMath::Max(TMath::Abs(dist0),TMath::Abs(dist1));
625 AliExternalTrackParam param0(*track0);
626 AliExternalTrackParam param1(*track1);
627 //
628 // Propagate using Magnetic field and correct fo material budget
629 //
9963b5e2 630 AliTracker::PropagateTrackTo(&param0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
631 AliTracker::PropagateTrackTo(&param1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
74235403 632 //
633 // Propagate rest to the 0,0 DCA - z should be ignored
634 //
635 //Bool_t b0 = ;
636 param0.PropagateToDCA(&vtx,bz,1000);
637 //Bool_t b1 =
638 param1.PropagateToDCA(&vtx,bz,1000);
74235403 639 param0.GetDZ(0,0,0,bz,dvertex0);
640 param1.GetDZ(0,0,0,bz,dvertex1);
dde68d36 641 Double_t xyz0[3];
642 Double_t xyz1[3];
74235403 643 param0.GetXYZ(xyz0);
644 param1.GetXYZ(xyz1);
645 Bool_t isPair = IsPair(&param0,&param1);
d3ce44cb 646 Bool_t isCross = IsCross(track0, track1);
dde68d36 647 Bool_t isSame = IsSame(track0, track1);
d3ce44cb 648
dde68d36 649 THnSparse* hist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
650 TString shortName=hist->ClassName();
651 shortName+="_MEAN_VDRIFT_COSMICS_";
652 delete hist;
653 hist=NULL;
654
3e55050f 655 if((isSame) || (isCross && isPair)){
656 if (track0->GetTPCNcls()+ track1->GetTPCNcls()> 80) {
74235403 657 fDz = param0.GetZ() - param1.GetZ();
cc65e4f5 658 Double_t sign=(nA0>nA1)? 1:-1;
659 fDz*=sign;
74235403 660 TTimeStamp tstamp(fTime);
661 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
662 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
dde68d36 663 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
664 THnSparse* curHist=NULL;
665 TString name="";
666
667 name=shortName;
668 name+=event->GetFiredTriggerClasses();
669 name.ToUpper();
670 curHist=(THnSparseF*)fArrayDz->FindObject(name);
74235403 671 if(!curHist){
dde68d36 672 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
673 fArrayDz->AddLast(curHist);
74235403 674 }
dde68d36 675// curHist=(THnSparseF*)(fMapDz->GetValue(event->GetFiredTriggerClasses()));
676// if(!curHist){
677// curHist=new THnSparseF(event->GetFiredTriggerClasses(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
678// fMapDz->Add(new TObjString(event->GetFiredTriggerClasses()),curHist);
679// }
680 curHist->Fill(vecDrift);
74235403 681
dde68d36 682 name=shortName;
683 name+="ALL";
684 name.ToUpper();
685 curHist=(THnSparseF*)fArrayDz->FindObject(name);
74235403 686 if(!curHist){
dde68d36 687 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
688 fArrayDz->AddLast(curHist);
74235403 689 }
dde68d36 690// curHist=(THnSparseF*)(fMapDz->GetValue("all"));
691// if(!curHist){
692// curHist=new THnSparseF("all","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
693// fMapDz->Add(new TObjString("all"),curHist);
694// }
695 curHist->Fill(vecDrift);
74235403 696 }
697 }
d3ce44cb 698 TTreeSRedirector *cstream = GetDebugStreamer();
699 if (fStreamLevel>0){
700 if (cstream){
701 (*cstream)<<"trackInfo"<<
702 "tr0.="<<track0<<
703 "tr1.="<<track1<<
704 "p0.="<<&param0<<
705 "p1.="<<&param1<<
cc65e4f5 706 "nAC="<<nAC<<
707 "nA0="<<nA0<<
708 "nA1="<<nA1<<
709 "nC0="<<nC0<<
710 "nC1="<<nC1<<
d3ce44cb 711 "isPair="<<isPair<<
712 "isCross="<<isCross<<
dde68d36 713 "isSame="<<isSame<<
d3ce44cb 714 "fDz="<<fDz<<
715 "fRun="<<fRun<<
716 "fTime="<<fTime<<
717 "\n";
718 }
719 }
c74c5f6c 720 } // end 2nd order loop
721 } // end 1st order loop
74235403 722
2be25cc9 723 if (fStreamLevel>0){
724 TTreeSRedirector *cstream = GetDebugStreamer();
725 if (cstream){
2be25cc9 726 (*cstream)<<"timeInfo"<<
727 "run="<<fRun<< // run number
728 "event="<<fEvent<< // event number
729 "time="<<fTime<< // time stamp of event
730 "trigger="<<fTrigger<< // trigger
731 "mag="<<fMagF<< // magnetic field
732 // Environment values
2be25cc9 733 //
734 // accumulated values
735 //
736 "fDz="<<fDz<< //! current delta z
2be25cc9 737 "trigger="<<event->GetFiredTriggerClasses()<<
738 "\n";
739 }
740 }
a390f11f 741 if (GetDebugLevel()>20) printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
c74c5f6c 742}
743
70036ef3 744void AliTPCcalibTime::ProcessBeam(const AliESDEvent *const event){
745 //
e0ef0ad7 746 // Process beam data - calculates vartex
747 // from A side and C side
748 // Histogram the differences
70036ef3 749 //
750 const Int_t kMinClusters =80;
b1fea14f 751 const Int_t kMinTracks =2; // minimal number of tracks to define the vertex
752 const Int_t kMinTracksVertex=30; // minimal number of tracks to define the cumulative vertex
70036ef3 753 const Double_t kMaxTgl =1.2; // maximal Tgl (z angle)
b1fea14f 754 const Double_t kMinPt =0.2; // minimal pt
755 const Double_t kMaxD0 =5.; // cut on distance to the primary vertex first guess
70036ef3 756 const Double_t kMaxZ0 =20;
b1fea14f 757 const Double_t kMaxD =2.5; // cut on distance to the primary vertex
758 const Double_t kMaxZ =4; // maximal z distance between tracks form the same side
759 const Double_t kMaxChi2 =15; // maximal chi2 of the TPCvertex
760 const Double_t kCumulCovarXY=0.003; //increase the error of cumul vertex 30 microns profile
761 const Double_t kCumulCovarZ=250.; //increase the error of cumul vertex
762 const Double_t kMaxDvertex = 1.0; // cut to accept the vertex;
763 //
764 Int_t flags=0;
765 const Int_t kBuffSize=100;
766 static Double_t deltaZ[kBuffSize]={0};
767 static Int_t counterZ=0;
768 static AliKFVertex cumulVertexA, cumulVertexC, cumulVertexAC; // cumulative vertex
769 AliKFVertex vertexA, vertexC;
70036ef3 770
771 Float_t dca0[2]={0,0};
772 Double_t dcaVertex[2]={0,0};
70036ef3 773 Int_t ntracks=event->GetNumberOfTracks();
774 if (ntracks==0) return;
012c4694 775 if (ntracks > fCutTracks) return;
70036ef3 776 //
777 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
778 //
779 // Divide tracks to A and C side tracks - using the cluster indexes
780 TObjArray tracksA(ntracks);
781 TObjArray tracksC(ntracks);
6ba68367 782 //
70036ef3 783 AliESDVertex *vertexSPD = (AliESDVertex *)event->GetPrimaryVertexSPD();
784 AliESDVertex *vertex = (AliESDVertex *)event->GetPrimaryVertex();
785 AliESDVertex *vertexTracks = (AliESDVertex *)event->GetPrimaryVertexTracks();
786 Double_t vertexZA[10000], vertexZC[10000];
6ba68367 787 //
70036ef3 788 Int_t ntracksA= 0;
789 Int_t ntracksC= 0;
790 //
791 for (Int_t itrack=0;itrack<ntracks;itrack++) {
792 AliESDtrack *track = event->GetTrack(itrack);
793 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(itrack);
794 if (!friendTrack) continue;
795 if (TMath::Abs(track->GetTgl())>kMaxTgl) continue;
796 if (TMath::Abs(track->Pt())<kMinPt) continue;
797 const AliExternalTrackParam * trackIn = track->GetInnerParam();
798 TObject *calibObject=0;
799 AliTPCseed *seed = 0;
800 Int_t nA=0, nC=0;
801 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
802 if (seed) {
803 for (Int_t irow=159;irow>0;irow--) {
804 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
805 if (!cl) continue;
806 if ((cl->GetDetector()%36)<18) nA++;
807 if ((cl->GetDetector()%36)>=18) nC++;
808 }
809 if ((nA>kMinClusters || nC>kMinClusters) && (nA*nC==0) ){
810 track->GetImpactParameters(dca0[0],dca0[1]);
811 if (TMath::Abs(dca0[0])>kMaxD0) continue;
812 if (TMath::Abs(dca0[1])>kMaxZ0) continue;
813 AliExternalTrackParam pTPCvertex(*trackIn);
814 if (!AliTracker::PropagateTrackToBxByBz(&pTPCvertex,4.+4.*TMath::Abs(dca0[0]),0.1,2,kTRUE)) continue;
815 pTPCvertex.PropagateToDCA(vertex,AliTracker::GetBz(), kMaxD, dcaVertex,0);
816 if (TMath::Abs(dcaVertex[0])>kMaxD) continue;
817 if (nA>kMinClusters &&nC==0) { tracksA.AddLast(pTPCvertex.Clone()); vertexZA[ntracksA++] = pTPCvertex.GetZ();}
818 if (nC>kMinClusters &&nA==0) {tracksC.AddLast(pTPCvertex.Clone()); vertexZC[ntracksC++] = pTPCvertex.GetZ();}
819 }
820 }
821 }
b1fea14f 822 Double_t medianZA=TMath::Median(ntracksA, vertexZA); // tracks median
823 Double_t medianZC=TMath::Median(ntracksC, vertexZC); // tracks median
70036ef3 824 //
825 ntracksA= tracksA.GetEntriesFast();
826 ntracksC= tracksC.GetEntriesFast();
827 if (ntracksA>kMinTracks && ntracksC>kMinTracks){
b1fea14f 828 deltaZ[counterZ%kBuffSize]=medianZA-medianZC;
70036ef3 829 counterZ+=1;
b1fea14f 830 Double_t medianDelta=(counterZ>=kBuffSize)? TMath::Median(kBuffSize, deltaZ): TMath::Median(counterZ, deltaZ);
70036ef3 831 if (TMath::Abs(medianDelta-(medianZA-medianZC))>kMaxZ) flags+=16;
832 // increse the error of cumulative vertex at the beginning of event
833 cumulVertexA.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
834 cumulVertexA.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
835 cumulVertexA.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
836 cumulVertexC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
837 cumulVertexC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
838 cumulVertexC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
839 cumulVertexAC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
840 cumulVertexAC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
841 cumulVertexAC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
842 //
843 for (Int_t iA=0; iA<ntracksA; iA++){
844 if (flags!=0) continue;
845 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksA.At(iA);
846 if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
847 AliKFParticle part(*aliTrack,211);
848 vertexA+=part;
70036ef3 849 }
850 for (Int_t iC=0; iC<ntracksC; iC++){
851 if (flags!=0) continue;
852 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksC.At(iC);
853 if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
854 AliKFParticle part(*aliTrack,211);
855 vertexC+=part;
70036ef3 856 }
b1fea14f 857 //
70036ef3 858 if (vertexA.GetNDF()<kMinTracks) flags+=32;
859 if (vertexC.GetNDF()<kMinTracks) flags+=32;
860 if (TMath::Abs(vertexA.Z()-medianZA)>kMaxZ) flags+=1; //apply cuts
861 if (TMath::Abs(vertexC.Z()-medianZC)>kMaxZ) flags+=2;
862 if (TMath::Abs(vertexA.GetChi2()/vertexA.GetNDF()+vertexC.GetChi2()/vertexC.GetNDF())> kMaxChi2) flags+=4;
b1fea14f 863 //
864 if (flags==0){
865 for (Int_t iA=0; iA<ntracksA; iA++){
866 if (flags!=0) continue;
867 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksA.At(iA);
868 if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
869 AliKFParticle part(*aliTrack,211);
870 cumulVertexA+=part;
871 cumulVertexAC+=part;
872 }
873 for (Int_t iC=0; iC<ntracksC; iC++){
874 if (flags!=0) continue;
875 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksC.At(iC);
876 if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
877 AliKFParticle part(*aliTrack,211);
878 cumulVertexC+=part;
879 cumulVertexAC+=part;
880 }
70036ef3 881 //
b1fea14f 882 if (TMath::Abs(cumulVertexA.X()-vertexA.X())>kMaxDvertex) flags+=64;
883 if (TMath::Abs(cumulVertexA.Y()-vertexA.Y())>kMaxDvertex) flags+=64;
884 if (TMath::Abs(cumulVertexA.Z()-vertexA.Z())>kMaxDvertex) flags+=64;
70036ef3 885 //
b1fea14f 886 if (TMath::Abs(cumulVertexC.X()-vertexC.X())>kMaxDvertex) flags+=64;
887 if (TMath::Abs(cumulVertexC.Y()-vertexC.Y())>kMaxDvertex) flags+=64;
888 if (TMath::Abs(cumulVertexC.Z()-vertexC.Z())>kMaxDvertex) flags+=64;
889
890
891 if ( flags==0 && cumulVertexC.GetNDF()>kMinTracksVertex&&cumulVertexA.GetNDF()>kMinTracksVertex){
892 Double_t cont[2]={0,fTime};
893 //
894 cont[0]= cumulVertexA.X();
895 fTPCVertex[0]->Fill(cont);
896 cont[0]= cumulVertexC.X();
897 fTPCVertex[1]->Fill(cont);
898 cont[0]= 0.5*(cumulVertexA.X()-cumulVertexC.X());
899 fTPCVertex[2]->Fill(cont);
900 cont[0]= 0.5*(cumulVertexA.X()+cumulVertexC.X())-vertexSPD->GetX();
901 fTPCVertex[3]->Fill(cont);
902 //
903 cont[0]= cumulVertexA.Y();
904 fTPCVertex[4]->Fill(cont);
905 cont[0]= cumulVertexC.Y();
906 fTPCVertex[5]->Fill(cont);
907 cont[0]= 0.5*(cumulVertexA.Y()-cumulVertexC.Y());
908 fTPCVertex[6]->Fill(cont);
909 cont[0]= 0.5*(cumulVertexA.Y()+cumulVertexC.Y())-vertexSPD->GetY();
910 fTPCVertex[7]->Fill(cont);
911 //
912 //
913 cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z());
914 fTPCVertex[8]->Fill(cont);
915 cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
916 fTPCVertex[9]->Fill(cont);
917 cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
918 fTPCVertex[10]->Fill(cont);
919 cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z())-vertexSPD->GetZ();
920 fTPCVertex[11]->Fill(cont);
921 //
922 Double_t correl[2]={0,0};
923 //
924 correl[0]=cumulVertexC.Z();
925 correl[1]=cumulVertexA.Z();
926 fTPCVertexCorrelation[0]->Fill(correl); // fill A side :TPC
927 correl[0]=cumulVertexA.Z();
928 correl[1]=cumulVertexC.Z();
929 fTPCVertexCorrelation[1]->Fill(correl); // fill C side :TPC
930 //
931 correl[0]=vertexSPD->GetZ();
932 correl[1]=cumulVertexA.Z()-correl[0];
933 fTPCVertexCorrelation[2]->Fill(correl); // fill A side :ITS
934 correl[1]=cumulVertexC.Z()-correl[0];
935 fTPCVertexCorrelation[3]->Fill(correl); // fill C side :ITS
936 correl[1]=0.5*(cumulVertexA.Z()+cumulVertexC.Z())-correl[0];
937 fTPCVertexCorrelation[4]->Fill(correl); // fill C side :ITS
938 }
939 }
70036ef3 940 TTreeSRedirector *cstream = GetDebugStreamer();
941 if (cstream){
942 /*
943 TCut cutChi2= "sqrt(vA.fChi2/vA.fNDF+vC.fChi2/vC.fNDF)<10"; // chi2 Cut e.g 10
944 TCut cutXY= "sqrt((vA.fP[0]-vC.fP[0])^2+(vA.fP[0]-vC.fP[1])^2)<5"; // vertex Cut
945 TCut cutZ= "abs(vA.fP[2]-mZA)<3&&abs(vC.fP[2]-mZC)<5"; // vertex Cut
946 tree->Draw("sqrt(vA.fChi2/vA.fNDF)","sqrt(vA.fChi2/vA.fNDF)<100","")
947
948 */
660f0543 949 //vertexA.Print();
950 //vertexC.Print();
70036ef3 951 (*cstream)<<"vertexTPC"<<
952 "flags="<<flags<< // rejection flags
953 "vSPD.="<<vertexSPD<< // SPD vertex
954 "vT.="<<vertexTracks<< // track vertex
955 "v.="<<vertex<< // esd vertex
956 "mZA="<<medianZA<< // median Z position at vertex A side
957 "mZC="<<medianZC<< // median Z position at vertex C side
958 "mDelta="<<medianDelta<< // median delta A side -C side
959 "counter="<<counterZ<< // counter Z
960 //
961 "vA.="<<&vertexA<< // vertex A side
962 "vC.="<<&vertexC<< // vertex C side
963 "cvA.="<<&cumulVertexA<< // cumulative vertex A side
964 "cvC.="<<&cumulVertexC<< // cumulative vertex C side
965 "cvAC.="<<&cumulVertexAC<< // cumulative vertex A+C side
966 "nA="<<ntracksA<< // contributors
967 "nC="<<ntracksC<< // contributors
968 "\n";
969 }
970 }
971 tracksA.Delete();
972 tracksC.Delete();
dde68d36 973}
974
6ba68367 975void AliTPCcalibTime::Analyze(){
976 //
977 // Special macro to analyze result of calibration and extract calibration entries
978 // Not yet ported to the Analyze function yet
979 //
980}
74235403 981
2811495d 982THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name) const
983{
6ba68367 984 //
985 // Get histogram for given trigger mask
986 //
dde68d36 987 TIterator* iterator = fArrayDz->MakeIterator();
988 iterator->Reset();
989 TString newName=name;
990 newName.ToUpper();
991 THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
992 THnSparse* addHist=NULL;
993 while((addHist=(THnSparseF*)iterator->Next())){
472f0066 994 // if(!addHist) continue;
dde68d36 995 TString histName=addHist->GetName();
996 if(!histName.Contains(newName)) continue;
997 addHist->Print();
998 newHist->Add(addHist);
74235403 999 }
dde68d36 1000 return newHist;
74235403 1001}
1002
2811495d 1003TObjArray* AliTPCcalibTime::GetHistoDrift() const
1004{
6ba68367 1005 //
1006 // return array of histograms
1007 //
dde68d36 1008 return fArrayDz;
74235403 1009}
1010
1011TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
6ba68367 1012 //
1013 // Make a drift velocity (delta Z) graph
1014 //
dde68d36 1015 THnSparse* histoDrift=GetHistoDrift(name);
1016 TGraphErrors* graphDrift=NULL;
1017 if(histoDrift){
1018 graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
1019 TString end=histoDrift->GetName();
1020 Int_t pos=end.Index("_");
1021 end=end(pos,end.Capacity()-pos);
1022 TString graphName=graphDrift->ClassName();
1023 graphName+=end;
1024 graphName.ToUpper();
1025 graphDrift->SetName(graphName);
74235403 1026 }
1027 return graphDrift;
1028}
1029
dde68d36 1030TObjArray* AliTPCcalibTime::GetGraphDrift(){
6ba68367 1031 //
1032 // make a array of drift graphs
1033 //
dde68d36 1034 TObjArray* arrayGraphDrift=new TObjArray();
1035 TIterator* iterator=fArrayDz->MakeIterator();
74235403 1036 iterator->Reset();
dde68d36 1037 THnSparse* addHist=NULL;
1038 while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
1039 return arrayGraphDrift;
74235403 1040}
1041
dde68d36 1042AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
6ba68367 1043 //
1044 // Make a fit AliSplinefit of drift velocity
1045 //
d3ce44cb 1046 TGraph* graphDrift=GetGraphDrift(name);
dde68d36 1047 AliSplineFit* fitDrift=NULL;
74235403 1048 if(graphDrift && graphDrift->GetN()){
dde68d36 1049 fitDrift=new AliSplineFit();
1050 fitDrift->SetGraph(graphDrift);
1051 fitDrift->SetMinPoints(graphDrift->GetN()+1);
1052 fitDrift->InitKnots(graphDrift,2,0,0.001);
1053 fitDrift->SplineFit(0);
1054 TString end=graphDrift->GetName();
1055 Int_t pos=end.Index("_");
1056 end=end(pos,end.Capacity()-pos);
1057 TString fitName=fitDrift->ClassName();
1058 fitName+=end;
1059 fitName.ToUpper();
1060 //fitDrift->SetName(fitName);
74235403 1061 delete graphDrift;
dde68d36 1062 graphDrift=NULL;
74235403 1063 }
1064 return fitDrift;
1065}
1066
74235403 1067
2811495d 1068Long64_t AliTPCcalibTime::Merge(TCollection *const li) {
6ba68367 1069 //
1070 // Object specific merging procedure
1071 //
c74c5f6c 1072 TIterator* iter = li->MakeIterator();
1073 AliTPCcalibTime* cal = 0;
660f0543 1074 //
c74c5f6c 1075 while ((cal = (AliTPCcalibTime*)iter->Next())) {
1076 if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
1077 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
1078 return -1;
1079 }
2be25cc9 1080 for (Int_t imeas=0; imeas<3; imeas++){
1081 if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
1082 fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
1083 fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
1084 }
1085 }
9112627b 1086 //
b1fea14f 1087 if (fTPCVertexCorrelation && cal->fTPCVertexCorrelation){
1088 for (Int_t imeas=0; imeas<5; imeas++){
1089 if (fTPCVertexCorrelation[imeas] && cal->fTPCVertexCorrelation[imeas]) fTPCVertexCorrelation[imeas]->Add(cal->fTPCVertexCorrelation[imeas]);
1090 }
1091 }
1092
472f0066 1093 if (fTPCVertex[0] && cal->fTPCVertex[0])
70036ef3 1094 for (Int_t imeas=0; imeas<12; imeas++){
660f0543 1095 if (fTPCVertex[imeas] && cal->fTPCVertex[imeas]) fTPCVertex[imeas]->Add(cal->fTPCVertex[imeas]);
70036ef3 1096 }
1097
660f0543 1098 if (fMemoryMode>0) for (Int_t imeas=0; imeas<5; imeas++){
1099 if (fMemoryMode>1){
1100 if ( cal->GetResHistoTPCCE(imeas) && cal->GetResHistoTPCCE(imeas)){
1101 fResHistoTPCCE[imeas]->Add(cal->fResHistoTPCCE[imeas]);
1102 }else{
1103 fResHistoTPCCE[imeas]=(THnSparse*)cal->fResHistoTPCCE[imeas]->Clone();
1104 }
b9908d0b 1105 }
660f0543 1106 //
1107 if ((fMemoryMode>0) &&cal->GetResHistoTPCITS(imeas) && cal->GetResHistoTPCITS(imeas)){
1108 if (fMemoryMode>1 || (imeas%2)==1) fResHistoTPCITS[imeas]->Add(cal->fResHistoTPCITS[imeas]);
1109 if (fMemoryMode>1) fResHistoTPCvertex[imeas]->Add(cal->fResHistoTPCvertex[imeas]);
33fc6b9e 1110 }
660f0543 1111 //
1112 if ((fMemoryMode>1) && cal->fResHistoTPCTRD[imeas]){
33fc6b9e 1113 if (fResHistoTPCTRD[imeas])
1114 fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
1115 else
1116 fResHistoTPCTRD[imeas]=(THnSparse*)cal->fResHistoTPCTRD[imeas]->Clone();
1117 }
660f0543 1118 //
1119 if ((fMemoryMode>1) && cal->fResHistoTPCTOF[imeas]){
33fc6b9e 1120 if (fResHistoTPCTOF[imeas])
1121 fResHistoTPCTOF[imeas]->Add(cal->fResHistoTPCTOF[imeas]);
1122 else
1123 fResHistoTPCTOF[imeas]=(THnSparse*)cal->fResHistoTPCTOF[imeas]->Clone();
9112627b 1124 }
660f0543 1125 //
97d3fea0 1126 if (cal->fArrayLaserA){
1127 fArrayLaserA->Expand(fArrayLaserA->GetEntriesFast()+cal->fArrayLaserA->GetEntriesFast());
1128 fArrayLaserC->Expand(fArrayLaserC->GetEntriesFast()+cal->fArrayLaserC->GetEntriesFast());
1129 for (Int_t ical=0; ical<cal->fArrayLaserA->GetEntriesFast(); ical++){
1130 if (cal->fArrayLaserA->UncheckedAt(ical)) fArrayLaserA->AddLast(cal->fArrayLaserA->UncheckedAt(ical)->Clone());
1131 if (cal->fArrayLaserC->UncheckedAt(ical)) fArrayLaserC->AddLast(cal->fArrayLaserC->UncheckedAt(ical)->Clone());
1132 }
1133 }
1134
9112627b 1135 }
dde68d36 1136 TObjArray* addArray=cal->GetHistoDrift();
1137 if(!addArray) return 0;
1138 TIterator* iterator = addArray->MakeIterator();
2be25cc9 1139 iterator->Reset();
dde68d36 1140 THnSparse* addHist=NULL;
660f0543 1141 if ((fMemoryMode>1)) while((addHist=(THnSparseF*)iterator->Next())){
472f0066 1142 // if(!addHist) continue;
2be25cc9 1143 addHist->Print();
dde68d36 1144 THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
2be25cc9 1145 if(!localHist){
1146 localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
dde68d36 1147 fArrayDz->AddLast(localHist);
2be25cc9 1148 }
1149 localHist->Add(addHist);
1150 }
391ffdb2 1151
d3ce44cb 1152 for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
3e55050f 1153 //
1154 // Merge alignment
1155 //
1156 for (Int_t itype=0; itype<3; itype++){
1157 //
1158 //
1159 TObjArray *arr0= 0;
1160 TObjArray *arr1= 0;
1161 if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
1162 if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
1163 if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
1164 if (!arr1) continue;
1165 if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
1166 if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
1167 arr0->Expand(arr1->GetEntriesFast());
1168 }
1169 for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
1170 AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
1171 AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
1172 if (!kalman1) continue;
1173 if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
1174 kalman0->SetRejectOutliers(kFALSE);
1175 kalman0->Merge(kalman1);
1176 }
1177 }
1178
c74c5f6c 1179 }
c74c5f6c 1180 return 0;
c74c5f6c 1181}
1182
e0ef0ad7 1183Bool_t AliTPCcalibTime::IsPair(const AliExternalTrackParam *tr0, const AliExternalTrackParam *tr1){
c74c5f6c 1184 /*
1185 // 0. Same direction - OPOSITE - cutDir +cutT
1186 TCut cutDir("cutDir","dir<-0.99")
1187 // 1.
1188 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
1189 //
1190 // 2. The same rphi
1191 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
1192 //
c74c5f6c 1193 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
1194 // 1/Pt diff cut
1195 */
1196 const Double_t *p0 = tr0->GetParameter();
1197 const Double_t *p1 = tr1->GetParameter();
74235403 1198 fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
1199 fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
d3ce44cb 1200 fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
74235403 1201 fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
1202 fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
1203
c74c5f6c 1204 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
1205 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
da6c0bc9 1206 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
c74c5f6c 1207 Double_t d0[3], d1[3];
1208 tr0->GetDirection(d0);
1209 tr1->GetDirection(d1);
1210 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
74235403 1211
d3ce44cb 1212 fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
1213 fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
1214 fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1215 fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
1216 fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
1217
c74c5f6c 1218 return kTRUE;
1219}
e0ef0ad7 1220Bool_t AliTPCcalibTime::IsCross(const AliESDtrack *const tr0, const AliESDtrack *const tr1){
cc65e4f5 1221 //
1222 // check if the cosmic pair of tracks crossed A/C side
1223 //
1224 Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
1225 if (result==kFALSE) return result;
1226 result=kTRUE;
1227 return result;
d3ce44cb 1228}
1229
e0ef0ad7 1230Bool_t AliTPCcalibTime::IsSame(const AliESDtrack *const tr0, const AliESDtrack *const tr1){
3e55050f 1231 //
1232 // track crossing the CE
1233 // 0. minimal number of clusters
1234 // 1. Same sector +-1
1235 // 2. Inner and outer track param on opposite side
1236 // 3. Outer and inner track parameter close each to other
1237 // 3.
1238 Bool_t result=kTRUE;
1239 //
1240 // inner and outer on opposite sides in z
1241 //
1242 const Int_t knclCut0 = 30;
1243 const Double_t kalphaCut = 0.4;
1244 //
1245 // 0. minimal number of clusters
1246 //
1247 if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
1248 if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
1249 //
1250 // 1. alpha cut - sector+-1
1251 //
1252 if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
1253 //
1254 // 2. Z crossing
1255 //
1256 if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
1257 if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
1258 if (result==kFALSE){
1259 return result;
1260 }
1261 //
1262 //
1263 const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
1264 const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
1265 const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1266 const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1267 //
1268 if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD) result&=kFALSE;
1269 if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1270 if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1271 if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1272 if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD) result&=kFALSE;
1273 if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1274 if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1275 if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1276 if (result==kTRUE){
1277 result=kTRUE; // just to put break point here
1278 }
1279 return result;
dde68d36 1280}
1281
31aa7c5c 1282
e0ef0ad7 1283void AliTPCcalibTime::ProcessSame(const AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
3e55050f 1284 //
1285 // Process TPC tracks crossing CE
1286 //
1287 // 0. Select only track crossing the CE
1288 // 1. Cut on the track length
b9908d0b 1289 // 2. Refit the the track on A and C side separatelly
3e55050f 1290 // 3. Fill time histograms
1291 const Int_t kMinNcl=100;
1292 const Int_t kMinNclS=25; // minimul number of clusters on the sides
b9908d0b 1293 const Double_t pimass=TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1294 const Double_t kMaxDy=1; // maximal distance in y
1295 const Double_t kMaxDsnp=0.05; // maximal distance in snp
1296 const Double_t kMaxDtheta=0.05; // maximal distance in theta
1297
3e55050f 1298 if (!friendTrack->GetTPCOut()) return;
1299 //
1300 // 0. Select only track crossing the CE
1301 //
1302 if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1303 //
1304 // 1. cut on track length
1305 //
1306 if (track->GetTPCNcls()<kMinNcl) return;
1307 //
1308 // 2. Refit track sepparatel on A and C side
1309 //
1310 TObject *calibObject;
1311 AliTPCseed *seed = 0;
1312 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1313 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1314 }
1315 if (!seed) return;
1316 //
1317 AliExternalTrackParam trackIn(*track->GetInnerParam());
1318 AliExternalTrackParam trackOut(*track->GetOuterParam());
1319 Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1320 Double_t xyz[3]={0,0.,0.0};
1321 Double_t bz =0;
1322 Int_t nclIn=0,nclOut=0;
70036ef3 1323 trackIn.ResetCovariance(1000.);
1324 trackOut.ResetCovariance(1000.);
3e55050f 1325 //
1326 //2.a Refit inner
1327 //
b9908d0b 1328 Int_t sideIn=0;
3e55050f 1329 for (Int_t irow=0;irow<159;irow++) {
1330 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1331 if (!cl) continue;
1332 if (cl->GetX()<80) continue;
b9908d0b 1333 if (sideIn==0){
1334 if (cl->GetDetector()%36<18) sideIn=1;
1335 if (cl->GetDetector()%36>=18) sideIn=-1;
1336 }
1337 if (sideIn== -1 && (cl->GetDetector()%36)<18) break;
1338 if (sideIn== 1 &&(cl->GetDetector()%36)>=18) break;
3e55050f 1339 Int_t sector = cl->GetDetector();
1340 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1341 if (TMath::Abs(dalpha)>0.01){
1342 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1343 }
1344 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1345 trackIn.GetXYZ(xyz);
1346 bz = AliTracker::GetBz(xyz);
b9908d0b 1347 AliTracker::PropagateTrackToBxByBz(&trackIn,r[0],1.,pimass,kFALSE);
3e55050f 1348 if (!trackIn.PropagateTo(r[0],bz)) break;
1349 nclIn++;
1350 trackIn.Update(&r[1],cov);
1351 }
1352 //
1353 //2.b Refit outer
b9908d0b 1354 //
1355 Int_t sideOut=0;
3e55050f 1356 for (Int_t irow=159;irow>0;irow--) {
1357 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1358 if (!cl) continue;
1359 if (cl->GetX()<80) continue;
b9908d0b 1360 if (sideOut==0){
1361 if (cl->GetDetector()%36<18) sideOut=1;
1362 if (cl->GetDetector()%36>=18) sideOut=-1;
1363 if (sideIn==sideOut) break;
1364 }
1365 if (sideOut== -1 && (cl->GetDetector()%36)<18) break;
1366 if (sideOut== 1 &&(cl->GetDetector()%36)>=18) break;
1367 //
3e55050f 1368 Int_t sector = cl->GetDetector();
1369 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1370 if (TMath::Abs(dalpha)>0.01){
1371 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1372 }
1373 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1374 trackOut.GetXYZ(xyz);
1375 bz = AliTracker::GetBz(xyz);
b9908d0b 1376 AliTracker::PropagateTrackToBxByBz(&trackOut,r[0],1.,pimass,kFALSE);
3e55050f 1377 if (!trackOut.PropagateTo(r[0],bz)) break;
1378 nclOut++;
1379 trackOut.Update(&r[1],cov);
1380 }
1381 trackOut.Rotate(trackIn.GetAlpha());
1382 Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1383 trackIn.PropagateTo(meanX,bz);
1384 trackOut.PropagateTo(meanX,bz);
b9908d0b 1385 if (TMath::Abs(trackIn.GetY()-trackOut.GetY())>kMaxDy) return;
1386 if (TMath::Abs(trackIn.GetSnp()-trackOut.GetSnp())>kMaxDsnp) return;
1387 if (TMath::Abs(trackIn.GetTgl()-trackOut.GetTgl())>kMaxDtheta) return;
1388 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1389 FillResHistoTPCCE(&trackIn,&trackOut);
1390 }
3e55050f 1391 TTreeSRedirector *cstream = GetDebugStreamer();
1392 if (cstream){
1393 TVectorD gxyz(3);
1394 trackIn.GetXYZ(gxyz.GetMatrixArray());
1395 TTimeStamp tstamp(fTime);
3e55050f 1396 (*cstream)<<"tpctpc"<<
1397 "run="<<fRun<< // run number
1398 "event="<<fEvent<< // event number
1399 "time="<<fTime<< // time stamp of event
1400 "trigger="<<fTrigger<< // trigger
1401 "mag="<<fMagF<< // magnetic field
3e55050f 1402 //
b9908d0b 1403 "sideIn="<<sideIn<< // side at inner part
1404 "sideOut="<<sideOut<< // side at puter part
3e55050f 1405 "xyz.="<<&gxyz<< // global position
1406 "tIn.="<<&trackIn<< // refitterd track in
1407 "tOut.="<<&trackOut<< // refitter track out
1408 "nclIn="<<nclIn<< //
1409 "nclOut="<<nclOut<< //
1410 "\n";
1411 }
1412 //
1413 // 3. Fill time histograms
1414 // Debug stremaer expression
1415 // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1416 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1417 fDz = trackOut.GetZ()-trackIn.GetZ();
1418 if (trackOut.GetTgl()<0) fDz*=-1.;
1419 TTimeStamp tstamp(fTime);
1420 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1421 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1422 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1423 //
1424 // fill histograms per trigger class and itegrated
1425 //
1426 THnSparse* curHist=NULL;
1427 for (Int_t itype=0; itype<2; itype++){
1428 TString name="MEAN_VDRIFT_CROSS_";
1429 if (itype==0){
1430 name+=event->GetFiredTriggerClasses();
1431 name.ToUpper();
1432 }else{
1433 name+="ALL";
1434 }
1435 curHist=(THnSparseF*)fArrayDz->FindObject(name);
1436 if(!curHist){
1437 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1438 fArrayDz->AddLast(curHist);
1439 }
1440 curHist->Fill(vecDrift);
1441 }
1442 }
1443
1444}
817766d5 1445
e0ef0ad7 1446void AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
817766d5 1447 //
3e55050f 1448 // Process track - Update TPC-ITS alignment
1449 // Updates:
1450 // 0. Apply standartd cuts
1451 // 1. Recalucluate the current statistic median/RMS
1452 // 2. Apply median+-rms cut
1453 // 3. Update kalman filter
1454 //
1455 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1456 const Int_t kMinITS = 3; // minimal number of ITS cluster
1457 const Double_t kMinZ = 10; // maximal dz distance
b96c3aef 1458 const Double_t kMaxDy = 2.; // maximal dy distance
b1fea14f 1459 const Double_t kMaxAngle= 0.07; // maximal angular distance
3e55050f 1460 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1461 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
97d3fea0 1462 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
3e55050f 1463 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1464 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
b96c3aef 1465 const Double_t kMinPt = 0.3; // minimal pt
b1fea14f 1466 const Double_t kMax1Pt=0.5; //maximal 1/pt distance
be67055b 1467 const Int_t kN=50; // deepnes of history
3e55050f 1468 static Int_t kglast=0;
1469 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
817766d5 1470 //
3e55050f 1471 // 0. Apply standard cuts
b1fea14f 1472 //
817766d5 1473 Int_t dummycl[1000];
817766d5 1474 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
b96c3aef 1475 if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
817766d5 1476 if (!track->GetInnerParam()) return;
1477 if (!track->GetOuterParam()) return;
b842d904 1478 if (track->GetInnerParam()->Pt()<kMinPt) return;
817766d5 1479 // exclude crossing track
1480 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
70036ef3 1481 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ/3.) return;
b96c3aef 1482 if (track->GetInnerParam()->GetX()>90) return;
817766d5 1483 //
1484 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
b1fea14f 1485 //
1486 AliExternalTrackParam pITS; // ITS standalone if possible
1487 AliExternalTrackParam pITS2; //TPC-ITS track
1488 if (friendTrack->GetITSOut()){
1489 pITS2=(*(friendTrack->GetITSOut())); //TPC-ITS track - snapshot ITS out
1490 pITS2.Rotate(pTPC.GetAlpha());
1491 AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1492 }
be67055b 1493
b842d904 1494 AliESDfriendTrack *itsfriendTrack=0;
1495 //
1496 // try to find standalone ITS track corresponing to the TPC if possible
1497 //
1498 Bool_t hasAlone=kFALSE;
1499 Int_t ntracks=event->GetNumberOfTracks();
1500 for (Int_t i=0; i<ntracks; i++){
b1fea14f 1501 AliESDtrack * trackITS = event->GetTrack(i);
1502 if (!trackITS) continue;
1503 if (trackITS->GetITSclusters(dummycl)<kMinITS) continue; // minimal amount of clusters
6ba68367 1504 itsfriendTrack = esdFriend->GetTrack(i);
b842d904 1505 if (!itsfriendTrack) continue;
1506 if (!itsfriendTrack->GetITSOut()) continue;
b1fea14f 1507
1508 if (TMath::Abs(pTPC.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1509 if (TMath::Abs(pTPC.GetSigned1Pt()-itsfriendTrack->GetITSOut()->GetSigned1Pt())> kMax1Pt) continue;
b842d904 1510 pITS=(*(itsfriendTrack->GetITSOut()));
1511 //
1512 pITS.Rotate(pTPC.GetAlpha());
be67055b 1513 AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
b1fea14f 1514 if (TMath::Abs(pTPC.GetY()-pITS.GetY())> kMaxDy) continue;
1515 if (TMath::Abs(pTPC.GetSnp()-pITS.GetSnp())> kMaxAngle) continue;
b842d904 1516 hasAlone=kTRUE;
1517 }
b1fea14f 1518 if (!hasAlone) {
1519 if (track->GetITSclusters(dummycl)<kMinITS) return;
1520 pITS=pITS2; // use combined track if it has ITS
1521 }
b842d904 1522 //
3e55050f 1523 if (TMath::Abs(pITS.GetY()-pTPC.GetY()) >kMaxDy) return;
1524 if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1525 if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1526 //
1527 // 1. Update median and RMS info
1528 //
b322e06a 1529 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
3e55050f 1530 TVectorD vecDeltaN(5);
1531 Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1532 vecDelta[4]=0;
1533 for (Int_t i=0;i<4;i++){
1534 vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1535 kgdP[i][kglast%kN]=vecDelta[i];
1536 }
1537 kglast=(kglast+1);
1538 Int_t entries=(kglast<kN)?kglast:kN;
1539 for (Int_t i=0;i<4;i++){
1540 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1541 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1542 vecDeltaN[i] = 0;
1543 if (vecRMS[i]>0.){
1544 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1545 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1546 }
1547 }
817766d5 1548 //
3e55050f 1549 // 2. Apply median+-rms cut
817766d5 1550 //
3e55050f 1551 if (kglast<3) return; //median and RMS to be defined
1552 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1553 //
1554 // 3. Update alignment
817766d5 1555 //
97d3fea0 1556 Int_t htime = (fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time bins number
4486a91f 1557 if (fAlignITSTPC->GetEntriesFast()<htime){
817766d5 1558 fAlignITSTPC->Expand(htime*2+20);
1559 }
1560 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1561 if (!align){
3e55050f 1562 // make Alignment object if doesn't exist
817766d5 1563 align=new AliRelAlignerKalman();
3e55050f 1564 align->SetRunNumber(fRun);
1565 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
97d3fea0 1566 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
3e55050f 1567 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1568 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1569 align->SetRejectOutliers(kFALSE);
1570
1571 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
817766d5 1572 align->SetMagField(fMagF);
1573 fAlignITSTPC->AddAt(align,htime);
1574 }
817766d5 1575 align->AddTrackParams(&pITS,&pTPC);
5647625c 1576 Double_t averageTime = fTime;
1577 if (align->GetTimeStamp()>0&&align->GetNUpdates()>0){
97d3fea0 1578 averageTime=((Double_t(align->GetTimeStamp())*Double_t(align->GetNUpdates())+Double_t(fTime)))/(Double_t(align->GetNUpdates())+1.);
5647625c 1579 }
1580 align->SetTimeStamp(Int_t(averageTime));
1581
3e55050f 1582 align->SetRunNumber(fRun );
391ffdb2 1583 Float_t dca[2],cov[3];
1584 track->GetImpactParameters(dca,cov);
882f5c06 1585 if (TMath::Abs(dca[0])<kMaxDy){
391ffdb2 1586 FillResHistoTPCITS(&pTPC,&pITS);
1587 FillResHistoTPC(track);
1588 }
3e55050f 1589 //
1590 Int_t nupdates=align->GetNUpdates();
1591 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1592 align->SetRejectOutliers(kFALSE);
817766d5 1593 TTreeSRedirector *cstream = GetDebugStreamer();
1594 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
817766d5 1595 TVectorD gpTPC(3), gdTPC(3);
1596 TVectorD gpITS(3), gdITS(3);
1597 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1598 pTPC.GetDirection(gdTPC.GetMatrixArray());
1599 pITS.GetXYZ(gpITS.GetMatrixArray());
1600 pITS.GetDirection(gdITS.GetMatrixArray());
1601 (*cstream)<<"itstpc"<<
1602 "run="<<fRun<< // run number
1603 "event="<<fEvent<< // event number
1604 "time="<<fTime<< // time stamp of event
1605 "trigger="<<fTrigger<< // trigger
1606 "mag="<<fMagF<< // magnetic field
817766d5 1607 //
b842d904 1608 "hasAlone="<<hasAlone<< // has ITS standalone ?
1609 "track.="<<track<< // track info
3e55050f 1610 "nmed="<<kglast<< // number of entries to define median and RMS
1611 "vMed.="<<&vecMedian<< // median of deltas
1612 "vRMS.="<<&vecRMS<< // rms of deltas
1613 "vDelta.="<<&vecDelta<< // delta in respect to median
1614 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
3e55050f 1615 "a.="<<align<< // current alignment
b842d904 1616 "pITS.="<<&pITS<< // track param ITS
1617 "pITS2.="<<&pITS2<< // track param ITS+TPC
3e55050f 1618 "pTPC.="<<&pTPC<< // track param TPC
1619 "gpTPC.="<<&gpTPC<< // global position TPC
1620 "gdTPC.="<<&gdTPC<< // global direction TPC
1621 "gpITS.="<<&gpITS<< // global position ITS
1622 "gdITS.="<<&gdITS<< // global position ITS
817766d5 1623 "\n";
1624 }
817766d5 1625}
3e55050f 1626
1627
1628
1629
e0ef0ad7 1630void AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack){
0dac7b3a 1631 //
3e55050f 1632 // Process track - Update TPC-TRD alignment
1633 // Updates:
1634 // 0. Apply standartd cuts
1635 // 1. Recalucluate the current statistic median/RMS
1636 // 2. Apply median+-rms cut
1637 // 3. Update kalman filter
1638 //
1639 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1640 const Int_t kMinTRD = 50; // minimal number of TRD cluster
1641 const Double_t kMinZ = 20; // maximal dz distance
9112627b 1642 const Double_t kMaxDy = 5.; // maximal dy distance
1643 const Double_t kMaxAngle= 0.1; // maximal angular distance
1644 const Double_t kSigmaCut= 10; // maximal sigma distance to median
3e55050f 1645 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
97d3fea0 1646 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
3e55050f 1647 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1648 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
be67055b 1649 const Double_t kRefX = 275; // reference X
1650 const Int_t kN=50; // deepnes of history
3e55050f 1651 static Int_t kglast=0;
1652 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
0dac7b3a 1653 //
3e55050f 1654 // 0. Apply standard cuts
0dac7b3a 1655 //
1656 Int_t dummycl[1000];
1657 if (track->GetTRDclusters(dummycl)<kMinTRD) return; // minimal amount of clusters
1658 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
3e55050f 1659 if (!friendTrack->GetTRDIn()) return;
be67055b 1660 if (!track->IsOn(AliESDtrack::kTRDrefit)) return;
1661 if (!track->IsOn(AliESDtrack::kTRDout)) return;
0dac7b3a 1662 if (!track->GetInnerParam()) return;
be67055b 1663 if (!friendTrack->GetTPCOut()) return;
0dac7b3a 1664 // exclude crossing track
be67055b 1665 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
3e55050f 1666 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
0dac7b3a 1667 //
be67055b 1668 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1669 AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
0dac7b3a 1670 AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
3e55050f 1671 pTRD.Rotate(pTPC.GetAlpha());
be67055b 1672 // pTRD.PropagateTo(pTPC.GetX(),fMagF);
1673 AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1674
3e55050f 1675 ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.; // increas sys errors
1676 ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1677
1678 if (TMath::Abs(pTRD.GetY()-pTPC.GetY()) >kMaxDy) return;
1679 if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
9112627b 1680 // if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
0dac7b3a 1681 //
3e55050f 1682 // 1. Update median and RMS info
0dac7b3a 1683 //
b322e06a 1684 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
3e55050f 1685 TVectorD vecDeltaN(5);
1686 Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1687 vecDelta[4]=0;
1688 for (Int_t i=0;i<4;i++){
1689 vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1690 kgdP[i][kglast%kN]=vecDelta[i];
1691 }
1692 kglast=(kglast+1);
1693 Int_t entries=(kglast<kN)?kglast:kN;
1694 for (Int_t i=0;i<4;i++){
1695 vecMedian[i] = TMath::Median(entries,kgdP[i]);
be67055b 1696
3e55050f 1697 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1698 vecDeltaN[i] = 0;
1699 if (vecRMS[i]>0.){
1700 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1701 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1702 }
1703 }
1704 //
1705 // 2. Apply median+-rms cut
1706 //
1707 if (kglast<3) return; //median and RMS to be defined
1708 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1709 //
1710 // 3. Update alignment
0dac7b3a 1711 //
5647625c 1712 //Int_t htime = fTime/3600; //time in hours
97d3fea0 1713 Int_t htime = (Int_t)(fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time in half hour
4486a91f 1714 if (fAlignTRDTPC->GetEntriesFast()<htime){
0dac7b3a 1715 fAlignTRDTPC->Expand(htime*2+20);
1716 }
1717 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1718 if (!align){
3e55050f 1719 // make Alignment object if doesn't exist
0dac7b3a 1720 align=new AliRelAlignerKalman();
3e55050f 1721 align->SetRunNumber(fRun);
1722 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
97d3fea0 1723 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
3e55050f 1724 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1725 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1726 align->SetRejectOutliers(kFALSE);
1727 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
0dac7b3a 1728 align->SetMagField(fMagF);
1729 fAlignTRDTPC->AddAt(align,htime);
1730 }
0dac7b3a 1731 align->AddTrackParams(&pTRD,&pTPC);
5647625c 1732 //align->SetTimeStamp(fTime);
1733 Double_t averageTime = fTime;
1734 if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1735 averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1736 //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1737 }
1738 align->SetTimeStamp((Int_t)averageTime);
1739
1740 //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1741
3e55050f 1742 align->SetRunNumber(fRun );
0b736a46 1743 Float_t dca[2],cov[3];
1744 track->GetImpactParameters(dca,cov);
1745 if (TMath::Abs(dca[0])<kMaxDy){
1746 FillResHistoTPCTRD(&pTPC,&pTRD); //only primaries
1747 }
3e55050f 1748 //
1749 Int_t nupdates=align->GetNUpdates();
1750 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1751 align->SetRejectOutliers(kFALSE);
0dac7b3a 1752 TTreeSRedirector *cstream = GetDebugStreamer();
1753 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
0dac7b3a 1754 TVectorD gpTPC(3), gdTPC(3);
1755 TVectorD gpTRD(3), gdTRD(3);
1756 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1757 pTPC.GetDirection(gdTPC.GetMatrixArray());
1758 pTRD.GetXYZ(gpTRD.GetMatrixArray());
1759 pTRD.GetDirection(gdTRD.GetMatrixArray());
3e55050f 1760 (*cstream)<<"trdtpc"<<
0dac7b3a 1761 "run="<<fRun<< // run number
1762 "event="<<fEvent<< // event number
1763 "time="<<fTime<< // time stamp of event
1764 "trigger="<<fTrigger<< // trigger
1765 "mag="<<fMagF<< // magnetic field
0dac7b3a 1766 //
3e55050f 1767 "nmed="<<kglast<< // number of entries to define median and RMS
1768 "vMed.="<<&vecMedian<< // median of deltas
1769 "vRMS.="<<&vecRMS<< // rms of deltas
1770 "vDelta.="<<&vecDelta<< // delta in respect to median
1771 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1772 "t.="<<track<< // ful track - find proper cuts
1773 "a.="<<align<< // current alignment
1774 "pTRD.="<<&pTRD<< // track param TRD
1775 "pTPC.="<<&pTPC<< // track param TPC
1776 "gpTPC.="<<&gpTPC<< // global position TPC
1777 "gdTPC.="<<&gdTPC<< // global direction TPC
1778 "gpTRD.="<<&gpTRD<< // global position TRD
1779 "gdTRD.="<<&gdTRD<< // global position TRD
0dac7b3a 1780 "\n";
1781 }
0dac7b3a 1782}
817766d5 1783
1784
e0ef0ad7 1785void AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack){
817766d5 1786 //
817766d5 1787 //
3e55050f 1788 // Process track - Update TPC-TOF alignment
1789 // Updates:
1790 // -1. Make a TOF "track"
1791 // 0. Apply standartd cuts
1792 // 1. Recalucluate the current statistic median/RMS
1793 // 2. Apply median+-rms cut
1794 // 3. Update kalman filter
1795 //
1796 const Int_t kMinTPC = 80; // minimal number of TPC cluster
b842d904 1797 // const Double_t kMinZ = 10; // maximal dz distance
3e55050f 1798 const Double_t kMaxDy = 5.; // maximal dy distance
70036ef3 1799 const Double_t kMaxAngle= 0.05; // maximal angular distance
3e55050f 1800 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1801 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
97d3fea0 1802 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
3e55050f 1803 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1804
1805 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
be67055b 1806 const Int_t kN=50; // deepnes of history
3e55050f 1807 static Int_t kglast=0;
1808 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
817766d5 1809 //
3e55050f 1810 // -1. Make a TOF track-
1811 // Clusters are not in friends - use alingment points
1812 //
1813 if (track->GetTOFsignal()<=0) return;
1814 if (!friendTrack->GetTPCOut()) return;
1815 if (!track->GetInnerParam()) return;
be67055b 1816 if (!friendTrack->GetTPCOut()) return;
817766d5 1817 const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1818 if (!points) return;
be67055b 1819 AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
3e55050f 1820 AliExternalTrackParam pTOF(pTPC);
1821 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
817766d5 1822 Int_t npoints = points->GetNPoints();
1823 AliTrackPoint point;
3e55050f 1824 Int_t naccept=0;
817766d5 1825 //
1826 for (Int_t ipoint=0;ipoint<npoints;ipoint++){
817766d5 1827 points->GetPoint(point,ipoint);
1828 Float_t xyz[3];
1829 point.GetXYZ(xyz);
1830 Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
3e55050f 1831 if (r<350) continue;
817766d5 1832 if (r>400) continue;
3e55050f 1833 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1834 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);
1835 AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1836 pTPC.PropagateTo(lpoint.GetX(),fMagF);
1837 pTOF=pTPC;
817766d5 1838 ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1839 ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
3e55050f 1840 ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1841 ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
0dac7b3a 1842 ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1843 ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
3e55050f 1844 naccept++;
1845 }
1846 if (naccept==0) return; // no tof match clusters
1847 //
1848 // 0. Apply standard cuts
1849 //
1850 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1851 // exclude crossing track
be67055b 1852 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
3e55050f 1853 //
1854 if (TMath::Abs(pTOF.GetY()-pTPC.GetY()) >kMaxDy) return;
1855 if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1856 if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1857 //
1858 // 1. Update median and RMS info
1859 //
b322e06a 1860 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
3e55050f 1861 TVectorD vecDeltaN(5);
1862 Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1863 vecDelta[4]=0;
1864 for (Int_t i=0;i<4;i++){
1865 vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1866 kgdP[i][kglast%kN]=vecDelta[i];
1867 }
1868 kglast=(kglast+1);
1869 Int_t entries=(kglast<kN)?kglast:kN;
1870 Bool_t isOK=kTRUE;
1871 for (Int_t i=0;i<4;i++){
1872 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1873 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1874 vecDeltaN[i] = 0;
1875 if (vecRMS[i]>0.){
1876 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1877 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1878 if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
817766d5 1879 }
3e55050f 1880 }
1881 //
1882 // 2. Apply median+-rms cut
1883 //
1884 if (kglast<10) return; //median and RMS to be defined
1885 if (!isOK) return;
1886 //
1887 // 3. Update alignment
1888 //
5647625c 1889 //Int_t htime = fTime/3600; //time in hours
97d3fea0 1890 Int_t htime = (Int_t)(fTime-fTimeKalmanBin)/fTimeKalmanBin; //time bin
4486a91f 1891 if (fAlignTOFTPC->GetEntriesFast()<htime){
3e55050f 1892 fAlignTOFTPC->Expand(htime*2+20);
1893 }
1894 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1895 if (!align){
1896 // make Alignment object if doesn't exist
1897 align=new AliRelAlignerKalman();
1898 align->SetRunNumber(fRun);
1899 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
97d3fea0 1900 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
3e55050f 1901 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1902 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1903 align->SetRejectOutliers(kFALSE);
1904 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1905 align->SetMagField(fMagF);
1906 fAlignTOFTPC->AddAt(align,htime);
1907 }
1908 align->AddTrackParams(&pTOF,&pTPC);
0b736a46 1909 Float_t dca[2],cov[3];
1910 track->GetImpactParameters(dca,cov);
1911 if (TMath::Abs(dca[0])<kMaxDy){
1912 FillResHistoTPCTOF(&pTPC,&pTOF);
1913 }
5647625c 1914 //align->SetTimeStamp(fTime);
1915 Double_t averageTime = fTime;
1916 if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1917 averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1918 //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1919 }
1920 align->SetTimeStamp((Int_t)averageTime);
1921
1922 //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1923
3e55050f 1924 align->SetRunNumber(fRun );
1925 //
1926 Int_t nupdates=align->GetNUpdates();
1927 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1928 align->SetRejectOutliers(kFALSE);
1929 TTreeSRedirector *cstream = GetDebugStreamer();
1930 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
3e55050f 1931 TVectorD gpTPC(3), gdTPC(3);
1932 TVectorD gpTOF(3), gdTOF(3);
1933 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1934 pTPC.GetDirection(gdTPC.GetMatrixArray());
1935 pTOF.GetXYZ(gpTOF.GetMatrixArray());
1936 pTOF.GetDirection(gdTOF.GetMatrixArray());
1937 (*cstream)<<"toftpc"<<
1938 "run="<<fRun<< // run number
1939 "event="<<fEvent<< // event number
1940 "time="<<fTime<< // time stamp of event
1941 "trigger="<<fTrigger<< // trigger
1942 "mag="<<fMagF<< // magnetic field
3e55050f 1943 //
1944 "nmed="<<kglast<< // number of entries to define median and RMS
1945 "vMed.="<<&vecMedian<< // median of deltas
1946 "vRMS.="<<&vecRMS<< // rms of deltas
1947 "vDelta.="<<&vecDelta<< // delta in respect to median
1948 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1949 "t.="<<track<< // ful track - find proper cuts
1950 "a.="<<align<< // current alignment
1951 "pTOF.="<<&pTOF<< // track param TOF
1952 "pTPC.="<<&pTPC<< // track param TPC
1953 "gpTPC.="<<&gpTPC<< // global position TPC
1954 "gdTPC.="<<&gdTPC<< // global direction TPC
1955 "gpTOF.="<<&gpTOF<< // global position TOF
1956 "gdTOF.="<<&gdTOF<< // global position TOF
1957 "\n";
817766d5 1958 }
817766d5 1959}
3e55050f 1960
1961
391ffdb2 1962void AliTPCcalibTime::BookDistortionMaps(){
1963 //
1964 // Book ndimensional histograms of distortions/residuals
1965 // Only primary tracks are selected for analysis
1966 //
1967
b9908d0b 1968 Double_t xminTrack[5], xmaxTrack[5];
1969 Int_t binsTrack[5];
1970 TString axisName[5];
1971 TString axisTitle[5];
391ffdb2 1972 //
882f5c06 1973 binsTrack[0] =50;
391ffdb2 1974 axisName[0] ="#Delta";
1975 axisTitle[0] ="#Delta";
1976 //
70036ef3 1977 binsTrack[1] =44;
1978 xminTrack[1] =-1.1; xmaxTrack[1]=1.1;
391ffdb2 1979 axisName[1] ="tanTheta";
1980 axisTitle[1] ="tan(#Theta)";
1981 //
70036ef3 1982 binsTrack[2] =180;
391ffdb2 1983 xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi();
1984 axisName[2] ="phi";
1985 axisTitle[2] ="#phi";
1986 //
1987 binsTrack[3] =20;
1988 xminTrack[3] =-1.; xmaxTrack[3]=1.; // 0.33 GeV cut
1989 axisName[3] ="snp";
b9908d0b 1990 axisTitle[3] ="snp";
1991 //
1992 binsTrack[4] =10;
1993 xminTrack[4] =120.; xmaxTrack[4]=215.; // crossing radius for CE only
1994 axisName[4] ="r";
1995 axisTitle[4] ="r(cm)";
391ffdb2 1996 //
1997 // delta y
be67055b 1998 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
b9908d0b 1999 fResHistoTPCCE[0] = new THnSparseS("TPCCE#Delta_{Y} (cm)","#Delta_{Y} (cm)", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2000 fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2001 fResHistoTPCvertex[0] = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
be67055b 2002 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
391ffdb2 2003 fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2004 xminTrack[0] =-5; xmaxTrack[0]=5; //
2005 fResHistoTPCTOF[0] = new THnSparseS("TPCTOF#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2006 //
2007 // delta z
9112627b 2008 xminTrack[0] =-3.; xmaxTrack[0]=3.; //
b9908d0b 2009 fResHistoTPCCE[1] = new THnSparseS("TPCCE#Delta_{Z} (cm)","#Delta_{Z} (cm)", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2010 fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2011 fResHistoTPCvertex[1] = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2012 fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2013 xminTrack[0] =-5.; xmaxTrack[0]=5.; //
2014 fResHistoTPCTOF[1] = new THnSparseS("TPCTOF#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2015 //
2016 // delta snp-P2
882f5c06 2017 xminTrack[0] =-0.015; xmaxTrack[0]=0.015; //
b9908d0b 2018 fResHistoTPCCE[2] = new THnSparseS("TPCCE#Delta_{#phi}","#Delta_{#phi}", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2019 fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
9112627b 2020 fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2021 fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2022 fResHistoTPCTOF[2] = new THnSparseS("TPCTOF#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2023 //
2024 // delta theta-P3
9112627b 2025 xminTrack[0] =-0.025; xmaxTrack[0]=0.025; //
b9908d0b 2026 fResHistoTPCCE[3] = new THnSparseS("TPCCE#Delta_{#theta}","#Delta_{#theta}", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2027 fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
9112627b 2028 fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2029 fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2030 fResHistoTPCTOF[3] = new THnSparseS("TPCTOF#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2031 //
2032 // delta theta-P4
9112627b 2033 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
b9908d0b 2034 fResHistoTPCCE[4] = new THnSparseS("TPCCE#Delta_{1/pt}","#Delta_{1/pt}", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2035 fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
9112627b 2036 fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2037 fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2038 fResHistoTPCTOF[4] = new THnSparseS("TPCTOF#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2039 //
2040 for (Int_t ivar=0;ivar<4;ivar++){
b9908d0b 2041 for (Int_t ivar2=0;ivar2<5;ivar2++){
2042 fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2043 fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2044 if (ivar2<4){
2045 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2046 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2047 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2048 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2049 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2050 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2051 }
391ffdb2 2052 }
2053 }
70036ef3 2054 //
2055 // Book vertex: time histograms
2056 //
2057 Int_t binsVertex[2]={500, fTimeBins};
2058 Double_t aminVertex[2]={-5,fTimeStart};
2059 Double_t amaxVertex[2]={5, fTimeEnd};
2060 const char* hnames[12]={"TPCXAside", "TPCXCside","TPCXACdiff","TPCXAPCdiff",
2061 "TPCYAside", "TPCYCside","TPCYACdiff","TPCYAPCdiff",
2062 "TPCZAPCside", "TPCZAMCside","TPCZACdiff","TPCZAPCdiff"};
2063 const char* anames[12]={"x (cm) - A side ", "x (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{x} (cm) - TPC-Common",
2064 "y (cm) - A side ", "y (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{y} (cm) - TPC-Common",
2065 "z (cm)", "#Delta_{Z} (cm) A-C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{Z} (cm) TPC-common"};
2066 for (Int_t ihis=0; ihis<12; ihis++) {
2067 if (ihis>=8) aminVertex[0]=-20.;
2068 if (ihis>=8) amaxVertex[0]=20.;
2069 fTPCVertex[ihis]=new THnSparseF(hnames[ihis],hnames[ihis],2,binsVertex,aminVertex,amaxVertex);
2070 fTPCVertex[ihis]->GetAxis(1)->SetTitle("Time");
2071 fTPCVertex[ihis]->GetAxis(0)->SetTitle(anames[ihis]);
2072 }
b1fea14f 2073
2074 Int_t binsVertexC[2]={40, 300};
2075 Double_t aminVertexC[2]={-20,-30};
2076 Double_t amaxVertexC[2]={20,30};
2077 const char* hnamesC[5]={"TPCA_TPC","TPCC_TPC","TPCA_ITS","TPCC_ITS","TPC_ITS"};
2078 for (Int_t ihis=0; ihis<5; ihis++) {
2079 fTPCVertexCorrelation[ihis]=new THnSparseF(hnamesC[ihis],hnamesC[ihis],2,binsVertexC,aminVertexC,amaxVertexC);
2080 fTPCVertexCorrelation[ihis]->GetAxis(1)->SetTitle("z (cm)");
2081 fTPCVertexCorrelation[ihis]->GetAxis(0)->SetTitle("z (cm)");
2082 }
391ffdb2 2083}
2084
2085
b9908d0b 2086void AliTPCcalibTime::FillResHistoTPCCE(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pTPCOut ){
2087 //
2088 // fill residual histograms pTPCOut-pTPCin - trac crossing CE
2089 // Histogram
2090 //
660f0543 2091 if (fMemoryMode<2) return;
b9908d0b 2092 Double_t histoX[5];
2093 Double_t xyz[3];
2094 pTPCIn->GetXYZ(xyz);
2095 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2096 histoX[1]= pTPCIn->GetTgl();
2097 histoX[2]= phi;
2098 histoX[3]= pTPCIn->GetSnp();
2099 histoX[4]= pTPCIn->GetX();
2100 AliExternalTrackParam lout(*pTPCOut);
2101 lout.Rotate(pTPCIn->GetAlpha());
2102 lout.PropagateTo(pTPCIn->GetX(),fMagF);
2103 //
2104 for (Int_t ihisto=0; ihisto<5; ihisto++){
2105 histoX[0]=lout.GetParameter()[ihisto]-pTPCIn->GetParameter()[ihisto];
2106 fResHistoTPCCE[ihisto]->Fill(histoX);
2107 }
2108}
391ffdb2 2109void AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
2110 //
2111 // fill residual histograms pTPCIn-pITSOut
2112 // Histogram is filled only for primary tracks
2113 //
2114 Double_t histoX[4];
2115 Double_t xyz[3];
2116 pTPCIn->GetXYZ(xyz);
2117 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2118 histoX[1]= pTPCIn->GetTgl();
2119 histoX[2]= phi;
2120 histoX[3]= pTPCIn->GetSnp();
9112627b 2121 AliExternalTrackParam lits(*pITSOut);
2122 lits.Rotate(pTPCIn->GetAlpha());
2123 lits.PropagateTo(pTPCIn->GetX(),fMagF);
391ffdb2 2124 //
2125 for (Int_t ihisto=0; ihisto<5; ihisto++){
9112627b 2126 histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
391ffdb2 2127 fResHistoTPCITS[ihisto]->Fill(histoX);
2128 }
2129}
2130
2131
2132void AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
2133 //
2134 // fill residual histograms pTPC - vertex
2135 // Histogram is filled only for primary tracks
2136 //
660f0543 2137 if (fMemoryMode<2) return;
391ffdb2 2138 Double_t histoX[4];
2139 const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
be67055b 2140 AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
9112627b 2141 //
2142 AliExternalTrackParam lits(*pTrack);
be67055b 2143 if (TMath::Abs(pTrack->GetY())>3) return; // beam pipe
2144 pTPCvertex.Rotate(lits.GetAlpha());
2145 //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
2146 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
2147 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
391ffdb2 2148 Double_t xyz[3];
2149 pTPCIn->GetXYZ(xyz);
2150 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2151 histoX[1]= pTPCIn->GetTgl();
2152 histoX[2]= phi;
2153 histoX[3]= pTPCIn->GetSnp();
2154 //
2155 Float_t dca[2], cov[3];
2156 pTrack->GetImpactParametersTPC(dca,cov);
882f5c06 2157 for (Int_t ihisto=0; ihisto<5; ihisto++){
be67055b 2158 histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2159 // if (ihisto<2) histoX[0]=dca[ihisto];
9112627b 2160 fResHistoTPCvertex[ihisto]->Fill(histoX);
391ffdb2 2161 }
2162}
2163
2164
2165void AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
2166 //
2167 // fill resuidual histogram TPCout-TRDin
2168 //
660f0543 2169 if (fMemoryMode<2) return;
391ffdb2 2170 Double_t histoX[4];
2171 Double_t xyz[3];
2172 pTPCOut->GetXYZ(xyz);
2173 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2174 histoX[1]= pTPCOut->GetTgl();
2175 histoX[2]= phi;
2176 histoX[3]= pTPCOut->GetSnp();
2177 //
882f5c06 2178 AliExternalTrackParam ltrd(*pTRDIn);
2179 ltrd.Rotate(pTPCOut->GetAlpha());
be67055b 2180 // ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
2181 AliTracker::PropagateTrackToBxByBz(&ltrd,pTPCOut->GetX(),0.1,0.1,kFALSE);
882f5c06 2182
391ffdb2 2183 for (Int_t ihisto=0; ihisto<5; ihisto++){
882f5c06 2184 histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
391ffdb2 2185 fResHistoTPCTRD[ihisto]->Fill(histoX);
2186 }
2187
2188}
0b736a46 2189
2190void AliTPCcalibTime::FillResHistoTPCTOF(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTOFIn ){
2191 //
2192 // fill resuidual histogram TPCout-TOFin
2193 // track propagated to the TOF position
660f0543 2194 if (fMemoryMode<2) return;
0b736a46 2195 Double_t histoX[4];
2196 Double_t xyz[3];
2197
2198 AliExternalTrackParam ltpc(*pTPCOut);
2199 ltpc.Rotate(pTOFIn->GetAlpha());
2200 AliTracker::PropagateTrackToBxByBz(&ltpc,pTOFIn->GetX(),0.1,0.1,kFALSE);
2201 //
2202 ltpc.GetXYZ(xyz);
2203 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2204 histoX[1]= ltpc.GetTgl();
2205 histoX[2]= phi;
2206 histoX[3]= ltpc.GetSnp();
2207 //
2208 for (Int_t ihisto=0; ihisto<2; ihisto++){
2209 histoX[0]=ltpc.GetParameter()[ihisto]-pTOFIn->GetParameter()[ihisto];
2210 fResHistoTPCTOF[ihisto]->Fill(histoX);
2211 }
2212
2213}