1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 Comments to be written here:
19 1. What do we calibrate.
21 Time dependence of gain and drift velocity in order to account for changes in: temperature, pressure, gas composition.
23 AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime",0, 1213.9e+06, 1213.96e+06, 0.04e+04, 0.04e+04);
25 2. How to interpret results
29 a) determine the required time range:
31 AliXRDPROOFtoolkit tool;
32 TChain * chain = tool.MakeChain("pass2.txt","esdTree",0,6000);
33 chain->Draw("GetTimeStamp()")
35 b) analyse calibration object on Proof in calibration train
37 AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime", StartTimeStamp, EndTimeStamp, IntegrationTimeVdrift);
41 gSystem->Load("libANALYSIS");
42 gSystem->Load("libTPCcalib");
44 TFile f("CalibObjectsTrain1.root");
45 AliTPCcalibTime *calib = (AliTPCcalibTime *)f->Get("calibTime");
46 calib->GetHistoDrift("all")->Projection(2,0)->Draw()
47 calib->GetFitDrift("all")->Draw("lp")
49 4. Analysis using debug streamers.
51 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
52 gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
53 AliXRDPROOFtoolkit tool;
54 TChain * chainTime = tool.MakeChainRandom("time.txt","trackInfo",0,10000);
56 AliXRDPROOFtoolkit::FilterList("timetpctpc.txt","* tpctpc",1)
57 AliXRDPROOFtoolkit::FilterList("timetoftpc.txt","* toftpc",1)
58 AliXRDPROOFtoolkit::FilterList("timeitstpc.txt","* itstpc",1)
59 AliXRDPROOFtoolkit::FilterList("timelaser.txt","* laserInfo",1)
60 TChain * chainTPCTPC = tool.MakeChainRandom("timetpctpc.txt.Good","tpctpc",0,10000);
61 TChain * chainTPCITS = tool.MakeChainRandom("timeitstpc.txt.Good","itstpc",0,10000);
62 TChain * chainTPCTOF = tool.MakeChainRandom("timetoftpc.txt.Good","toftpc",0,10000);
63 TChain * chainLaser = tool.MakeChainRandom("timelaser.txt.Good","laserInfo",0,10000);
68 #include "Riostream.h"
69 #include "TDatabasePDG.h"
70 #include "TGraphErrors.h"
72 #include "THnSparse.h"
75 #include "TTimeStamp.h"
81 #include "AliDCSSensor.h"
82 #include "AliDCSSensorArray.h"
83 #include "AliESDEvent.h"
84 #include "AliESDInputHandler.h"
85 #include "AliESDVertex.h"
86 #include "AliESDfriend.h"
88 #include "AliRelAlignerKalman.h"
89 #include "AliTPCCalROC.h"
90 #include "AliTPCParam.h"
91 #include "AliTPCTracklet.h"
92 #include "AliTPCcalibDB.h"
93 #include "AliTPCcalibLaser.h"
94 #include "AliTPCcalibTime.h"
95 #include "AliTPCclusterMI.h"
96 #include "AliTPCseed.h"
97 #include "AliTrackPointArray.h"
98 #include "AliTracker.h"
100 ClassImp(AliTPCcalibTime)
103 AliTPCcalibTime::AliTPCcalibTime()
105 fLaser(0), // pointer to laser calibration
106 fDz(0), // current delta z
107 fCutMaxD(3), // maximal distance in rfi ditection
108 fCutMaxDz(25), // maximal distance in rfi ditection
109 fCutTheta(0.03), // maximal distan theta
110 fCutMinDir(-0.99), // direction vector products
112 fArrayDz(0), //NEW! Tmap of V drifts for different triggers
113 fAlignITSTPC(0), //alignemnt array ITS TPC match
114 fAlignTRDTPC(0), //alignemnt array TRD TPC match
115 fAlignTOFTPC(0), //alignemnt array TOF TPC match
130 // default constructor
132 AliInfo("Default Constructor");
133 for (Int_t i=0;i<3;i++) {
134 fHistVdriftLaserA[i]=0;
135 fHistVdriftLaserC[i]=0;
137 for (Int_t i=0;i<10;i++) {
138 fCosmiMatchingHisto[i]=0;
141 for (Int_t i=0;i<5;i++) {
142 fResHistoTPCITS[i]=0;
143 fResHistoTPCTRD[i]=0;
144 fResHistoTPCTOF[i]=0;
145 fResHistoTPCvertex[i]=0;
150 AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeVdrift)
152 fLaser(0), // pointer to laser calibration
153 fDz(0), // current delta z
154 fCutMaxD(5*0.5356), // maximal distance in rfi ditection
155 fCutMaxDz(40), // maximal distance in rfi ditection
156 fCutTheta(5*0.004644),// maximal distan theta
157 fCutMinDir(-0.99), // direction vector products
159 fArrayDz(0), //Tmap of V drifts for different triggers
160 fAlignITSTPC(0), //alignemnt array ITS TPC match
161 fAlignTRDTPC(0), //alignemnt array TRD TPC match
162 fAlignTOFTPC(0), //alignemnt array TOF TPC match
177 // Non deafaul constructor - to be used in the Calibration setups
182 for (Int_t i=0;i<3;i++) {
183 fHistVdriftLaserA[i]=0;
184 fHistVdriftLaserC[i]=0;
187 for (Int_t i=0;i<5;i++) {
188 fResHistoTPCITS[i]=0;
189 fResHistoTPCTRD[i]=0;
190 fResHistoTPCTOF[i]=0;
191 fResHistoTPCvertex[i]=0;
195 AliInfo("Non Default Constructor");
196 fTimeBins =(EndTime-StartTime)/deltaIntegrationTimeVdrift;
197 fTimeStart =StartTime; //(((TObjString*)(mapGRP->GetValue("fAliceStartTime")))->GetString()).Atoi();
198 fTimeEnd =EndTime; //(((TObjString*)(mapGRP->GetValue("fAliceStopTime")))->GetString()).Atoi();
209 Int_t binsVdriftLaser[4] = {fTimeBins , fPtBins , fVdriftBins*20, fRunBins };
210 Double_t xminVdriftLaser[4] = {fTimeStart, fPtStart, fVdriftStart , fRunStart};
211 Double_t xmaxVdriftLaser[4] = {fTimeEnd , fPtEnd , fVdriftEnd , fRunEnd };
212 TString axisTitle[4]={
218 TString histoName[3]={
225 for (Int_t i=0;i<3;i++) {
226 fHistVdriftLaserA[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
227 fHistVdriftLaserC[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
228 fHistVdriftLaserA[i]->SetName(histoName[i]);
229 fHistVdriftLaserC[i]->SetName(histoName[i]);
230 for (Int_t iaxis=0; iaxis<4;iaxis++){
231 fHistVdriftLaserA[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
232 fHistVdriftLaserC[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
235 fBinsVdrift[0] = fTimeBins;
236 fBinsVdrift[1] = fPtBins;
237 fBinsVdrift[2] = fVdriftBins;
238 fBinsVdrift[3] = fRunBins;
239 fXminVdrift[0] = fTimeStart;
240 fXminVdrift[1] = fPtStart;
241 fXminVdrift[2] = fVdriftStart;
242 fXminVdrift[3] = fRunStart;
243 fXmaxVdrift[0] = fTimeEnd;
244 fXmaxVdrift[1] = fPtEnd;
245 fXmaxVdrift[2] = fVdriftEnd;
246 fXmaxVdrift[3] = fRunEnd;
248 fArrayDz=new TObjArray();
249 fAlignITSTPC = new TObjArray; //alignemnt array ITS TPC match
250 fAlignTRDTPC = new TObjArray; //alignemnt array ITS TPC match
251 fAlignTOFTPC = new TObjArray; //alignemnt array ITS TPC match
252 fAlignITSTPC->SetOwner(kTRUE);
253 fAlignTRDTPC->SetOwner(kTRUE);
254 fAlignTOFTPC->SetOwner(kTRUE);
257 fCosmiMatchingHisto[0]=new TH1F("Cosmics matching","p0-all" ,100,-10*0.5356 ,10*0.5356 );
258 fCosmiMatchingHisto[1]=new TH1F("Cosmics matching","p1-all" ,100,-10*4.541 ,10*4.541 );
259 fCosmiMatchingHisto[2]=new TH1F("Cosmics matching","p2-all" ,100,-10*0.01134 ,10*0.01134 );
260 fCosmiMatchingHisto[3]=new TH1F("Cosmics matching","p3-all" ,100,-10*0.004644,10*0.004644);
261 fCosmiMatchingHisto[4]=new TH1F("Cosmics matching","p4-all" ,100,-10*0.03773 ,10*0.03773 );
262 fCosmiMatchingHisto[5]=new TH1F("Cosmics matching","p0-isPair",100,-10*0.5356 ,10*0.5356 );
263 fCosmiMatchingHisto[6]=new TH1F("Cosmics matching","p1-isPair",100,-10*4.541 ,10*4.541 );
264 fCosmiMatchingHisto[7]=new TH1F("Cosmics matching","p2-isPair",100,-10*0.01134 ,10*0.01134 );
265 fCosmiMatchingHisto[8]=new TH1F("Cosmics matching","p3-isPair",100,-10*0.004644,10*0.004644);
266 fCosmiMatchingHisto[9]=new TH1F("Cosmics matching","p4-isPair",100,-10*0.03773 ,10*0.03773 );
267 BookDistortionMaps();
270 AliTPCcalibTime::~AliTPCcalibTime(){
272 // Virtual Destructor
274 for(Int_t i=0;i<3;i++){
275 if(fHistVdriftLaserA[i]){
276 delete fHistVdriftLaserA[i];
277 fHistVdriftLaserA[i]=NULL;
279 if(fHistVdriftLaserC[i]){
280 delete fHistVdriftLaserC[i];
281 fHistVdriftLaserC[i]=NULL;
285 fArrayDz->SetOwner();
290 for(Int_t i=0;i<5;i++){
291 if(fCosmiMatchingHisto[i]){
292 delete fCosmiMatchingHisto[i];
293 fCosmiMatchingHisto[i]=NULL;
297 for (Int_t i=0;i<5;i++) {
298 delete fResHistoTPCITS[i];
299 delete fResHistoTPCTRD[i];
300 delete fResHistoTPCTOF[i];
301 delete fResHistoTPCvertex[i];
302 fResHistoTPCITS[i]=0;
303 fResHistoTPCTRD[i]=0;
304 fResHistoTPCTOF[i]=0;
305 fResHistoTPCvertex[i]=0;
309 fAlignITSTPC->SetOwner(kTRUE);
310 fAlignTRDTPC->SetOwner(kTRUE);
311 fAlignTOFTPC->SetOwner(kTRUE);
313 fAlignITSTPC->Delete();
314 fAlignTRDTPC->Delete();
315 fAlignTOFTPC->Delete();
321 Bool_t AliTPCcalibTime::IsLaser(const AliESDEvent *const /*event*/){
323 // Indicator is laser event not yet implemented - to be done using trigger info or event specie
325 return kTRUE; //More accurate creteria to be added
327 Bool_t AliTPCcalibTime::IsCosmics(const AliESDEvent *const /*event*/){
329 // Indicator is cosmic event not yet implemented - to be done using trigger info or event specie
332 return kTRUE; //More accurate creteria to be added
334 Bool_t AliTPCcalibTime::IsBeam(const AliESDEvent *const /*event*/){
336 // Indicator is physic event not yet implemented - to be done using trigger info or event specie
339 return kTRUE; //More accurate creteria to be added
341 void AliTPCcalibTime::ResetCurrent(){
342 fDz=0; //Reset current dz
347 void AliTPCcalibTime::Process(AliESDEvent *event){
349 // main function to make calibration
352 if (event->GetNumberOfTracks()<2) return;
354 if(IsLaser (event)) ProcessLaser (event);
355 if(IsCosmics(event)) ProcessCosmic(event);
356 if(IsBeam (event)) ProcessBeam (event);
359 void AliTPCcalibTime::ProcessLaser(AliESDEvent *event){
361 // Fit drift velocity using laser
364 const Int_t kMinTracks = 40; // minimal number of laser tracks
365 const Int_t kMinTracksSide = 20; // minimal number of tracks per side
366 const Float_t kMaxDeltaZ = 30.; // maximal trigger delay
367 const Float_t kMaxDeltaV = 0.05; // maximal deltaV
368 const Float_t kMaxRMS = 0.1; // maximal RMS of tracks
371 TCut cutRMS("sqrt(laserA.fElements[4])<0.1&&sqrt(laserC.fElements[4])<0.1");
372 TCut cutZ("abs(laserA.fElements[0]-laserC.fElements[0])<3");
373 TCut cutV("abs(laserA.fElements[1]-laserC.fElements[1])<0.01");
374 TCut cutY("abs(laserA.fElements[2]-laserC.fElements[2])<2");
375 TCut cutAll = cutRMS+cutZ+cutV+cutY;
377 if (event->GetNumberOfTracks()<kMinTracks) return;
379 if(!fLaser) fLaser = new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
380 fLaser->Process(event);
381 if (fLaser->GetNtracks()<kMinTracks) return; // small amount of tracks cut
382 if (fLaser->fFitAside->GetNrows()==0 && fLaser->fFitCside->GetNrows()==0) return; // no fit neither a or C side
384 // debug streamer - activate stream level
385 // Use it for tuning of the cuts
387 // cuts to be applied
389 Int_t isReject[2]={0,0};
392 if (TMath::Abs((*fLaser->fFitAside)[3]) < kMinTracksSide) isReject[0]|=1;
393 if (TMath::Abs((*fLaser->fFitCside)[3]) < kMinTracksSide) isReject[1]|=1;
394 // unreasonable z offset
395 if (TMath::Abs((*fLaser->fFitAside)[0])>kMaxDeltaZ) isReject[0]|=2;
396 if (TMath::Abs((*fLaser->fFitCside)[0])>kMaxDeltaZ) isReject[1]|=2;
397 // unreasonable drift velocity
398 if (TMath::Abs((*fLaser->fFitAside)[1]-1)>kMaxDeltaV) isReject[0]|=4;
399 if (TMath::Abs((*fLaser->fFitCside)[1]-1)>kMaxDeltaV) isReject[1]|=4;
401 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitAside)[4]))>kMaxRMS ) isReject[0]|=8;
402 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitCside)[4]))>kMaxRMS ) isReject[1]|=8;
408 printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
410 TTreeSRedirector *cstream = GetDebugStreamer();
412 TTimeStamp tstamp(fTime);
413 (*cstream)<<"laserInfo"<<
414 "run="<<fRun<< // run number
415 "event="<<fEvent<< // event number
416 "time="<<fTime<< // time stamp of event
417 "trigger="<<fTrigger<< // trigger
418 "mag="<<fMagF<< // magnetic field
420 "rejectA="<<isReject[0]<<
421 "rejectC="<<isReject[1]<<
422 "laserA.="<<fLaser->fFitAside<<
423 "laserC.="<<fLaser->fFitCside<<
424 "laserAC.="<<fLaser->fFitACside<<
425 "trigger="<<event->GetFiredTriggerClasses()<<
432 TVectorD vdriftA(5), vdriftC(5),vdriftAC(5);
433 vdriftA=*(fLaser->fFitAside);
434 vdriftC=*(fLaser->fFitCside);
435 vdriftAC=*(fLaser->fFitACside);
436 Int_t npointsA=0, npointsC=0;
437 Float_t chi2A=0, chi2C=0;
438 npointsA= TMath::Nint(vdriftA[3]);
440 npointsC= TMath::Nint(vdriftC[3]);
443 TTimeStamp tstamp(fTime);
444 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
445 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
446 Double_t driftA=0, driftC=0;
447 if (vdriftA[1]>1.-kMaxDeltaV) driftA = 1./vdriftA[1]-1.;
448 if (vdriftC[1]>1.-kMaxDeltaV) driftC = 1./vdriftC[1]-1.;
450 Double_t vecDriftLaserA[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftA,event->GetRunNumber()};
451 Double_t vecDriftLaserC[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftC,event->GetRunNumber()};
452 // Double_t vecDrift[4] ={fTime,(ptrelative0+ptrelative1)/2.0,1./((*(fLaser->fFitACside))[1])-1,event->GetRunNumber()};
454 for (Int_t icalib=0;icalib<3;icalib++){
455 if (icalib==0){ //z0 shift
456 vecDriftLaserA[2]=vdriftA[0]/250.;
457 vecDriftLaserC[2]=vdriftC[0]/250.;
459 if (icalib==1){ //vdrel shift
460 vecDriftLaserA[2]=driftA;
461 vecDriftLaserC[2]=driftC;
463 if (icalib==2){ //gy shift - full gy - full drift
464 vecDriftLaserA[2]=vdriftA[2]/250.;
465 vecDriftLaserC[2]=vdriftC[2]/250.;
467 //if (isReject[0]==0) fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
468 //if (isReject[1]==0) fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
469 fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
470 fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
473 // THnSparse* curHist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
474 // TString shortName=curHist->ClassName();
475 // shortName+="_MEAN_DRIFT_LASER_";
481 // name+=event->GetFiredTriggerClasses();
483 // curHist=(THnSparseF*)fArrayDz->FindObject(name);
485 // curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
486 // fArrayDz->AddLast(curHist);
488 // curHist->Fill(vecDrift);
493 // curHist=(THnSparseF*)fArrayDz->FindObject(name);
495 // curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
496 // fArrayDz->AddLast(curHist);
498 // curHist->Fill(vecDrift);
501 void AliTPCcalibTime::ProcessCosmic(const AliESDEvent *const event){
503 // process Cosmic event - track matching A side C side
506 Printf("ERROR: ESD not available");
509 if (event->GetTimeStamp() == 0 ) {
510 Printf("no time stamp!");
517 // Track0 is choosen in upper TPC part
518 // Track1 is choosen in lower TPC part
520 const Int_t kMinClustersCross =30;
521 const Int_t kMinClusters =80;
522 Int_t ntracks=event->GetNumberOfTracks();
523 if (ntracks==0) return;
524 if (ntracks > fCutTracks) return;
526 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
527 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
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);
536 TArrayI clusterSideA(ntracks);
537 TArrayI clusterSideC(ntracks);
538 for (Int_t i=0;i<ntracks;++i) {
541 AliESDtrack *track = event->GetTrack(i);
543 const AliExternalTrackParam * trackIn = track->GetInnerParam();
544 const AliExternalTrackParam * trackOut = track->GetOuterParam();
545 if (!trackIn) continue;
546 if (!trackOut) continue;
548 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
549 if (!friendTrack) continue;
550 if (friendTrack) ProcessSame(track,friendTrack,event);
551 if (friendTrack) ProcessAlignITS(track,friendTrack,event,esdFriend);
552 if (friendTrack) ProcessAlignTRD(track,friendTrack);
553 if (friendTrack) ProcessAlignTOF(track,friendTrack);
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;
558 tpcSeeds.AddAt(seed,i);
560 for (Int_t irow=159;irow>0;irow--) {
561 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
563 if ((cl->GetDetector()%36)<18) nA++;
564 if ((cl->GetDetector()%36)>=18) nC++;
570 if (ntracks<2) return;
575 for (Int_t i=0;i<ntracks;++i) {
576 AliESDtrack *track0 = event->GetTrack(i);
577 // track0 - choosen upper part
578 if (!track0) continue;
579 if (!track0->GetOuterParam()) continue;
580 if (track0->GetOuterParam()->GetAlpha()<0) continue;
582 track0->GetDirection(d1);
583 for (Int_t j=0;j<ntracks;++j) {
585 AliESDtrack *track1 = event->GetTrack(j);
587 if (!track1) continue;
588 if (!track1->GetOuterParam()) continue;
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];
597 // if (track1->GetOuterParam()->GetAlpha()>0) continue;
600 track1->GetDirection(d2);
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);
610 // conservative cuts - convergence to be guarantied
611 // applying before track propagation
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
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;
624 Float_t dmax = TMath::Max(TMath::Abs(dist0),TMath::Abs(dist1));
625 AliExternalTrackParam param0(*track0);
626 AliExternalTrackParam param1(*track1);
628 // Propagate using Magnetic field and correct fo material budget
630 AliTracker::PropagateTrackTo(¶m0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
631 AliTracker::PropagateTrackTo(¶m1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
633 // Propagate rest to the 0,0 DCA - z should be ignored
636 param0.PropagateToDCA(&vtx,bz,1000);
638 param1.PropagateToDCA(&vtx,bz,1000);
639 param0.GetDZ(0,0,0,bz,dvertex0);
640 param1.GetDZ(0,0,0,bz,dvertex1);
645 Bool_t isPair = IsPair(¶m0,¶m1);
646 Bool_t isCross = IsCross(track0, track1);
647 Bool_t isSame = IsSame(track0, track1);
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_";
655 if((isSame) || (isCross && isPair)){
656 if (track0->GetTPCNcls()+ track1->GetTPCNcls()> 80) {
657 fDz = param0.GetZ() - param1.GetZ();
658 Double_t sign=(nA0>nA1)? 1:-1;
660 TTimeStamp tstamp(fTime);
661 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
662 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
663 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
664 THnSparse* curHist=NULL;
668 name+=event->GetFiredTriggerClasses();
670 curHist=(THnSparseF*)fArrayDz->FindObject(name);
672 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
673 fArrayDz->AddLast(curHist);
675 // curHist=(THnSparseF*)(fMapDz->GetValue(event->GetFiredTriggerClasses()));
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);
680 curHist->Fill(vecDrift);
685 curHist=(THnSparseF*)fArrayDz->FindObject(name);
687 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
688 fArrayDz->AddLast(curHist);
690 // curHist=(THnSparseF*)(fMapDz->GetValue("all"));
692 // curHist=new THnSparseF("all","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
693 // fMapDz->Add(new TObjString("all"),curHist);
695 curHist->Fill(vecDrift);
698 TTreeSRedirector *cstream = GetDebugStreamer();
701 (*cstream)<<"trackInfo"<<
712 "isCross="<<isCross<<
720 } // end 2nd order loop
721 } // end 1st order loop
724 TTreeSRedirector *cstream = GetDebugStreamer();
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
734 // accumulated values
736 "fDz="<<fDz<< //! current delta z
737 "trigger="<<event->GetFiredTriggerClasses()<<
741 if (GetDebugLevel()>20) printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
744 void AliTPCcalibTime::ProcessBeam(const AliESDEvent *const /*event*/){
746 // Not special treatment yet - the same for cosmic and physic event
750 void AliTPCcalibTime::Analyze(){
752 // Special macro to analyze result of calibration and extract calibration entries
753 // Not yet ported to the Analyze function yet
757 THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name) const
760 // Get histogram for given trigger mask
762 TIterator* iterator = fArrayDz->MakeIterator();
764 TString newName=name;
766 THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
767 THnSparse* addHist=NULL;
768 while((addHist=(THnSparseF*)iterator->Next())){
769 if(!addHist) continue;
770 TString histName=addHist->GetName();
771 if(!histName.Contains(newName)) continue;
773 newHist->Add(addHist);
778 TObjArray* AliTPCcalibTime::GetHistoDrift() const
781 // return array of histograms
786 TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
788 // Make a drift velocity (delta Z) graph
790 THnSparse* histoDrift=GetHistoDrift(name);
791 TGraphErrors* graphDrift=NULL;
793 graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
794 TString end=histoDrift->GetName();
795 Int_t pos=end.Index("_");
796 end=end(pos,end.Capacity()-pos);
797 TString graphName=graphDrift->ClassName();
800 graphDrift->SetName(graphName);
805 TObjArray* AliTPCcalibTime::GetGraphDrift(){
807 // make a array of drift graphs
809 TObjArray* arrayGraphDrift=new TObjArray();
810 TIterator* iterator=fArrayDz->MakeIterator();
812 THnSparse* addHist=NULL;
813 while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
814 return arrayGraphDrift;
817 AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
819 // Make a fit AliSplinefit of drift velocity
821 TGraph* graphDrift=GetGraphDrift(name);
822 AliSplineFit* fitDrift=NULL;
823 if(graphDrift && graphDrift->GetN()){
824 fitDrift=new AliSplineFit();
825 fitDrift->SetGraph(graphDrift);
826 fitDrift->SetMinPoints(graphDrift->GetN()+1);
827 fitDrift->InitKnots(graphDrift,2,0,0.001);
828 fitDrift->SplineFit(0);
829 TString end=graphDrift->GetName();
830 Int_t pos=end.Index("_");
831 end=end(pos,end.Capacity()-pos);
832 TString fitName=fitDrift->ClassName();
835 //fitDrift->SetName(fitName);
843 Long64_t AliTPCcalibTime::Merge(TCollection *const li) {
845 // Object specific merging procedure
847 TIterator* iter = li->MakeIterator();
848 AliTPCcalibTime* cal = 0;
850 while ((cal = (AliTPCcalibTime*)iter->Next())) {
851 if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
852 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
855 for (Int_t imeas=0; imeas<3; imeas++){
856 if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
857 fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
858 fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
862 for (Int_t imeas=0; imeas<5; imeas++){
863 if (cal->GetResHistoTPCITS(imeas) && cal->GetResHistoTPCITS(imeas)){
864 fResHistoTPCITS[imeas]->Add(cal->fResHistoTPCITS[imeas]);
865 fResHistoTPCvertex[imeas]->Add(cal->fResHistoTPCvertex[imeas]);
867 if (cal->fResHistoTPCTRD[imeas]){
868 if (fResHistoTPCTRD[imeas])
869 fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
871 fResHistoTPCTRD[imeas]=(THnSparse*)cal->fResHistoTPCTRD[imeas]->Clone();
873 if (cal->fResHistoTPCTOF[imeas]){
874 if (fResHistoTPCTOF[imeas])
875 fResHistoTPCTOF[imeas]->Add(cal->fResHistoTPCTOF[imeas]);
877 fResHistoTPCTOF[imeas]=(THnSparse*)cal->fResHistoTPCTOF[imeas]->Clone();
880 TObjArray* addArray=cal->GetHistoDrift();
881 if(!addArray) return 0;
882 TIterator* iterator = addArray->MakeIterator();
884 THnSparse* addHist=NULL;
885 while((addHist=(THnSparseF*)iterator->Next())){
886 if(!addHist) continue;
888 THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
890 localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
891 fArrayDz->AddLast(localHist);
893 localHist->Add(addHist);
896 for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
900 for (Int_t itype=0; itype<3; itype++){
905 if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
906 if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
907 if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
909 if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
910 if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
911 arr0->Expand(arr1->GetEntriesFast());
913 for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
914 AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
915 AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
916 if (!kalman1) continue;
917 if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
918 kalman0->SetRejectOutliers(kFALSE);
919 kalman0->Merge(kalman1);
927 Bool_t AliTPCcalibTime::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
929 // 0. Same direction - OPOSITE - cutDir +cutT
930 TCut cutDir("cutDir","dir<-0.99")
932 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
935 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
937 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
940 const Double_t *p0 = tr0->GetParameter();
941 const Double_t *p1 = tr1->GetParameter();
942 fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
943 fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
944 fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
945 fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
946 fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
948 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
949 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
950 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
951 Double_t d0[3], d1[3];
952 tr0->GetDirection(d0);
953 tr1->GetDirection(d1);
954 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
956 fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
957 fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
958 fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
959 fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
960 fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
964 Bool_t AliTPCcalibTime::IsCross(AliESDtrack *const tr0, AliESDtrack *const tr1){
966 // check if the cosmic pair of tracks crossed A/C side
968 Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
969 if (result==kFALSE) return result;
974 Bool_t AliTPCcalibTime::IsSame(AliESDtrack *const tr0, AliESDtrack *const tr1){
976 // track crossing the CE
977 // 0. minimal number of clusters
978 // 1. Same sector +-1
979 // 2. Inner and outer track param on opposite side
980 // 3. Outer and inner track parameter close each to other
984 // inner and outer on opposite sides in z
986 const Int_t knclCut0 = 30;
987 const Double_t kalphaCut = 0.4;
989 // 0. minimal number of clusters
991 if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
992 if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
994 // 1. alpha cut - sector+-1
996 if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
1000 if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
1001 if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
1002 if (result==kFALSE){
1007 const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
1008 const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
1009 const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1010 const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1012 if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD) result&=kFALSE;
1013 if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1014 if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1015 if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1016 if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD) result&=kFALSE;
1017 if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1018 if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1019 if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1021 result=kTRUE; // just to put break point here
1027 void AliTPCcalibTime::ProcessSame(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
1029 // Process TPC tracks crossing CE
1031 // 0. Select only track crossing the CE
1032 // 1. Cut on the track length
1033 // 2. Refit the terack on A and C side separatelly
1034 // 3. Fill time histograms
1035 const Int_t kMinNcl=100;
1036 const Int_t kMinNclS=25; // minimul number of clusters on the sides
1037 if (!friendTrack->GetTPCOut()) return;
1039 // 0. Select only track crossing the CE
1041 if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1043 // 1. cut on track length
1045 if (track->GetTPCNcls()<kMinNcl) return;
1047 // 2. Refit track sepparatel on A and C side
1049 TObject *calibObject;
1050 AliTPCseed *seed = 0;
1051 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1052 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1056 AliExternalTrackParam trackIn(*track->GetInnerParam());
1057 AliExternalTrackParam trackOut(*track->GetOuterParam());
1058 Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1059 Double_t xyz[3]={0,0.,0.0};
1061 Int_t nclIn=0,nclOut=0;
1062 trackIn.ResetCovariance(30.);
1063 trackOut.ResetCovariance(30.);
1067 for (Int_t irow=0;irow<159;irow++) {
1068 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1070 if (cl->GetX()<80) continue;
1071 if (track->GetInnerParam()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1072 if (track->GetInnerParam()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1073 Int_t sector = cl->GetDetector();
1074 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1075 if (TMath::Abs(dalpha)>0.01){
1076 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1078 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1079 trackIn.GetXYZ(xyz);
1080 bz = AliTracker::GetBz(xyz);
1081 if (!trackIn.PropagateTo(r[0],bz)) break;
1083 trackIn.Update(&r[1],cov);
1088 for (Int_t irow=159;irow>0;irow--) {
1089 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1091 if (cl->GetX()<80) continue;
1092 if (cl->GetZ()*track->GetOuterParam()->GetZ()<0) break;
1093 if (friendTrack->GetTPCOut()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1094 if (friendTrack->GetTPCOut()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1095 Int_t sector = cl->GetDetector();
1096 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1097 if (TMath::Abs(dalpha)>0.01){
1098 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1100 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1101 trackOut.GetXYZ(xyz);
1102 bz = AliTracker::GetBz(xyz);
1103 if (!trackOut.PropagateTo(r[0],bz)) break;
1105 trackOut.Update(&r[1],cov);
1107 trackOut.Rotate(trackIn.GetAlpha());
1108 Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1109 trackIn.PropagateTo(meanX,bz);
1110 trackOut.PropagateTo(meanX,bz);
1111 TTreeSRedirector *cstream = GetDebugStreamer();
1114 trackIn.GetXYZ(gxyz.GetMatrixArray());
1115 TTimeStamp tstamp(fTime);
1116 (*cstream)<<"tpctpc"<<
1117 "run="<<fRun<< // run number
1118 "event="<<fEvent<< // event number
1119 "time="<<fTime<< // time stamp of event
1120 "trigger="<<fTrigger<< // trigger
1121 "mag="<<fMagF<< // magnetic field
1123 "xyz.="<<&gxyz<< // global position
1124 "tIn.="<<&trackIn<< // refitterd track in
1125 "tOut.="<<&trackOut<< // refitter track out
1126 "nclIn="<<nclIn<< //
1127 "nclOut="<<nclOut<< //
1131 // 3. Fill time histograms
1132 // Debug stremaer expression
1133 // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1134 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1135 fDz = trackOut.GetZ()-trackIn.GetZ();
1136 if (trackOut.GetTgl()<0) fDz*=-1.;
1137 TTimeStamp tstamp(fTime);
1138 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1139 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1140 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1142 // fill histograms per trigger class and itegrated
1144 THnSparse* curHist=NULL;
1145 for (Int_t itype=0; itype<2; itype++){
1146 TString name="MEAN_VDRIFT_CROSS_";
1148 name+=event->GetFiredTriggerClasses();
1153 curHist=(THnSparseF*)fArrayDz->FindObject(name);
1155 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1156 fArrayDz->AddLast(curHist);
1158 curHist->Fill(vecDrift);
1164 void AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
1166 // Process track - Update TPC-ITS alignment
1168 // 0. Apply standartd cuts
1169 // 1. Recalucluate the current statistic median/RMS
1170 // 2. Apply median+-rms cut
1171 // 3. Update kalman filter
1173 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1174 const Int_t kMinITS = 3; // minimal number of ITS cluster
1175 const Double_t kMinZ = 10; // maximal dz distance
1176 const Double_t kMaxDy = 2.; // maximal dy distance
1177 const Double_t kMaxAngle= 0.015; // maximal angular distance
1178 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1179 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1180 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1181 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1182 const Double_t kMinPt = 0.3; // minimal pt
1183 const Int_t kN=50; // deepnes of history
1184 static Int_t kglast=0;
1185 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1188 TCut cut="abs(pTPC.fP[2]-pITS.fP[2])<0.01&&abs(pTPC.fP[3]-pITS.fP[3])<0.01&&abs(pTPC.fP[2]-pITS.fP[2])<1";
1191 // 0. Apply standard cuts
1193 Int_t dummycl[1000];
1194 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1195 if (track->GetITSclusters(dummycl)<kMinITS) return; // minimal amount of clusters
1196 if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
1197 if (!friendTrack->GetITSOut()) return;
1198 if (!track->GetInnerParam()) return;
1199 if (!track->GetOuterParam()) return;
1200 if (track->GetInnerParam()->Pt()<kMinPt) return;
1201 // exclude crossing track
1202 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1203 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1204 if (track->GetInnerParam()->GetX()>90) return;
1206 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
1207 AliExternalTrackParam pITS(*(friendTrack->GetITSOut())); // ITS standalone if possible
1208 AliExternalTrackParam pITS2(*(friendTrack->GetITSOut())); //TPC-ITS track
1209 pITS2.Rotate(pTPC.GetAlpha());
1210 // pITS2.PropagateTo(pTPC.GetX(),fMagF);
1211 AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1213 AliESDfriendTrack *itsfriendTrack=0;
1215 // try to find standalone ITS track corresponing to the TPC if possible
1217 Bool_t hasAlone=kFALSE;
1218 Int_t ntracks=event->GetNumberOfTracks();
1219 for (Int_t i=0; i<ntracks; i++){
1220 AliESDtrack *trackS = event->GetTrack(i);
1221 if (trackS->GetTPCNcls()>0) continue; //continue if has TPC info
1222 itsfriendTrack = esdFriend->GetTrack(i);
1223 if (!itsfriendTrack) continue;
1224 if (!itsfriendTrack->GetITSOut()) continue;
1225 if (TMath::Abs(pITS2.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1226 pITS=(*(itsfriendTrack->GetITSOut()));
1228 pITS.Rotate(pTPC.GetAlpha());
1229 //pITS.PropagateTo(pTPC.GetX(),fMagF);
1230 AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
1231 if (TMath::Abs(pITS2.GetY()-pITS.GetY())> kMaxDy) continue;
1234 if (!hasAlone) pITS=pITS2;
1236 if (TMath::Abs(pITS.GetY()-pTPC.GetY()) >kMaxDy) return;
1237 if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1238 if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1240 // 1. Update median and RMS info
1242 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1243 TVectorD vecDeltaN(5);
1244 Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1246 for (Int_t i=0;i<4;i++){
1247 vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1248 kgdP[i][kglast%kN]=vecDelta[i];
1251 Int_t entries=(kglast<kN)?kglast:kN;
1252 for (Int_t i=0;i<4;i++){
1253 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1254 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1257 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1258 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1262 // 2. Apply median+-rms cut
1264 if (kglast<3) return; //median and RMS to be defined
1265 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1267 // 3. Update alignment
1269 Int_t htime = fTime/3600; //time in hours
1270 if (fAlignITSTPC->GetEntriesFast()<htime){
1271 fAlignITSTPC->Expand(htime*2+20);
1273 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1275 // make Alignment object if doesn't exist
1276 align=new AliRelAlignerKalman();
1277 align->SetRunNumber(fRun);
1278 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1279 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1280 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1281 align->SetRejectOutliers(kFALSE);
1283 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1284 align->SetMagField(fMagF);
1285 fAlignITSTPC->AddAt(align,htime);
1287 align->AddTrackParams(&pITS,&pTPC);
1288 align->SetTimeStamp(fTime);
1289 align->SetRunNumber(fRun );
1290 Float_t dca[2],cov[3];
1291 track->GetImpactParameters(dca,cov);
1292 if (TMath::Abs(dca[0])<kMaxDy){
1293 FillResHistoTPCITS(&pTPC,&pITS);
1294 FillResHistoTPC(track);
1297 Int_t nupdates=align->GetNUpdates();
1298 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1299 align->SetRejectOutliers(kFALSE);
1300 TTreeSRedirector *cstream = GetDebugStreamer();
1301 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1302 TVectorD gpTPC(3), gdTPC(3);
1303 TVectorD gpITS(3), gdITS(3);
1304 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1305 pTPC.GetDirection(gdTPC.GetMatrixArray());
1306 pITS.GetXYZ(gpITS.GetMatrixArray());
1307 pITS.GetDirection(gdITS.GetMatrixArray());
1308 (*cstream)<<"itstpc"<<
1309 "run="<<fRun<< // run number
1310 "event="<<fEvent<< // event number
1311 "time="<<fTime<< // time stamp of event
1312 "trigger="<<fTrigger<< // trigger
1313 "mag="<<fMagF<< // magnetic field
1315 "hasAlone="<<hasAlone<< // has ITS standalone ?
1316 "track.="<<track<< // track info
1317 "nmed="<<kglast<< // number of entries to define median and RMS
1318 "vMed.="<<&vecMedian<< // median of deltas
1319 "vRMS.="<<&vecRMS<< // rms of deltas
1320 "vDelta.="<<&vecDelta<< // delta in respect to median
1321 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1322 "t.="<<track<< // ful track - find proper cuts
1323 "a.="<<align<< // current alignment
1324 "pITS.="<<&pITS<< // track param ITS
1325 "pITS2.="<<&pITS2<< // track param ITS+TPC
1326 "pTPC.="<<&pTPC<< // track param TPC
1327 "gpTPC.="<<&gpTPC<< // global position TPC
1328 "gdTPC.="<<&gdTPC<< // global direction TPC
1329 "gpITS.="<<&gpITS<< // global position ITS
1330 "gdITS.="<<&gdITS<< // global position ITS
1338 void AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1340 // Process track - Update TPC-TRD alignment
1342 // 0. Apply standartd cuts
1343 // 1. Recalucluate the current statistic median/RMS
1344 // 2. Apply median+-rms cut
1345 // 3. Update kalman filter
1347 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1348 const Int_t kMinTRD = 50; // minimal number of TRD cluster
1349 const Double_t kMinZ = 20; // maximal dz distance
1350 const Double_t kMaxDy = 5.; // maximal dy distance
1351 const Double_t kMaxAngle= 0.1; // maximal angular distance
1352 const Double_t kSigmaCut= 10; // maximal sigma distance to median
1353 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1354 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1355 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1356 const Double_t kRefX = 275; // reference X
1357 const Int_t kN=50; // deepnes of history
1358 static Int_t kglast=0;
1359 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1361 // 0. Apply standard cuts
1363 Int_t dummycl[1000];
1364 if (track->GetTRDclusters(dummycl)<kMinTRD) return; // minimal amount of clusters
1365 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1366 if (!friendTrack->GetTRDIn()) return;
1367 if (!track->IsOn(AliESDtrack::kTRDrefit)) return;
1368 if (!track->IsOn(AliESDtrack::kTRDout)) return;
1369 if (!track->GetInnerParam()) return;
1370 if (!friendTrack->GetTPCOut()) return;
1371 // exclude crossing track
1372 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1373 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1375 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1376 AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
1377 AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
1378 pTRD.Rotate(pTPC.GetAlpha());
1379 // pTRD.PropagateTo(pTPC.GetX(),fMagF);
1380 AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1382 ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.; // increas sys errors
1383 ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1385 if (TMath::Abs(pTRD.GetY()-pTPC.GetY()) >kMaxDy) return;
1386 if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1387 // if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1389 // 1. Update median and RMS info
1391 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1392 TVectorD vecDeltaN(5);
1393 Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1395 for (Int_t i=0;i<4;i++){
1396 vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1397 kgdP[i][kglast%kN]=vecDelta[i];
1400 Int_t entries=(kglast<kN)?kglast:kN;
1401 for (Int_t i=0;i<4;i++){
1402 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1404 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1407 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1408 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1412 // 2. Apply median+-rms cut
1414 if (kglast<3) return; //median and RMS to be defined
1415 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1417 // 3. Update alignment
1419 Int_t htime = fTime/3600; //time in hours
1420 if (fAlignTRDTPC->GetEntriesFast()<htime){
1421 fAlignTRDTPC->Expand(htime*2+20);
1423 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1425 // make Alignment object if doesn't exist
1426 align=new AliRelAlignerKalman();
1427 align->SetRunNumber(fRun);
1428 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1429 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1430 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1431 align->SetRejectOutliers(kFALSE);
1432 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1433 align->SetMagField(fMagF);
1434 fAlignTRDTPC->AddAt(align,htime);
1436 align->AddTrackParams(&pTRD,&pTPC);
1437 align->SetTimeStamp(fTime);
1438 align->SetRunNumber(fRun );
1439 Float_t dca[2],cov[3];
1440 track->GetImpactParameters(dca,cov);
1441 if (TMath::Abs(dca[0])<kMaxDy){
1442 FillResHistoTPCTRD(&pTPC,&pTRD); //only primaries
1445 Int_t nupdates=align->GetNUpdates();
1446 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1447 align->SetRejectOutliers(kFALSE);
1448 TTreeSRedirector *cstream = GetDebugStreamer();
1449 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1450 TVectorD gpTPC(3), gdTPC(3);
1451 TVectorD gpTRD(3), gdTRD(3);
1452 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1453 pTPC.GetDirection(gdTPC.GetMatrixArray());
1454 pTRD.GetXYZ(gpTRD.GetMatrixArray());
1455 pTRD.GetDirection(gdTRD.GetMatrixArray());
1456 (*cstream)<<"trdtpc"<<
1457 "run="<<fRun<< // run number
1458 "event="<<fEvent<< // event number
1459 "time="<<fTime<< // time stamp of event
1460 "trigger="<<fTrigger<< // trigger
1461 "mag="<<fMagF<< // magnetic field
1463 "nmed="<<kglast<< // number of entries to define median and RMS
1464 "vMed.="<<&vecMedian<< // median of deltas
1465 "vRMS.="<<&vecRMS<< // rms of deltas
1466 "vDelta.="<<&vecDelta<< // delta in respect to median
1467 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1468 "t.="<<track<< // ful track - find proper cuts
1469 "a.="<<align<< // current alignment
1470 "pTRD.="<<&pTRD<< // track param TRD
1471 "pTPC.="<<&pTPC<< // track param TPC
1472 "gpTPC.="<<&gpTPC<< // global position TPC
1473 "gdTPC.="<<&gdTPC<< // global direction TPC
1474 "gpTRD.="<<&gpTRD<< // global position TRD
1475 "gdTRD.="<<&gdTRD<< // global position TRD
1481 void AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1484 // Process track - Update TPC-TOF alignment
1486 // -1. Make a TOF "track"
1487 // 0. Apply standartd cuts
1488 // 1. Recalucluate the current statistic median/RMS
1489 // 2. Apply median+-rms cut
1490 // 3. Update kalman filter
1492 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1493 // const Double_t kMinZ = 10; // maximal dz distance
1494 const Double_t kMaxDy = 5.; // maximal dy distance
1495 const Double_t kMaxAngle= 0.015; // maximal angular distance
1496 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1497 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1498 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1500 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1501 const Int_t kN=50; // deepnes of history
1502 static Int_t kglast=0;
1503 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1505 // -1. Make a TOF track-
1506 // Clusters are not in friends - use alingment points
1508 if (track->GetTOFsignal()<=0) return;
1509 if (!friendTrack->GetTPCOut()) return;
1510 if (!track->GetInnerParam()) return;
1511 if (!friendTrack->GetTPCOut()) return;
1512 const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1513 if (!points) return;
1514 AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
1515 AliExternalTrackParam pTOF(pTPC);
1516 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1517 Int_t npoints = points->GetNPoints();
1518 AliTrackPoint point;
1521 for (Int_t ipoint=0;ipoint<npoints;ipoint++){
1522 points->GetPoint(point,ipoint);
1525 Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1526 if (r<350) continue;
1527 if (r>400) continue;
1528 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1529 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);
1530 AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1531 pTPC.PropagateTo(lpoint.GetX(),fMagF);
1533 ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1534 ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
1535 ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1536 ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
1537 ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1538 ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
1541 if (naccept==0) return; // no tof match clusters
1543 // 0. Apply standard cuts
1545 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1546 // exclude crossing track
1547 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1549 if (TMath::Abs(pTOF.GetY()-pTPC.GetY()) >kMaxDy) return;
1550 if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1551 if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1553 // 1. Update median and RMS info
1555 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1556 TVectorD vecDeltaN(5);
1557 Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1559 for (Int_t i=0;i<4;i++){
1560 vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1561 kgdP[i][kglast%kN]=vecDelta[i];
1564 Int_t entries=(kglast<kN)?kglast:kN;
1566 for (Int_t i=0;i<4;i++){
1567 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1568 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1571 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1572 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1573 if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
1577 // 2. Apply median+-rms cut
1579 if (kglast<10) return; //median and RMS to be defined
1582 // 3. Update alignment
1584 Int_t htime = fTime/3600; //time in hours
1585 if (fAlignTOFTPC->GetEntriesFast()<htime){
1586 fAlignTOFTPC->Expand(htime*2+20);
1588 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1590 // make Alignment object if doesn't exist
1591 align=new AliRelAlignerKalman();
1592 align->SetRunNumber(fRun);
1593 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1594 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1595 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1596 align->SetRejectOutliers(kFALSE);
1597 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1598 align->SetMagField(fMagF);
1599 fAlignTOFTPC->AddAt(align,htime);
1601 align->AddTrackParams(&pTOF,&pTPC);
1602 Float_t dca[2],cov[3];
1603 track->GetImpactParameters(dca,cov);
1604 if (TMath::Abs(dca[0])<kMaxDy){
1605 FillResHistoTPCTOF(&pTPC,&pTOF);
1607 align->SetTimeStamp(fTime);
1608 align->SetRunNumber(fRun );
1610 Int_t nupdates=align->GetNUpdates();
1611 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1612 align->SetRejectOutliers(kFALSE);
1613 TTreeSRedirector *cstream = GetDebugStreamer();
1614 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1615 TVectorD gpTPC(3), gdTPC(3);
1616 TVectorD gpTOF(3), gdTOF(3);
1617 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1618 pTPC.GetDirection(gdTPC.GetMatrixArray());
1619 pTOF.GetXYZ(gpTOF.GetMatrixArray());
1620 pTOF.GetDirection(gdTOF.GetMatrixArray());
1621 (*cstream)<<"toftpc"<<
1622 "run="<<fRun<< // run number
1623 "event="<<fEvent<< // event number
1624 "time="<<fTime<< // time stamp of event
1625 "trigger="<<fTrigger<< // trigger
1626 "mag="<<fMagF<< // magnetic field
1628 "nmed="<<kglast<< // number of entries to define median and RMS
1629 "vMed.="<<&vecMedian<< // median of deltas
1630 "vRMS.="<<&vecRMS<< // rms of deltas
1631 "vDelta.="<<&vecDelta<< // delta in respect to median
1632 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1633 "t.="<<track<< // ful track - find proper cuts
1634 "a.="<<align<< // current alignment
1635 "pTOF.="<<&pTOF<< // track param TOF
1636 "pTPC.="<<&pTPC<< // track param TPC
1637 "gpTPC.="<<&gpTPC<< // global position TPC
1638 "gdTPC.="<<&gdTPC<< // global direction TPC
1639 "gpTOF.="<<&gpTOF<< // global position TOF
1640 "gdTOF.="<<&gdTOF<< // global position TOF
1646 void AliTPCcalibTime::BookDistortionMaps(){
1648 // Book ndimensional histograms of distortions/residuals
1649 // Only primary tracks are selected for analysis
1652 Double_t xminTrack[4], xmaxTrack[4];
1654 TString axisName[4];
1655 TString axisTitle[4];
1658 axisName[0] ="#Delta";
1659 axisTitle[0] ="#Delta";
1662 xminTrack[1] =-1.5; xmaxTrack[1]=1.5;
1663 axisName[1] ="tanTheta";
1664 axisTitle[1] ="tan(#Theta)";
1667 xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi();
1669 axisTitle[2] ="#phi";
1672 xminTrack[3] =-1.; xmaxTrack[3]=1.; // 0.33 GeV cut
1676 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
1677 fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1678 fResHistoTPCvertex[0] = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1679 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
1680 fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1681 xminTrack[0] =-5; xmaxTrack[0]=5; //
1682 fResHistoTPCTOF[0] = new THnSparseS("TPCTOF#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1685 xminTrack[0] =-3.; xmaxTrack[0]=3.; //
1686 fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1687 fResHistoTPCvertex[1] = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1688 fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1689 xminTrack[0] =-5.; xmaxTrack[0]=5.; //
1690 fResHistoTPCTOF[1] = new THnSparseS("TPCTOF#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1693 xminTrack[0] =-0.015; xmaxTrack[0]=0.015; //
1694 fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1695 fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1696 fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1697 fResHistoTPCTOF[2] = new THnSparseS("TPCTOF#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1700 xminTrack[0] =-0.025; xmaxTrack[0]=0.025; //
1701 fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1702 fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1703 fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1704 fResHistoTPCTOF[3] = new THnSparseS("TPCTOF#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1707 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
1708 fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1709 fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1710 fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1711 fResHistoTPCTOF[4] = new THnSparseS("TPCTOF#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1713 for (Int_t ivar=0;ivar<4;ivar++){
1714 for (Int_t ivar2=0;ivar2<4;ivar2++){
1715 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1716 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1717 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1718 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1719 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1720 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1726 void AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
1728 // fill residual histograms pTPCIn-pITSOut
1729 // Histogram is filled only for primary tracks
1733 pTPCIn->GetXYZ(xyz);
1734 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1735 histoX[1]= pTPCIn->GetTgl();
1737 histoX[3]= pTPCIn->GetSnp();
1738 AliExternalTrackParam lits(*pITSOut);
1739 lits.Rotate(pTPCIn->GetAlpha());
1740 lits.PropagateTo(pTPCIn->GetX(),fMagF);
1742 for (Int_t ihisto=0; ihisto<5; ihisto++){
1743 histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
1744 fResHistoTPCITS[ihisto]->Fill(histoX);
1749 void AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
1751 // fill residual histograms pTPC - vertex
1752 // Histogram is filled only for primary tracks
1755 const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
1756 AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
1758 AliExternalTrackParam lits(*pTrack);
1759 if (TMath::Abs(pTrack->GetY())>3) return; // beam pipe
1760 pTPCvertex.Rotate(lits.GetAlpha());
1761 //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
1762 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
1763 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
1765 pTPCIn->GetXYZ(xyz);
1766 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1767 histoX[1]= pTPCIn->GetTgl();
1769 histoX[3]= pTPCIn->GetSnp();
1771 Float_t dca[2], cov[3];
1772 pTrack->GetImpactParametersTPC(dca,cov);
1773 for (Int_t ihisto=0; ihisto<5; ihisto++){
1774 histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
1775 // if (ihisto<2) histoX[0]=dca[ihisto];
1776 fResHistoTPCvertex[ihisto]->Fill(histoX);
1781 void AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
1783 // fill resuidual histogram TPCout-TRDin
1787 pTPCOut->GetXYZ(xyz);
1788 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1789 histoX[1]= pTPCOut->GetTgl();
1791 histoX[3]= pTPCOut->GetSnp();
1793 AliExternalTrackParam ltrd(*pTRDIn);
1794 ltrd.Rotate(pTPCOut->GetAlpha());
1795 // ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
1796 AliTracker::PropagateTrackToBxByBz(<rd,pTPCOut->GetX(),0.1,0.1,kFALSE);
1798 for (Int_t ihisto=0; ihisto<5; ihisto++){
1799 histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
1800 fResHistoTPCTRD[ihisto]->Fill(histoX);
1805 void AliTPCcalibTime::FillResHistoTPCTOF(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTOFIn ){
1807 // fill resuidual histogram TPCout-TOFin
1808 // track propagated to the TOF position
1812 AliExternalTrackParam ltpc(*pTPCOut);
1813 ltpc.Rotate(pTOFIn->GetAlpha());
1814 AliTracker::PropagateTrackToBxByBz(<pc,pTOFIn->GetX(),0.1,0.1,kFALSE);
1817 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1818 histoX[1]= ltpc.GetTgl();
1820 histoX[3]= ltpc.GetSnp();
1822 for (Int_t ihisto=0; ihisto<2; ihisto++){
1823 histoX[0]=ltpc.GetParameter()[ihisto]-pTOFIn->GetParameter()[ihisto];
1824 fResHistoTPCTOF[ihisto]->Fill(histoX);