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 fResHistoTPCvertex[i]=0;
149 AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeVdrift)
151 fLaser(0), // pointer to laser calibration
152 fDz(0), // current delta z
153 fCutMaxD(5*0.5356), // maximal distance in rfi ditection
154 fCutMaxDz(40), // maximal distance in rfi ditection
155 fCutTheta(5*0.004644),// maximal distan theta
156 fCutMinDir(-0.99), // direction vector products
158 fArrayDz(0), //Tmap of V drifts for different triggers
159 fAlignITSTPC(0), //alignemnt array ITS TPC match
160 fAlignTRDTPC(0), //alignemnt array TRD TPC match
161 fAlignTOFTPC(0), //alignemnt array TOF TPC match
176 // Non deafaul constructor - to be used in the Calibration setups
181 for (Int_t i=0;i<3;i++) {
182 fHistVdriftLaserA[i]=0;
183 fHistVdriftLaserC[i]=0;
186 for (Int_t i=0;i<5;i++) {
187 fResHistoTPCITS[i]=0;
188 fResHistoTPCTRD[i]=0;
189 fResHistoTPCvertex[i]=0;
193 AliInfo("Non Default Constructor");
194 fTimeBins =(EndTime-StartTime)/deltaIntegrationTimeVdrift;
195 fTimeStart =StartTime; //(((TObjString*)(mapGRP->GetValue("fAliceStartTime")))->GetString()).Atoi();
196 fTimeEnd =EndTime; //(((TObjString*)(mapGRP->GetValue("fAliceStopTime")))->GetString()).Atoi();
207 Int_t binsVdriftLaser[4] = {fTimeBins , fPtBins , fVdriftBins*20, fRunBins };
208 Double_t xminVdriftLaser[4] = {fTimeStart, fPtStart, fVdriftStart , fRunStart};
209 Double_t xmaxVdriftLaser[4] = {fTimeEnd , fPtEnd , fVdriftEnd , fRunEnd };
210 TString axisTitle[4]={
216 TString histoName[3]={
223 for (Int_t i=0;i<3;i++) {
224 fHistVdriftLaserA[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
225 fHistVdriftLaserC[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
226 fHistVdriftLaserA[i]->SetName(histoName[i]);
227 fHistVdriftLaserC[i]->SetName(histoName[i]);
228 for (Int_t iaxis=0; iaxis<4;iaxis++){
229 fHistVdriftLaserA[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
230 fHistVdriftLaserC[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
233 fBinsVdrift[0] = fTimeBins;
234 fBinsVdrift[1] = fPtBins;
235 fBinsVdrift[2] = fVdriftBins;
236 fBinsVdrift[3] = fRunBins;
237 fXminVdrift[0] = fTimeStart;
238 fXminVdrift[1] = fPtStart;
239 fXminVdrift[2] = fVdriftStart;
240 fXminVdrift[3] = fRunStart;
241 fXmaxVdrift[0] = fTimeEnd;
242 fXmaxVdrift[1] = fPtEnd;
243 fXmaxVdrift[2] = fVdriftEnd;
244 fXmaxVdrift[3] = fRunEnd;
246 fArrayDz=new TObjArray();
247 fAlignITSTPC = new TObjArray; //alignemnt array ITS TPC match
248 fAlignTRDTPC = new TObjArray; //alignemnt array ITS TPC match
249 fAlignTOFTPC = new TObjArray; //alignemnt array ITS TPC match
250 fAlignITSTPC->SetOwner(kTRUE);
251 fAlignTRDTPC->SetOwner(kTRUE);
252 fAlignTOFTPC->SetOwner(kTRUE);
255 fCosmiMatchingHisto[0]=new TH1F("Cosmics matching","p0-all" ,100,-10*0.5356 ,10*0.5356 );
256 fCosmiMatchingHisto[1]=new TH1F("Cosmics matching","p1-all" ,100,-10*4.541 ,10*4.541 );
257 fCosmiMatchingHisto[2]=new TH1F("Cosmics matching","p2-all" ,100,-10*0.01134 ,10*0.01134 );
258 fCosmiMatchingHisto[3]=new TH1F("Cosmics matching","p3-all" ,100,-10*0.004644,10*0.004644);
259 fCosmiMatchingHisto[4]=new TH1F("Cosmics matching","p4-all" ,100,-10*0.03773 ,10*0.03773 );
260 fCosmiMatchingHisto[5]=new TH1F("Cosmics matching","p0-isPair",100,-10*0.5356 ,10*0.5356 );
261 fCosmiMatchingHisto[6]=new TH1F("Cosmics matching","p1-isPair",100,-10*4.541 ,10*4.541 );
262 fCosmiMatchingHisto[7]=new TH1F("Cosmics matching","p2-isPair",100,-10*0.01134 ,10*0.01134 );
263 fCosmiMatchingHisto[8]=new TH1F("Cosmics matching","p3-isPair",100,-10*0.004644,10*0.004644);
264 fCosmiMatchingHisto[9]=new TH1F("Cosmics matching","p4-isPair",100,-10*0.03773 ,10*0.03773 );
265 BookDistortionMaps();
268 AliTPCcalibTime::~AliTPCcalibTime(){
270 // Virtual Destructor
272 for(Int_t i=0;i<3;i++){
273 if(fHistVdriftLaserA[i]){
274 delete fHistVdriftLaserA[i];
275 fHistVdriftLaserA[i]=NULL;
277 if(fHistVdriftLaserC[i]){
278 delete fHistVdriftLaserC[i];
279 fHistVdriftLaserC[i]=NULL;
283 fArrayDz->SetOwner();
288 for(Int_t i=0;i<5;i++){
289 if(fCosmiMatchingHisto[i]){
290 delete fCosmiMatchingHisto[i];
291 fCosmiMatchingHisto[i]=NULL;
295 for (Int_t i=0;i<5;i++) {
296 delete fResHistoTPCITS[i];
297 delete fResHistoTPCTRD[i];
298 delete fResHistoTPCvertex[i];
299 fResHistoTPCITS[i]=0;
300 fResHistoTPCTRD[i]=0;
301 fResHistoTPCvertex[i]=0;
305 fAlignITSTPC->SetOwner(kTRUE);
306 fAlignTRDTPC->SetOwner(kTRUE);
307 fAlignTOFTPC->SetOwner(kTRUE);
309 fAlignITSTPC->Delete();
310 fAlignTRDTPC->Delete();
311 fAlignTOFTPC->Delete();
317 Bool_t AliTPCcalibTime::IsLaser(const AliESDEvent *const /*event*/){
319 // Indicator is laser event not yet implemented - to be done using trigger info or event specie
321 return kTRUE; //More accurate creteria to be added
323 Bool_t AliTPCcalibTime::IsCosmics(const AliESDEvent *const /*event*/){
325 // Indicator is cosmic event not yet implemented - to be done using trigger info or event specie
328 return kTRUE; //More accurate creteria to be added
330 Bool_t AliTPCcalibTime::IsBeam(const AliESDEvent *const /*event*/){
332 // Indicator is physic event not yet implemented - to be done using trigger info or event specie
335 return kTRUE; //More accurate creteria to be added
337 void AliTPCcalibTime::ResetCurrent(){
338 fDz=0; //Reset current dz
343 void AliTPCcalibTime::Process(AliESDEvent *event){
345 // main function to make calibration
348 if (event->GetNumberOfTracks()<2) return;
350 if(IsLaser (event)) ProcessLaser (event);
351 if(IsCosmics(event)) ProcessCosmic(event);
352 if(IsBeam (event)) ProcessBeam (event);
355 void AliTPCcalibTime::ProcessLaser(AliESDEvent *event){
357 // Fit drift velocity using laser
360 const Int_t kMinTracks = 40; // minimal number of laser tracks
361 const Int_t kMinTracksSide = 20; // minimal number of tracks per side
362 const Float_t kMaxDeltaZ = 30.; // maximal trigger delay
363 const Float_t kMaxDeltaV = 0.05; // maximal deltaV
364 const Float_t kMaxRMS = 0.1; // maximal RMS of tracks
367 TCut cutRMS("sqrt(laserA.fElements[4])<0.1&&sqrt(laserC.fElements[4])<0.1");
368 TCut cutZ("abs(laserA.fElements[0]-laserC.fElements[0])<3");
369 TCut cutV("abs(laserA.fElements[1]-laserC.fElements[1])<0.01");
370 TCut cutY("abs(laserA.fElements[2]-laserC.fElements[2])<2");
371 TCut cutAll = cutRMS+cutZ+cutV+cutY;
373 if (event->GetNumberOfTracks()<kMinTracks) return;
375 if(!fLaser) fLaser = new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
376 fLaser->Process(event);
377 if (fLaser->GetNtracks()<kMinTracks) return; // small amount of tracks cut
378 if (fLaser->fFitAside->GetNrows()==0 && fLaser->fFitCside->GetNrows()==0) return; // no fit neither a or C side
380 // debug streamer - activate stream level
381 // Use it for tuning of the cuts
383 // cuts to be applied
385 Int_t isReject[2]={0,0};
388 if (TMath::Abs((*fLaser->fFitAside)[3]) < kMinTracksSide) isReject[0]|=1;
389 if (TMath::Abs((*fLaser->fFitCside)[3]) < kMinTracksSide) isReject[1]|=1;
390 // unreasonable z offset
391 if (TMath::Abs((*fLaser->fFitAside)[0])>kMaxDeltaZ) isReject[0]|=2;
392 if (TMath::Abs((*fLaser->fFitCside)[0])>kMaxDeltaZ) isReject[1]|=2;
393 // unreasonable drift velocity
394 if (TMath::Abs((*fLaser->fFitAside)[1]-1)>kMaxDeltaV) isReject[0]|=4;
395 if (TMath::Abs((*fLaser->fFitCside)[1]-1)>kMaxDeltaV) isReject[1]|=4;
397 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitAside)[4]))>kMaxRMS ) isReject[0]|=8;
398 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitCside)[4]))>kMaxRMS ) isReject[1]|=8;
404 printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
406 TTreeSRedirector *cstream = GetDebugStreamer();
408 TTimeStamp tstamp(fTime);
409 (*cstream)<<"laserInfo"<<
410 "run="<<fRun<< // run number
411 "event="<<fEvent<< // event number
412 "time="<<fTime<< // time stamp of event
413 "trigger="<<fTrigger<< // trigger
414 "mag="<<fMagF<< // magnetic field
416 "rejectA="<<isReject[0]<<
417 "rejectC="<<isReject[1]<<
418 "laserA.="<<fLaser->fFitAside<<
419 "laserC.="<<fLaser->fFitCside<<
420 "laserAC.="<<fLaser->fFitACside<<
421 "trigger="<<event->GetFiredTriggerClasses()<<
428 TVectorD vdriftA(5), vdriftC(5),vdriftAC(5);
429 vdriftA=*(fLaser->fFitAside);
430 vdriftC=*(fLaser->fFitCside);
431 vdriftAC=*(fLaser->fFitACside);
432 Int_t npointsA=0, npointsC=0;
433 Float_t chi2A=0, chi2C=0;
434 npointsA= TMath::Nint(vdriftA[3]);
436 npointsC= TMath::Nint(vdriftC[3]);
439 TTimeStamp tstamp(fTime);
440 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
441 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
442 Double_t driftA=0, driftC=0;
443 if (vdriftA[1]>1.-kMaxDeltaV) driftA = 1./vdriftA[1]-1.;
444 if (vdriftC[1]>1.-kMaxDeltaV) driftC = 1./vdriftC[1]-1.;
446 Double_t vecDriftLaserA[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftA,event->GetRunNumber()};
447 Double_t vecDriftLaserC[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftC,event->GetRunNumber()};
448 // Double_t vecDrift[4] ={fTime,(ptrelative0+ptrelative1)/2.0,1./((*(fLaser->fFitACside))[1])-1,event->GetRunNumber()};
450 for (Int_t icalib=0;icalib<3;icalib++){
451 if (icalib==0){ //z0 shift
452 vecDriftLaserA[2]=vdriftA[0]/250.;
453 vecDriftLaserC[2]=vdriftC[0]/250.;
455 if (icalib==1){ //vdrel shift
456 vecDriftLaserA[2]=driftA;
457 vecDriftLaserC[2]=driftC;
459 if (icalib==2){ //gy shift - full gy - full drift
460 vecDriftLaserA[2]=vdriftA[2]/250.;
461 vecDriftLaserC[2]=vdriftC[2]/250.;
463 //if (isReject[0]==0) fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
464 //if (isReject[1]==0) fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
465 fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
466 fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
469 // THnSparse* curHist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
470 // TString shortName=curHist->ClassName();
471 // shortName+="_MEAN_DRIFT_LASER_";
477 // name+=event->GetFiredTriggerClasses();
479 // curHist=(THnSparseF*)fArrayDz->FindObject(name);
481 // curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
482 // fArrayDz->AddLast(curHist);
484 // curHist->Fill(vecDrift);
489 // curHist=(THnSparseF*)fArrayDz->FindObject(name);
491 // curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
492 // fArrayDz->AddLast(curHist);
494 // curHist->Fill(vecDrift);
497 void AliTPCcalibTime::ProcessCosmic(const AliESDEvent *const event){
499 // process Cosmic event - track matching A side C side
502 Printf("ERROR: ESD not available");
505 if (event->GetTimeStamp() == 0 ) {
506 Printf("no time stamp!");
513 // Track0 is choosen in upper TPC part
514 // Track1 is choosen in lower TPC part
516 const Int_t kMinClustersCross =30;
517 const Int_t kMinClusters =80;
518 Int_t ntracks=event->GetNumberOfTracks();
519 if (ntracks==0) return;
520 if (ntracks > fCutTracks) return;
522 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
523 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
525 TObjArray tpcSeeds(ntracks);
526 Double_t vtxx[3]={0,0,0};
527 Double_t svtxx[3]={0.000001,0.000001,100.};
528 AliESDVertex vtx(vtxx,svtxx);
532 TArrayI clusterSideA(ntracks);
533 TArrayI clusterSideC(ntracks);
534 for (Int_t i=0;i<ntracks;++i) {
537 AliESDtrack *track = event->GetTrack(i);
539 const AliExternalTrackParam * trackIn = track->GetInnerParam();
540 const AliExternalTrackParam * trackOut = track->GetOuterParam();
541 if (!trackIn) continue;
542 if (!trackOut) continue;
544 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
545 if (!friendTrack) continue;
546 if (friendTrack) ProcessSame(track,friendTrack,event);
547 if (friendTrack) ProcessAlignITS(track,friendTrack,event,esdFriend);
548 if (friendTrack) ProcessAlignTRD(track,friendTrack);
549 if (friendTrack) ProcessAlignTOF(track,friendTrack);
550 TObject *calibObject;
551 AliTPCseed *seed = 0;
552 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
554 tpcSeeds.AddAt(seed,i);
556 for (Int_t irow=159;irow>0;irow--) {
557 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
559 if ((cl->GetDetector()%36)<18) nA++;
560 if ((cl->GetDetector()%36)>=18) nC++;
566 if (ntracks<2) return;
571 for (Int_t i=0;i<ntracks;++i) {
572 AliESDtrack *track0 = event->GetTrack(i);
573 // track0 - choosen upper part
574 if (!track0) continue;
575 if (!track0->GetOuterParam()) continue;
576 if (track0->GetOuterParam()->GetAlpha()<0) continue;
578 track0->GetDirection(d1);
579 for (Int_t j=0;j<ntracks;++j) {
581 AliESDtrack *track1 = event->GetTrack(j);
583 if (!track1) continue;
584 if (!track1->GetOuterParam()) continue;
585 if (track0->GetTPCNcls()+ track1->GetTPCNcls()< kMinClusters) continue;
586 Int_t nAC = TMath::Max( TMath::Min(clusterSideA[i], clusterSideC[j]),
587 TMath::Min(clusterSideC[i], clusterSideA[j]));
588 if (nAC<kMinClustersCross) continue;
589 Int_t nA0=clusterSideA[i];
590 Int_t nC0=clusterSideC[i];
591 Int_t nA1=clusterSideA[j];
592 Int_t nC1=clusterSideC[j];
593 // if (track1->GetOuterParam()->GetAlpha()>0) continue;
596 track1->GetDirection(d2);
598 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
599 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
600 if (! seed0) continue;
601 if (! seed1) continue;
602 Float_t dir = (d1[0]*d2[0] + d1[1]*d2[1] + d1[2]*d2[2]);
603 Float_t dist0 = track0->GetLinearD(0,0);
604 Float_t dist1 = track1->GetLinearD(0,0);
606 // conservative cuts - convergence to be guarantied
607 // applying before track propagation
608 if (TMath::Abs(TMath::Abs(dist0)-TMath::Abs(dist1))>fCutMaxD) continue; // distance to the 0,0
609 if (TMath::Abs(dir)<TMath::Abs(fCutMinDir)) continue; // direction vector product
610 Float_t bz = AliTracker::GetBz();
611 Float_t dvertex0[2]; //distance to 0,0
612 Float_t dvertex1[2]; //distance to 0,0
613 track0->GetDZ(0,0,0,bz,dvertex0);
614 track1->GetDZ(0,0,0,bz,dvertex1);
615 if (TMath::Abs(dvertex0[1])>250) continue;
616 if (TMath::Abs(dvertex1[1])>250) continue;
620 Float_t dmax = TMath::Max(TMath::Abs(dist0),TMath::Abs(dist1));
621 AliExternalTrackParam param0(*track0);
622 AliExternalTrackParam param1(*track1);
624 // Propagate using Magnetic field and correct fo material budget
626 AliTracker::PropagateTrackTo(¶m0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
627 AliTracker::PropagateTrackTo(¶m1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
629 // Propagate rest to the 0,0 DCA - z should be ignored
632 param0.PropagateToDCA(&vtx,bz,1000);
634 param1.PropagateToDCA(&vtx,bz,1000);
635 param0.GetDZ(0,0,0,bz,dvertex0);
636 param1.GetDZ(0,0,0,bz,dvertex1);
641 Bool_t isPair = IsPair(¶m0,¶m1);
642 Bool_t isCross = IsCross(track0, track1);
643 Bool_t isSame = IsSame(track0, track1);
645 THnSparse* hist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
646 TString shortName=hist->ClassName();
647 shortName+="_MEAN_VDRIFT_COSMICS_";
651 if((isSame) || (isCross && isPair)){
652 if (track0->GetTPCNcls()+ track1->GetTPCNcls()> 80) {
653 fDz = param0.GetZ() - param1.GetZ();
654 Double_t sign=(nA0>nA1)? 1:-1;
656 TTimeStamp tstamp(fTime);
657 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
658 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
659 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
660 THnSparse* curHist=NULL;
664 name+=event->GetFiredTriggerClasses();
666 curHist=(THnSparseF*)fArrayDz->FindObject(name);
668 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
669 fArrayDz->AddLast(curHist);
671 // curHist=(THnSparseF*)(fMapDz->GetValue(event->GetFiredTriggerClasses()));
673 // curHist=new THnSparseF(event->GetFiredTriggerClasses(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
674 // fMapDz->Add(new TObjString(event->GetFiredTriggerClasses()),curHist);
676 curHist->Fill(vecDrift);
681 curHist=(THnSparseF*)fArrayDz->FindObject(name);
683 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
684 fArrayDz->AddLast(curHist);
686 // curHist=(THnSparseF*)(fMapDz->GetValue("all"));
688 // curHist=new THnSparseF("all","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
689 // fMapDz->Add(new TObjString("all"),curHist);
691 curHist->Fill(vecDrift);
694 TTreeSRedirector *cstream = GetDebugStreamer();
697 (*cstream)<<"trackInfo"<<
708 "isCross="<<isCross<<
716 } // end 2nd order loop
717 } // end 1st order loop
720 TTreeSRedirector *cstream = GetDebugStreamer();
722 (*cstream)<<"timeInfo"<<
723 "run="<<fRun<< // run number
724 "event="<<fEvent<< // event number
725 "time="<<fTime<< // time stamp of event
726 "trigger="<<fTrigger<< // trigger
727 "mag="<<fMagF<< // magnetic field
728 // Environment values
730 // accumulated values
732 "fDz="<<fDz<< //! current delta z
733 "trigger="<<event->GetFiredTriggerClasses()<<
737 if (GetDebugLevel()>20) printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
740 void AliTPCcalibTime::ProcessBeam(const AliESDEvent *const /*event*/){
742 // Not special treatment yet - the same for cosmic and physic event
746 void AliTPCcalibTime::Analyze(){
748 // Special macro to analyze result of calibration and extract calibration entries
749 // Not yet ported to the Analyze function yet
753 THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name) const
756 // Get histogram for given trigger mask
758 TIterator* iterator = fArrayDz->MakeIterator();
760 TString newName=name;
762 THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
763 THnSparse* addHist=NULL;
764 while((addHist=(THnSparseF*)iterator->Next())){
765 if(!addHist) continue;
766 TString histName=addHist->GetName();
767 if(!histName.Contains(newName)) continue;
769 newHist->Add(addHist);
774 TObjArray* AliTPCcalibTime::GetHistoDrift() const
777 // return array of histograms
782 TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
784 // Make a drift velocity (delta Z) graph
786 THnSparse* histoDrift=GetHistoDrift(name);
787 TGraphErrors* graphDrift=NULL;
789 graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
790 TString end=histoDrift->GetName();
791 Int_t pos=end.Index("_");
792 end=end(pos,end.Capacity()-pos);
793 TString graphName=graphDrift->ClassName();
796 graphDrift->SetName(graphName);
801 TObjArray* AliTPCcalibTime::GetGraphDrift(){
803 // make a array of drift graphs
805 TObjArray* arrayGraphDrift=new TObjArray();
806 TIterator* iterator=fArrayDz->MakeIterator();
808 THnSparse* addHist=NULL;
809 while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
810 return arrayGraphDrift;
813 AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
815 // Make a fit AliSplinefit of drift velocity
817 TGraph* graphDrift=GetGraphDrift(name);
818 AliSplineFit* fitDrift=NULL;
819 if(graphDrift && graphDrift->GetN()){
820 fitDrift=new AliSplineFit();
821 fitDrift->SetGraph(graphDrift);
822 fitDrift->SetMinPoints(graphDrift->GetN()+1);
823 fitDrift->InitKnots(graphDrift,2,0,0.001);
824 fitDrift->SplineFit(0);
825 TString end=graphDrift->GetName();
826 Int_t pos=end.Index("_");
827 end=end(pos,end.Capacity()-pos);
828 TString fitName=fitDrift->ClassName();
831 //fitDrift->SetName(fitName);
839 Long64_t AliTPCcalibTime::Merge(TCollection *const li) {
841 // Object specific merging procedure
843 TIterator* iter = li->MakeIterator();
844 AliTPCcalibTime* cal = 0;
846 while ((cal = (AliTPCcalibTime*)iter->Next())) {
847 if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
848 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
851 for (Int_t imeas=0; imeas<3; imeas++){
852 if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
853 fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
854 fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
858 for (Int_t imeas=0; imeas<5; imeas++){
859 if (cal->GetResHistoTPCITS(imeas) && cal->GetResHistoTPCITS(imeas)){
860 fResHistoTPCITS[imeas]->Add(cal->fResHistoTPCITS[imeas]);
861 fResHistoTPCvertex[imeas]->Add(cal->fResHistoTPCvertex[imeas]);
862 fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
865 TObjArray* addArray=cal->GetHistoDrift();
866 if(!addArray) return 0;
867 TIterator* iterator = addArray->MakeIterator();
869 THnSparse* addHist=NULL;
870 while((addHist=(THnSparseF*)iterator->Next())){
871 if(!addHist) continue;
873 THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
875 localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
876 fArrayDz->AddLast(localHist);
878 localHist->Add(addHist);
881 for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
885 for (Int_t itype=0; itype<3; itype++){
890 if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
891 if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
892 if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
894 if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
895 if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
896 arr0->Expand(arr1->GetEntriesFast());
898 for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
899 AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
900 AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
901 if (!kalman1) continue;
902 if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
903 kalman0->SetRejectOutliers(kFALSE);
904 kalman0->Merge(kalman1);
912 Bool_t AliTPCcalibTime::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
914 // 0. Same direction - OPOSITE - cutDir +cutT
915 TCut cutDir("cutDir","dir<-0.99")
917 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
920 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
922 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
925 const Double_t *p0 = tr0->GetParameter();
926 const Double_t *p1 = tr1->GetParameter();
927 fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
928 fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
929 fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
930 fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
931 fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
933 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
934 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
935 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
936 Double_t d0[3], d1[3];
937 tr0->GetDirection(d0);
938 tr1->GetDirection(d1);
939 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
941 fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
942 fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
943 fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
944 fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
945 fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
949 Bool_t AliTPCcalibTime::IsCross(AliESDtrack *const tr0, AliESDtrack *const tr1){
951 // check if the cosmic pair of tracks crossed A/C side
953 Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
954 if (result==kFALSE) return result;
959 Bool_t AliTPCcalibTime::IsSame(AliESDtrack *const tr0, AliESDtrack *const tr1){
961 // track crossing the CE
962 // 0. minimal number of clusters
963 // 1. Same sector +-1
964 // 2. Inner and outer track param on opposite side
965 // 3. Outer and inner track parameter close each to other
969 // inner and outer on opposite sides in z
971 const Int_t knclCut0 = 30;
972 const Double_t kalphaCut = 0.4;
974 // 0. minimal number of clusters
976 if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
977 if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
979 // 1. alpha cut - sector+-1
981 if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
985 if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
986 if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
992 const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
993 const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
994 const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
995 const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
997 if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD) result&=kFALSE;
998 if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
999 if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1000 if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1001 if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD) result&=kFALSE;
1002 if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1003 if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1004 if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1006 result=kTRUE; // just to put break point here
1012 void AliTPCcalibTime::ProcessSame(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
1014 // Process TPC tracks crossing CE
1016 // 0. Select only track crossing the CE
1017 // 1. Cut on the track length
1018 // 2. Refit the terack on A and C side separatelly
1019 // 3. Fill time histograms
1020 const Int_t kMinNcl=100;
1021 const Int_t kMinNclS=25; // minimul number of clusters on the sides
1022 if (!friendTrack->GetTPCOut()) return;
1024 // 0. Select only track crossing the CE
1026 if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1028 // 1. cut on track length
1030 if (track->GetTPCNcls()<kMinNcl) return;
1032 // 2. Refit track sepparatel on A and C side
1034 TObject *calibObject;
1035 AliTPCseed *seed = 0;
1036 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1037 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1041 AliExternalTrackParam trackIn(*track->GetInnerParam());
1042 AliExternalTrackParam trackOut(*track->GetOuterParam());
1043 Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1044 Double_t xyz[3]={0,0.,0.0};
1046 Int_t nclIn=0,nclOut=0;
1047 trackIn.ResetCovariance(30.);
1048 trackOut.ResetCovariance(30.);
1052 for (Int_t irow=0;irow<159;irow++) {
1053 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1055 if (cl->GetX()<80) continue;
1056 if (track->GetInnerParam()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1057 if (track->GetInnerParam()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1058 Int_t sector = cl->GetDetector();
1059 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1060 if (TMath::Abs(dalpha)>0.01){
1061 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1063 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1064 trackIn.GetXYZ(xyz);
1065 bz = AliTracker::GetBz(xyz);
1066 if (!trackIn.PropagateTo(r[0],bz)) break;
1068 trackIn.Update(&r[1],cov);
1073 for (Int_t irow=159;irow>0;irow--) {
1074 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1076 if (cl->GetX()<80) continue;
1077 if (cl->GetZ()*track->GetOuterParam()->GetZ()<0) break;
1078 if (friendTrack->GetTPCOut()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1079 if (friendTrack->GetTPCOut()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1080 Int_t sector = cl->GetDetector();
1081 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1082 if (TMath::Abs(dalpha)>0.01){
1083 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1085 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1086 trackOut.GetXYZ(xyz);
1087 bz = AliTracker::GetBz(xyz);
1088 if (!trackOut.PropagateTo(r[0],bz)) break;
1090 trackOut.Update(&r[1],cov);
1092 trackOut.Rotate(trackIn.GetAlpha());
1093 Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1094 trackIn.PropagateTo(meanX,bz);
1095 trackOut.PropagateTo(meanX,bz);
1096 TTreeSRedirector *cstream = GetDebugStreamer();
1099 trackIn.GetXYZ(gxyz.GetMatrixArray());
1100 TTimeStamp tstamp(fTime);
1101 (*cstream)<<"tpctpc"<<
1102 "run="<<fRun<< // run number
1103 "event="<<fEvent<< // event number
1104 "time="<<fTime<< // time stamp of event
1105 "trigger="<<fTrigger<< // trigger
1106 "mag="<<fMagF<< // magnetic field
1108 "xyz.="<<&gxyz<< // global position
1109 "tIn.="<<&trackIn<< // refitterd track in
1110 "tOut.="<<&trackOut<< // refitter track out
1111 "nclIn="<<nclIn<< //
1112 "nclOut="<<nclOut<< //
1116 // 3. Fill time histograms
1117 // Debug stremaer expression
1118 // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1119 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1120 fDz = trackOut.GetZ()-trackIn.GetZ();
1121 if (trackOut.GetTgl()<0) fDz*=-1.;
1122 TTimeStamp tstamp(fTime);
1123 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1124 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1125 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1127 // fill histograms per trigger class and itegrated
1129 THnSparse* curHist=NULL;
1130 for (Int_t itype=0; itype<2; itype++){
1131 TString name="MEAN_VDRIFT_CROSS_";
1133 name+=event->GetFiredTriggerClasses();
1138 curHist=(THnSparseF*)fArrayDz->FindObject(name);
1140 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1141 fArrayDz->AddLast(curHist);
1143 curHist->Fill(vecDrift);
1149 void AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
1151 // Process track - Update TPC-ITS alignment
1153 // 0. Apply standartd cuts
1154 // 1. Recalucluate the current statistic median/RMS
1155 // 2. Apply median+-rms cut
1156 // 3. Update kalman filter
1158 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1159 const Int_t kMinITS = 3; // minimal number of ITS cluster
1160 const Double_t kMinZ = 10; // maximal dz distance
1161 const Double_t kMaxDy = 2.; // maximal dy distance
1162 const Double_t kMaxAngle= 0.015; // maximal angular distance
1163 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1164 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1165 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1166 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1167 const Double_t kMinPt = 0.3; // minimal pt
1168 const Int_t kN=50; // deepnes of history
1169 static Int_t kglast=0;
1170 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1173 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";
1176 // 0. Apply standard cuts
1178 Int_t dummycl[1000];
1179 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1180 if (track->GetITSclusters(dummycl)<kMinITS) return; // minimal amount of clusters
1181 if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
1182 if (!friendTrack->GetITSOut()) return;
1183 if (!track->GetInnerParam()) return;
1184 if (!track->GetOuterParam()) return;
1185 if (track->GetInnerParam()->Pt()<kMinPt) return;
1186 // exclude crossing track
1187 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1188 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1189 if (track->GetInnerParam()->GetX()>90) return;
1191 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
1192 AliExternalTrackParam pITS(*(friendTrack->GetITSOut())); // ITS standalone if possible
1193 AliExternalTrackParam pITS2(*(friendTrack->GetITSOut())); //TPC-ITS track
1194 pITS2.Rotate(pTPC.GetAlpha());
1195 // pITS2.PropagateTo(pTPC.GetX(),fMagF);
1196 AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1198 AliESDfriendTrack *itsfriendTrack=0;
1200 // try to find standalone ITS track corresponing to the TPC if possible
1202 Bool_t hasAlone=kFALSE;
1203 Int_t ntracks=event->GetNumberOfTracks();
1204 for (Int_t i=0; i<ntracks; i++){
1205 AliESDtrack *trackS = event->GetTrack(i);
1206 if (trackS->GetTPCNcls()>0) continue; //continue if has TPC info
1207 itsfriendTrack = esdFriend->GetTrack(i);
1208 if (!itsfriendTrack) continue;
1209 if (!itsfriendTrack->GetITSOut()) continue;
1210 if (TMath::Abs(pITS2.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1211 pITS=(*(itsfriendTrack->GetITSOut()));
1213 pITS.Rotate(pTPC.GetAlpha());
1214 //pITS.PropagateTo(pTPC.GetX(),fMagF);
1215 AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
1216 if (TMath::Abs(pITS2.GetY()-pITS.GetY())> kMaxDy) continue;
1219 if (!hasAlone) pITS=pITS2;
1221 if (TMath::Abs(pITS.GetY()-pTPC.GetY()) >kMaxDy) return;
1222 if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1223 if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1225 // 1. Update median and RMS info
1227 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1228 TVectorD vecDeltaN(5);
1229 Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1231 for (Int_t i=0;i<4;i++){
1232 vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1233 kgdP[i][kglast%kN]=vecDelta[i];
1236 Int_t entries=(kglast<kN)?kglast:kN;
1237 for (Int_t i=0;i<4;i++){
1238 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1239 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1242 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1243 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1247 // 2. Apply median+-rms cut
1249 if (kglast<3) return; //median and RMS to be defined
1250 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1252 // 3. Update alignment
1254 Int_t htime = fTime/3600; //time in hours
1255 if (fAlignITSTPC->GetEntriesFast()<htime){
1256 fAlignITSTPC->Expand(htime*2+20);
1258 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1260 // make Alignment object if doesn't exist
1261 align=new AliRelAlignerKalman();
1262 align->SetRunNumber(fRun);
1263 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1264 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1265 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1266 align->SetRejectOutliers(kFALSE);
1268 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1269 align->SetMagField(fMagF);
1270 fAlignITSTPC->AddAt(align,htime);
1272 align->AddTrackParams(&pITS,&pTPC);
1273 align->SetTimeStamp(fTime);
1274 align->SetRunNumber(fRun );
1275 Float_t dca[2],cov[3];
1276 track->GetImpactParameters(dca,cov);
1277 if (TMath::Abs(dca[0])<kMaxDy){
1278 FillResHistoTPCITS(&pTPC,&pITS);
1279 FillResHistoTPC(track);
1282 Int_t nupdates=align->GetNUpdates();
1283 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1284 align->SetRejectOutliers(kFALSE);
1285 TTreeSRedirector *cstream = GetDebugStreamer();
1286 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1287 TVectorD gpTPC(3), gdTPC(3);
1288 TVectorD gpITS(3), gdITS(3);
1289 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1290 pTPC.GetDirection(gdTPC.GetMatrixArray());
1291 pITS.GetXYZ(gpITS.GetMatrixArray());
1292 pITS.GetDirection(gdITS.GetMatrixArray());
1293 (*cstream)<<"itstpc"<<
1294 "run="<<fRun<< // run number
1295 "event="<<fEvent<< // event number
1296 "time="<<fTime<< // time stamp of event
1297 "trigger="<<fTrigger<< // trigger
1298 "mag="<<fMagF<< // magnetic field
1300 "hasAlone="<<hasAlone<< // has ITS standalone ?
1301 "track.="<<track<< // track info
1302 "nmed="<<kglast<< // number of entries to define median and RMS
1303 "vMed.="<<&vecMedian<< // median of deltas
1304 "vRMS.="<<&vecRMS<< // rms of deltas
1305 "vDelta.="<<&vecDelta<< // delta in respect to median
1306 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1307 "t.="<<track<< // ful track - find proper cuts
1308 "a.="<<align<< // current alignment
1309 "pITS.="<<&pITS<< // track param ITS
1310 "pITS2.="<<&pITS2<< // track param ITS+TPC
1311 "pTPC.="<<&pTPC<< // track param TPC
1312 "gpTPC.="<<&gpTPC<< // global position TPC
1313 "gdTPC.="<<&gdTPC<< // global direction TPC
1314 "gpITS.="<<&gpITS<< // global position ITS
1315 "gdITS.="<<&gdITS<< // global position ITS
1323 void AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1325 // Process track - Update TPC-TRD alignment
1327 // 0. Apply standartd cuts
1328 // 1. Recalucluate the current statistic median/RMS
1329 // 2. Apply median+-rms cut
1330 // 3. Update kalman filter
1332 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1333 const Int_t kMinTRD = 50; // minimal number of TRD cluster
1334 const Double_t kMinZ = 20; // maximal dz distance
1335 const Double_t kMaxDy = 5.; // maximal dy distance
1336 const Double_t kMaxAngle= 0.1; // maximal angular distance
1337 const Double_t kSigmaCut= 10; // maximal sigma distance to median
1338 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1339 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1340 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1341 const Double_t kRefX = 275; // reference X
1342 const Int_t kN=50; // deepnes of history
1343 static Int_t kglast=0;
1344 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1346 // 0. Apply standard cuts
1348 Int_t dummycl[1000];
1349 if (track->GetTRDclusters(dummycl)<kMinTRD) return; // minimal amount of clusters
1350 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1351 if (!friendTrack->GetTRDIn()) return;
1352 if (!track->IsOn(AliESDtrack::kTRDrefit)) return;
1353 if (!track->IsOn(AliESDtrack::kTRDout)) return;
1354 if (!track->GetInnerParam()) return;
1355 if (!friendTrack->GetTPCOut()) return;
1356 // exclude crossing track
1357 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1358 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1360 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1361 AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
1362 AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
1363 pTRD.Rotate(pTPC.GetAlpha());
1364 // pTRD.PropagateTo(pTPC.GetX(),fMagF);
1365 AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1367 ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.; // increas sys errors
1368 ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1370 if (TMath::Abs(pTRD.GetY()-pTPC.GetY()) >kMaxDy) return;
1371 if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1372 // if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1374 // 1. Update median and RMS info
1376 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1377 TVectorD vecDeltaN(5);
1378 Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1380 for (Int_t i=0;i<4;i++){
1381 vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1382 kgdP[i][kglast%kN]=vecDelta[i];
1385 Int_t entries=(kglast<kN)?kglast:kN;
1386 for (Int_t i=0;i<4;i++){
1387 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1389 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1392 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1393 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1397 // 2. Apply median+-rms cut
1399 if (kglast<3) return; //median and RMS to be defined
1400 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1402 // 3. Update alignment
1404 Int_t htime = fTime/3600; //time in hours
1405 if (fAlignTRDTPC->GetEntriesFast()<htime){
1406 fAlignTRDTPC->Expand(htime*2+20);
1408 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1410 // make Alignment object if doesn't exist
1411 align=new AliRelAlignerKalman();
1412 align->SetRunNumber(fRun);
1413 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1414 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1415 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1416 align->SetRejectOutliers(kFALSE);
1417 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1418 align->SetMagField(fMagF);
1419 fAlignTRDTPC->AddAt(align,htime);
1421 align->AddTrackParams(&pTRD,&pTPC);
1422 align->SetTimeStamp(fTime);
1423 align->SetRunNumber(fRun );
1424 FillResHistoTPCTRD(&pTPC,&pTRD);
1426 Int_t nupdates=align->GetNUpdates();
1427 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1428 align->SetRejectOutliers(kFALSE);
1429 TTreeSRedirector *cstream = GetDebugStreamer();
1430 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1431 TVectorD gpTPC(3), gdTPC(3);
1432 TVectorD gpTRD(3), gdTRD(3);
1433 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1434 pTPC.GetDirection(gdTPC.GetMatrixArray());
1435 pTRD.GetXYZ(gpTRD.GetMatrixArray());
1436 pTRD.GetDirection(gdTRD.GetMatrixArray());
1437 (*cstream)<<"trdtpc"<<
1438 "run="<<fRun<< // run number
1439 "event="<<fEvent<< // event number
1440 "time="<<fTime<< // time stamp of event
1441 "trigger="<<fTrigger<< // trigger
1442 "mag="<<fMagF<< // magnetic field
1444 "nmed="<<kglast<< // number of entries to define median and RMS
1445 "vMed.="<<&vecMedian<< // median of deltas
1446 "vRMS.="<<&vecRMS<< // rms of deltas
1447 "vDelta.="<<&vecDelta<< // delta in respect to median
1448 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1449 "t.="<<track<< // ful track - find proper cuts
1450 "a.="<<align<< // current alignment
1451 "pTRD.="<<&pTRD<< // track param TRD
1452 "pTPC.="<<&pTPC<< // track param TPC
1453 "gpTPC.="<<&gpTPC<< // global position TPC
1454 "gdTPC.="<<&gdTPC<< // global direction TPC
1455 "gpTRD.="<<&gpTRD<< // global position TRD
1456 "gdTRD.="<<&gdTRD<< // global position TRD
1462 void AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1465 // Process track - Update TPC-TOF alignment
1467 // -1. Make a TOF "track"
1468 // 0. Apply standartd cuts
1469 // 1. Recalucluate the current statistic median/RMS
1470 // 2. Apply median+-rms cut
1471 // 3. Update kalman filter
1473 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1474 // const Double_t kMinZ = 10; // maximal dz distance
1475 const Double_t kMaxDy = 5.; // maximal dy distance
1476 const Double_t kMaxAngle= 0.015; // maximal angular distance
1477 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1478 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1479 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1481 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1482 const Int_t kN=50; // deepnes of history
1483 static Int_t kglast=0;
1484 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1486 // -1. Make a TOF track-
1487 // Clusters are not in friends - use alingment points
1489 if (track->GetTOFsignal()<=0) return;
1490 if (!friendTrack->GetTPCOut()) return;
1491 if (!track->GetInnerParam()) return;
1492 if (!friendTrack->GetTPCOut()) return;
1493 const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1494 if (!points) return;
1495 AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
1496 AliExternalTrackParam pTOF(pTPC);
1497 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1498 Int_t npoints = points->GetNPoints();
1499 AliTrackPoint point;
1502 for (Int_t ipoint=0;ipoint<npoints;ipoint++){
1503 points->GetPoint(point,ipoint);
1506 Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1507 if (r<350) continue;
1508 if (r>400) continue;
1509 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1510 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);
1511 AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1512 pTPC.PropagateTo(lpoint.GetX(),fMagF);
1514 ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1515 ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
1516 ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1517 ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
1518 ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1519 ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
1522 if (naccept==0) return; // no tof match clusters
1524 // 0. Apply standard cuts
1526 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1527 // exclude crossing track
1528 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1530 if (TMath::Abs(pTOF.GetY()-pTPC.GetY()) >kMaxDy) return;
1531 if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1532 if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1534 // 1. Update median and RMS info
1536 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1537 TVectorD vecDeltaN(5);
1538 Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1540 for (Int_t i=0;i<4;i++){
1541 vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1542 kgdP[i][kglast%kN]=vecDelta[i];
1545 Int_t entries=(kglast<kN)?kglast:kN;
1547 for (Int_t i=0;i<4;i++){
1548 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1549 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1552 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1553 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1554 if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
1558 // 2. Apply median+-rms cut
1560 if (kglast<10) return; //median and RMS to be defined
1563 // 3. Update alignment
1565 Int_t htime = fTime/3600; //time in hours
1566 if (fAlignTOFTPC->GetEntriesFast()<htime){
1567 fAlignTOFTPC->Expand(htime*2+20);
1569 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1571 // make Alignment object if doesn't exist
1572 align=new AliRelAlignerKalman();
1573 align->SetRunNumber(fRun);
1574 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1575 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1576 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1577 align->SetRejectOutliers(kFALSE);
1578 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1579 align->SetMagField(fMagF);
1580 fAlignTOFTPC->AddAt(align,htime);
1582 align->AddTrackParams(&pTOF,&pTPC);
1583 align->SetTimeStamp(fTime);
1584 align->SetRunNumber(fRun );
1586 Int_t nupdates=align->GetNUpdates();
1587 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1588 align->SetRejectOutliers(kFALSE);
1589 TTreeSRedirector *cstream = GetDebugStreamer();
1590 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1591 TVectorD gpTPC(3), gdTPC(3);
1592 TVectorD gpTOF(3), gdTOF(3);
1593 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1594 pTPC.GetDirection(gdTPC.GetMatrixArray());
1595 pTOF.GetXYZ(gpTOF.GetMatrixArray());
1596 pTOF.GetDirection(gdTOF.GetMatrixArray());
1597 (*cstream)<<"toftpc"<<
1598 "run="<<fRun<< // run number
1599 "event="<<fEvent<< // event number
1600 "time="<<fTime<< // time stamp of event
1601 "trigger="<<fTrigger<< // trigger
1602 "mag="<<fMagF<< // magnetic field
1604 "nmed="<<kglast<< // number of entries to define median and RMS
1605 "vMed.="<<&vecMedian<< // median of deltas
1606 "vRMS.="<<&vecRMS<< // rms of deltas
1607 "vDelta.="<<&vecDelta<< // delta in respect to median
1608 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1609 "t.="<<track<< // ful track - find proper cuts
1610 "a.="<<align<< // current alignment
1611 "pTOF.="<<&pTOF<< // track param TOF
1612 "pTPC.="<<&pTPC<< // track param TPC
1613 "gpTPC.="<<&gpTPC<< // global position TPC
1614 "gdTPC.="<<&gdTPC<< // global direction TPC
1615 "gpTOF.="<<&gpTOF<< // global position TOF
1616 "gdTOF.="<<&gdTOF<< // global position TOF
1622 void AliTPCcalibTime::BookDistortionMaps(){
1624 // Book ndimensional histograms of distortions/residuals
1625 // Only primary tracks are selected for analysis
1628 Double_t xminTrack[4], xmaxTrack[4];
1630 TString axisName[4];
1631 TString axisTitle[4];
1634 axisName[0] ="#Delta";
1635 axisTitle[0] ="#Delta";
1638 xminTrack[1] =-1.5; xmaxTrack[1]=1.5;
1639 axisName[1] ="tanTheta";
1640 axisTitle[1] ="tan(#Theta)";
1643 xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi();
1645 axisTitle[2] ="#phi";
1648 xminTrack[3] =-1.; xmaxTrack[3]=1.; // 0.33 GeV cut
1652 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
1653 fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1654 fResHistoTPCvertex[0] = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1655 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
1656 fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1659 xminTrack[0] =-3.; xmaxTrack[0]=3.; //
1660 fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1661 fResHistoTPCvertex[1] = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1662 fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1665 xminTrack[0] =-0.015; xmaxTrack[0]=0.015; //
1666 fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1667 fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1668 fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1671 xminTrack[0] =-0.025; xmaxTrack[0]=0.025; //
1672 fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1673 fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1674 fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1677 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
1678 fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1679 fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1680 fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1682 for (Int_t ivar=0;ivar<4;ivar++){
1683 for (Int_t ivar2=0;ivar2<4;ivar2++){
1684 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1685 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1686 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1687 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1688 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1689 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1695 void AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
1697 // fill residual histograms pTPCIn-pITSOut
1698 // Histogram is filled only for primary tracks
1702 pTPCIn->GetXYZ(xyz);
1703 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1704 histoX[1]= pTPCIn->GetTgl();
1706 histoX[3]= pTPCIn->GetSnp();
1707 AliExternalTrackParam lits(*pITSOut);
1708 lits.Rotate(pTPCIn->GetAlpha());
1709 lits.PropagateTo(pTPCIn->GetX(),fMagF);
1711 for (Int_t ihisto=0; ihisto<5; ihisto++){
1712 histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
1713 fResHistoTPCITS[ihisto]->Fill(histoX);
1718 void AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
1720 // fill residual histograms pTPC - vertex
1721 // Histogram is filled only for primary tracks
1724 const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
1725 AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
1727 AliExternalTrackParam lits(*pTrack);
1728 if (TMath::Abs(pTrack->GetY())>3) return; // beam pipe
1729 pTPCvertex.Rotate(lits.GetAlpha());
1730 //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
1731 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
1732 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
1734 pTPCIn->GetXYZ(xyz);
1735 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1736 histoX[1]= pTPCIn->GetTgl();
1738 histoX[3]= pTPCIn->GetSnp();
1740 Float_t dca[2], cov[3];
1741 pTrack->GetImpactParametersTPC(dca,cov);
1742 for (Int_t ihisto=0; ihisto<5; ihisto++){
1743 histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
1744 // if (ihisto<2) histoX[0]=dca[ihisto];
1745 fResHistoTPCvertex[ihisto]->Fill(histoX);
1750 void AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
1752 // fill resuidual histogram TPCout-TRDin
1756 pTPCOut->GetXYZ(xyz);
1757 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1758 histoX[1]= pTPCOut->GetTgl();
1760 histoX[3]= pTPCOut->GetSnp();
1762 AliExternalTrackParam ltrd(*pTRDIn);
1763 ltrd.Rotate(pTPCOut->GetAlpha());
1764 // ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
1765 AliTracker::PropagateTrackToBxByBz(<rd,pTPCOut->GetX(),0.1,0.1,kFALSE);
1767 for (Int_t ihisto=0; ihisto<5; ihisto++){
1768 histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
1769 fResHistoTPCTRD[ihisto]->Fill(histoX);