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]);
866 fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
867 fResHistoTPCTOF[imeas]->Add(cal->fResHistoTPCTOF[imeas]);
870 TObjArray* addArray=cal->GetHistoDrift();
871 if(!addArray) return 0;
872 TIterator* iterator = addArray->MakeIterator();
874 THnSparse* addHist=NULL;
875 while((addHist=(THnSparseF*)iterator->Next())){
876 if(!addHist) continue;
878 THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
880 localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
881 fArrayDz->AddLast(localHist);
883 localHist->Add(addHist);
886 for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
890 for (Int_t itype=0; itype<3; itype++){
895 if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
896 if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
897 if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
899 if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
900 if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
901 arr0->Expand(arr1->GetEntriesFast());
903 for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
904 AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
905 AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
906 if (!kalman1) continue;
907 if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
908 kalman0->SetRejectOutliers(kFALSE);
909 kalman0->Merge(kalman1);
917 Bool_t AliTPCcalibTime::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
919 // 0. Same direction - OPOSITE - cutDir +cutT
920 TCut cutDir("cutDir","dir<-0.99")
922 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
925 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
927 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
930 const Double_t *p0 = tr0->GetParameter();
931 const Double_t *p1 = tr1->GetParameter();
932 fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
933 fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
934 fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
935 fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
936 fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
938 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
939 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
940 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
941 Double_t d0[3], d1[3];
942 tr0->GetDirection(d0);
943 tr1->GetDirection(d1);
944 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
946 fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
947 fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
948 fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
949 fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
950 fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
954 Bool_t AliTPCcalibTime::IsCross(AliESDtrack *const tr0, AliESDtrack *const tr1){
956 // check if the cosmic pair of tracks crossed A/C side
958 Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
959 if (result==kFALSE) return result;
964 Bool_t AliTPCcalibTime::IsSame(AliESDtrack *const tr0, AliESDtrack *const tr1){
966 // track crossing the CE
967 // 0. minimal number of clusters
968 // 1. Same sector +-1
969 // 2. Inner and outer track param on opposite side
970 // 3. Outer and inner track parameter close each to other
974 // inner and outer on opposite sides in z
976 const Int_t knclCut0 = 30;
977 const Double_t kalphaCut = 0.4;
979 // 0. minimal number of clusters
981 if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
982 if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
984 // 1. alpha cut - sector+-1
986 if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
990 if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
991 if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
997 const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
998 const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
999 const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1000 const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1002 if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD) result&=kFALSE;
1003 if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1004 if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1005 if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1006 if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD) result&=kFALSE;
1007 if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1008 if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1009 if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1011 result=kTRUE; // just to put break point here
1017 void AliTPCcalibTime::ProcessSame(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
1019 // Process TPC tracks crossing CE
1021 // 0. Select only track crossing the CE
1022 // 1. Cut on the track length
1023 // 2. Refit the terack on A and C side separatelly
1024 // 3. Fill time histograms
1025 const Int_t kMinNcl=100;
1026 const Int_t kMinNclS=25; // minimul number of clusters on the sides
1027 if (!friendTrack->GetTPCOut()) return;
1029 // 0. Select only track crossing the CE
1031 if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1033 // 1. cut on track length
1035 if (track->GetTPCNcls()<kMinNcl) return;
1037 // 2. Refit track sepparatel on A and C side
1039 TObject *calibObject;
1040 AliTPCseed *seed = 0;
1041 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1042 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1046 AliExternalTrackParam trackIn(*track->GetInnerParam());
1047 AliExternalTrackParam trackOut(*track->GetOuterParam());
1048 Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1049 Double_t xyz[3]={0,0.,0.0};
1051 Int_t nclIn=0,nclOut=0;
1052 trackIn.ResetCovariance(30.);
1053 trackOut.ResetCovariance(30.);
1057 for (Int_t irow=0;irow<159;irow++) {
1058 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1060 if (cl->GetX()<80) continue;
1061 if (track->GetInnerParam()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1062 if (track->GetInnerParam()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1063 Int_t sector = cl->GetDetector();
1064 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1065 if (TMath::Abs(dalpha)>0.01){
1066 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1068 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1069 trackIn.GetXYZ(xyz);
1070 bz = AliTracker::GetBz(xyz);
1071 if (!trackIn.PropagateTo(r[0],bz)) break;
1073 trackIn.Update(&r[1],cov);
1078 for (Int_t irow=159;irow>0;irow--) {
1079 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1081 if (cl->GetX()<80) continue;
1082 if (cl->GetZ()*track->GetOuterParam()->GetZ()<0) break;
1083 if (friendTrack->GetTPCOut()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1084 if (friendTrack->GetTPCOut()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1085 Int_t sector = cl->GetDetector();
1086 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1087 if (TMath::Abs(dalpha)>0.01){
1088 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1090 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1091 trackOut.GetXYZ(xyz);
1092 bz = AliTracker::GetBz(xyz);
1093 if (!trackOut.PropagateTo(r[0],bz)) break;
1095 trackOut.Update(&r[1],cov);
1097 trackOut.Rotate(trackIn.GetAlpha());
1098 Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1099 trackIn.PropagateTo(meanX,bz);
1100 trackOut.PropagateTo(meanX,bz);
1101 TTreeSRedirector *cstream = GetDebugStreamer();
1104 trackIn.GetXYZ(gxyz.GetMatrixArray());
1105 TTimeStamp tstamp(fTime);
1106 (*cstream)<<"tpctpc"<<
1107 "run="<<fRun<< // run number
1108 "event="<<fEvent<< // event number
1109 "time="<<fTime<< // time stamp of event
1110 "trigger="<<fTrigger<< // trigger
1111 "mag="<<fMagF<< // magnetic field
1113 "xyz.="<<&gxyz<< // global position
1114 "tIn.="<<&trackIn<< // refitterd track in
1115 "tOut.="<<&trackOut<< // refitter track out
1116 "nclIn="<<nclIn<< //
1117 "nclOut="<<nclOut<< //
1121 // 3. Fill time histograms
1122 // Debug stremaer expression
1123 // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1124 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1125 fDz = trackOut.GetZ()-trackIn.GetZ();
1126 if (trackOut.GetTgl()<0) fDz*=-1.;
1127 TTimeStamp tstamp(fTime);
1128 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1129 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1130 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1132 // fill histograms per trigger class and itegrated
1134 THnSparse* curHist=NULL;
1135 for (Int_t itype=0; itype<2; itype++){
1136 TString name="MEAN_VDRIFT_CROSS_";
1138 name+=event->GetFiredTriggerClasses();
1143 curHist=(THnSparseF*)fArrayDz->FindObject(name);
1145 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1146 fArrayDz->AddLast(curHist);
1148 curHist->Fill(vecDrift);
1154 void AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
1156 // Process track - Update TPC-ITS alignment
1158 // 0. Apply standartd cuts
1159 // 1. Recalucluate the current statistic median/RMS
1160 // 2. Apply median+-rms cut
1161 // 3. Update kalman filter
1163 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1164 const Int_t kMinITS = 3; // minimal number of ITS cluster
1165 const Double_t kMinZ = 10; // maximal dz distance
1166 const Double_t kMaxDy = 2.; // maximal dy distance
1167 const Double_t kMaxAngle= 0.015; // maximal angular distance
1168 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1169 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1170 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1171 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1172 const Double_t kMinPt = 0.3; // minimal pt
1173 const Int_t kN=50; // deepnes of history
1174 static Int_t kglast=0;
1175 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1178 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";
1181 // 0. Apply standard cuts
1183 Int_t dummycl[1000];
1184 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1185 if (track->GetITSclusters(dummycl)<kMinITS) return; // minimal amount of clusters
1186 if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
1187 if (!friendTrack->GetITSOut()) return;
1188 if (!track->GetInnerParam()) return;
1189 if (!track->GetOuterParam()) return;
1190 if (track->GetInnerParam()->Pt()<kMinPt) return;
1191 // exclude crossing track
1192 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1193 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1194 if (track->GetInnerParam()->GetX()>90) return;
1196 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
1197 AliExternalTrackParam pITS(*(friendTrack->GetITSOut())); // ITS standalone if possible
1198 AliExternalTrackParam pITS2(*(friendTrack->GetITSOut())); //TPC-ITS track
1199 pITS2.Rotate(pTPC.GetAlpha());
1200 // pITS2.PropagateTo(pTPC.GetX(),fMagF);
1201 AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1203 AliESDfriendTrack *itsfriendTrack=0;
1205 // try to find standalone ITS track corresponing to the TPC if possible
1207 Bool_t hasAlone=kFALSE;
1208 Int_t ntracks=event->GetNumberOfTracks();
1209 for (Int_t i=0; i<ntracks; i++){
1210 AliESDtrack *trackS = event->GetTrack(i);
1211 if (trackS->GetTPCNcls()>0) continue; //continue if has TPC info
1212 itsfriendTrack = esdFriend->GetTrack(i);
1213 if (!itsfriendTrack) continue;
1214 if (!itsfriendTrack->GetITSOut()) continue;
1215 if (TMath::Abs(pITS2.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1216 pITS=(*(itsfriendTrack->GetITSOut()));
1218 pITS.Rotate(pTPC.GetAlpha());
1219 //pITS.PropagateTo(pTPC.GetX(),fMagF);
1220 AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
1221 if (TMath::Abs(pITS2.GetY()-pITS.GetY())> kMaxDy) continue;
1224 if (!hasAlone) pITS=pITS2;
1226 if (TMath::Abs(pITS.GetY()-pTPC.GetY()) >kMaxDy) return;
1227 if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1228 if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1230 // 1. Update median and RMS info
1232 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1233 TVectorD vecDeltaN(5);
1234 Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1236 for (Int_t i=0;i<4;i++){
1237 vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1238 kgdP[i][kglast%kN]=vecDelta[i];
1241 Int_t entries=(kglast<kN)?kglast:kN;
1242 for (Int_t i=0;i<4;i++){
1243 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1244 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1247 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1248 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1252 // 2. Apply median+-rms cut
1254 if (kglast<3) return; //median and RMS to be defined
1255 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1257 // 3. Update alignment
1259 Int_t htime = fTime/3600; //time in hours
1260 if (fAlignITSTPC->GetEntriesFast()<htime){
1261 fAlignITSTPC->Expand(htime*2+20);
1263 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1265 // make Alignment object if doesn't exist
1266 align=new AliRelAlignerKalman();
1267 align->SetRunNumber(fRun);
1268 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1269 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1270 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1271 align->SetRejectOutliers(kFALSE);
1273 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1274 align->SetMagField(fMagF);
1275 fAlignITSTPC->AddAt(align,htime);
1277 align->AddTrackParams(&pITS,&pTPC);
1278 align->SetTimeStamp(fTime);
1279 align->SetRunNumber(fRun );
1280 Float_t dca[2],cov[3];
1281 track->GetImpactParameters(dca,cov);
1282 if (TMath::Abs(dca[0])<kMaxDy){
1283 FillResHistoTPCITS(&pTPC,&pITS);
1284 FillResHistoTPC(track);
1287 Int_t nupdates=align->GetNUpdates();
1288 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1289 align->SetRejectOutliers(kFALSE);
1290 TTreeSRedirector *cstream = GetDebugStreamer();
1291 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1292 TVectorD gpTPC(3), gdTPC(3);
1293 TVectorD gpITS(3), gdITS(3);
1294 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1295 pTPC.GetDirection(gdTPC.GetMatrixArray());
1296 pITS.GetXYZ(gpITS.GetMatrixArray());
1297 pITS.GetDirection(gdITS.GetMatrixArray());
1298 (*cstream)<<"itstpc"<<
1299 "run="<<fRun<< // run number
1300 "event="<<fEvent<< // event number
1301 "time="<<fTime<< // time stamp of event
1302 "trigger="<<fTrigger<< // trigger
1303 "mag="<<fMagF<< // magnetic field
1305 "hasAlone="<<hasAlone<< // has ITS standalone ?
1306 "track.="<<track<< // track info
1307 "nmed="<<kglast<< // number of entries to define median and RMS
1308 "vMed.="<<&vecMedian<< // median of deltas
1309 "vRMS.="<<&vecRMS<< // rms of deltas
1310 "vDelta.="<<&vecDelta<< // delta in respect to median
1311 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1312 "t.="<<track<< // ful track - find proper cuts
1313 "a.="<<align<< // current alignment
1314 "pITS.="<<&pITS<< // track param ITS
1315 "pITS2.="<<&pITS2<< // track param ITS+TPC
1316 "pTPC.="<<&pTPC<< // track param TPC
1317 "gpTPC.="<<&gpTPC<< // global position TPC
1318 "gdTPC.="<<&gdTPC<< // global direction TPC
1319 "gpITS.="<<&gpITS<< // global position ITS
1320 "gdITS.="<<&gdITS<< // global position ITS
1328 void AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1330 // Process track - Update TPC-TRD alignment
1332 // 0. Apply standartd cuts
1333 // 1. Recalucluate the current statistic median/RMS
1334 // 2. Apply median+-rms cut
1335 // 3. Update kalman filter
1337 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1338 const Int_t kMinTRD = 50; // minimal number of TRD cluster
1339 const Double_t kMinZ = 20; // maximal dz distance
1340 const Double_t kMaxDy = 5.; // maximal dy distance
1341 const Double_t kMaxAngle= 0.1; // maximal angular distance
1342 const Double_t kSigmaCut= 10; // maximal sigma distance to median
1343 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1344 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1345 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1346 const Double_t kRefX = 275; // reference X
1347 const Int_t kN=50; // deepnes of history
1348 static Int_t kglast=0;
1349 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1351 // 0. Apply standard cuts
1353 Int_t dummycl[1000];
1354 if (track->GetTRDclusters(dummycl)<kMinTRD) return; // minimal amount of clusters
1355 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1356 if (!friendTrack->GetTRDIn()) return;
1357 if (!track->IsOn(AliESDtrack::kTRDrefit)) return;
1358 if (!track->IsOn(AliESDtrack::kTRDout)) return;
1359 if (!track->GetInnerParam()) return;
1360 if (!friendTrack->GetTPCOut()) return;
1361 // exclude crossing track
1362 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1363 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1365 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1366 AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
1367 AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
1368 pTRD.Rotate(pTPC.GetAlpha());
1369 // pTRD.PropagateTo(pTPC.GetX(),fMagF);
1370 AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1372 ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.; // increas sys errors
1373 ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1375 if (TMath::Abs(pTRD.GetY()-pTPC.GetY()) >kMaxDy) return;
1376 if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1377 // if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1379 // 1. Update median and RMS info
1381 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1382 TVectorD vecDeltaN(5);
1383 Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1385 for (Int_t i=0;i<4;i++){
1386 vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1387 kgdP[i][kglast%kN]=vecDelta[i];
1390 Int_t entries=(kglast<kN)?kglast:kN;
1391 for (Int_t i=0;i<4;i++){
1392 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1394 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1397 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1398 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1402 // 2. Apply median+-rms cut
1404 if (kglast<3) return; //median and RMS to be defined
1405 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1407 // 3. Update alignment
1409 Int_t htime = fTime/3600; //time in hours
1410 if (fAlignTRDTPC->GetEntriesFast()<htime){
1411 fAlignTRDTPC->Expand(htime*2+20);
1413 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1415 // make Alignment object if doesn't exist
1416 align=new AliRelAlignerKalman();
1417 align->SetRunNumber(fRun);
1418 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1419 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1420 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1421 align->SetRejectOutliers(kFALSE);
1422 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1423 align->SetMagField(fMagF);
1424 fAlignTRDTPC->AddAt(align,htime);
1426 align->AddTrackParams(&pTRD,&pTPC);
1427 align->SetTimeStamp(fTime);
1428 align->SetRunNumber(fRun );
1429 Float_t dca[2],cov[3];
1430 track->GetImpactParameters(dca,cov);
1431 if (TMath::Abs(dca[0])<kMaxDy){
1432 FillResHistoTPCTRD(&pTPC,&pTRD); //only primaries
1435 Int_t nupdates=align->GetNUpdates();
1436 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1437 align->SetRejectOutliers(kFALSE);
1438 TTreeSRedirector *cstream = GetDebugStreamer();
1439 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1440 TVectorD gpTPC(3), gdTPC(3);
1441 TVectorD gpTRD(3), gdTRD(3);
1442 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1443 pTPC.GetDirection(gdTPC.GetMatrixArray());
1444 pTRD.GetXYZ(gpTRD.GetMatrixArray());
1445 pTRD.GetDirection(gdTRD.GetMatrixArray());
1446 (*cstream)<<"trdtpc"<<
1447 "run="<<fRun<< // run number
1448 "event="<<fEvent<< // event number
1449 "time="<<fTime<< // time stamp of event
1450 "trigger="<<fTrigger<< // trigger
1451 "mag="<<fMagF<< // magnetic field
1453 "nmed="<<kglast<< // number of entries to define median and RMS
1454 "vMed.="<<&vecMedian<< // median of deltas
1455 "vRMS.="<<&vecRMS<< // rms of deltas
1456 "vDelta.="<<&vecDelta<< // delta in respect to median
1457 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1458 "t.="<<track<< // ful track - find proper cuts
1459 "a.="<<align<< // current alignment
1460 "pTRD.="<<&pTRD<< // track param TRD
1461 "pTPC.="<<&pTPC<< // track param TPC
1462 "gpTPC.="<<&gpTPC<< // global position TPC
1463 "gdTPC.="<<&gdTPC<< // global direction TPC
1464 "gpTRD.="<<&gpTRD<< // global position TRD
1465 "gdTRD.="<<&gdTRD<< // global position TRD
1471 void AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1474 // Process track - Update TPC-TOF alignment
1476 // -1. Make a TOF "track"
1477 // 0. Apply standartd cuts
1478 // 1. Recalucluate the current statistic median/RMS
1479 // 2. Apply median+-rms cut
1480 // 3. Update kalman filter
1482 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1483 // const Double_t kMinZ = 10; // maximal dz distance
1484 const Double_t kMaxDy = 5.; // maximal dy distance
1485 const Double_t kMaxAngle= 0.015; // maximal angular distance
1486 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1487 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1488 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1490 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1491 const Int_t kN=50; // deepnes of history
1492 static Int_t kglast=0;
1493 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1495 // -1. Make a TOF track-
1496 // Clusters are not in friends - use alingment points
1498 if (track->GetTOFsignal()<=0) return;
1499 if (!friendTrack->GetTPCOut()) return;
1500 if (!track->GetInnerParam()) return;
1501 if (!friendTrack->GetTPCOut()) return;
1502 const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1503 if (!points) return;
1504 AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
1505 AliExternalTrackParam pTOF(pTPC);
1506 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1507 Int_t npoints = points->GetNPoints();
1508 AliTrackPoint point;
1511 for (Int_t ipoint=0;ipoint<npoints;ipoint++){
1512 points->GetPoint(point,ipoint);
1515 Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1516 if (r<350) continue;
1517 if (r>400) continue;
1518 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1519 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);
1520 AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1521 pTPC.PropagateTo(lpoint.GetX(),fMagF);
1523 ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1524 ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
1525 ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1526 ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
1527 ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1528 ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
1531 if (naccept==0) return; // no tof match clusters
1533 // 0. Apply standard cuts
1535 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1536 // exclude crossing track
1537 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1539 if (TMath::Abs(pTOF.GetY()-pTPC.GetY()) >kMaxDy) return;
1540 if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1541 if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1543 // 1. Update median and RMS info
1545 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1546 TVectorD vecDeltaN(5);
1547 Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1549 for (Int_t i=0;i<4;i++){
1550 vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1551 kgdP[i][kglast%kN]=vecDelta[i];
1554 Int_t entries=(kglast<kN)?kglast:kN;
1556 for (Int_t i=0;i<4;i++){
1557 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1558 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1561 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1562 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1563 if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
1567 // 2. Apply median+-rms cut
1569 if (kglast<10) return; //median and RMS to be defined
1572 // 3. Update alignment
1574 Int_t htime = fTime/3600; //time in hours
1575 if (fAlignTOFTPC->GetEntriesFast()<htime){
1576 fAlignTOFTPC->Expand(htime*2+20);
1578 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1580 // make Alignment object if doesn't exist
1581 align=new AliRelAlignerKalman();
1582 align->SetRunNumber(fRun);
1583 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1584 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1585 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1586 align->SetRejectOutliers(kFALSE);
1587 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1588 align->SetMagField(fMagF);
1589 fAlignTOFTPC->AddAt(align,htime);
1591 align->AddTrackParams(&pTOF,&pTPC);
1592 Float_t dca[2],cov[3];
1593 track->GetImpactParameters(dca,cov);
1594 if (TMath::Abs(dca[0])<kMaxDy){
1595 FillResHistoTPCTOF(&pTPC,&pTOF);
1597 align->SetTimeStamp(fTime);
1598 align->SetRunNumber(fRun );
1600 Int_t nupdates=align->GetNUpdates();
1601 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1602 align->SetRejectOutliers(kFALSE);
1603 TTreeSRedirector *cstream = GetDebugStreamer();
1604 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1605 TVectorD gpTPC(3), gdTPC(3);
1606 TVectorD gpTOF(3), gdTOF(3);
1607 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1608 pTPC.GetDirection(gdTPC.GetMatrixArray());
1609 pTOF.GetXYZ(gpTOF.GetMatrixArray());
1610 pTOF.GetDirection(gdTOF.GetMatrixArray());
1611 (*cstream)<<"toftpc"<<
1612 "run="<<fRun<< // run number
1613 "event="<<fEvent<< // event number
1614 "time="<<fTime<< // time stamp of event
1615 "trigger="<<fTrigger<< // trigger
1616 "mag="<<fMagF<< // magnetic field
1618 "nmed="<<kglast<< // number of entries to define median and RMS
1619 "vMed.="<<&vecMedian<< // median of deltas
1620 "vRMS.="<<&vecRMS<< // rms of deltas
1621 "vDelta.="<<&vecDelta<< // delta in respect to median
1622 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1623 "t.="<<track<< // ful track - find proper cuts
1624 "a.="<<align<< // current alignment
1625 "pTOF.="<<&pTOF<< // track param TOF
1626 "pTPC.="<<&pTPC<< // track param TPC
1627 "gpTPC.="<<&gpTPC<< // global position TPC
1628 "gdTPC.="<<&gdTPC<< // global direction TPC
1629 "gpTOF.="<<&gpTOF<< // global position TOF
1630 "gdTOF.="<<&gdTOF<< // global position TOF
1636 void AliTPCcalibTime::BookDistortionMaps(){
1638 // Book ndimensional histograms of distortions/residuals
1639 // Only primary tracks are selected for analysis
1642 Double_t xminTrack[4], xmaxTrack[4];
1644 TString axisName[4];
1645 TString axisTitle[4];
1648 axisName[0] ="#Delta";
1649 axisTitle[0] ="#Delta";
1652 xminTrack[1] =-1.5; xmaxTrack[1]=1.5;
1653 axisName[1] ="tanTheta";
1654 axisTitle[1] ="tan(#Theta)";
1657 xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi();
1659 axisTitle[2] ="#phi";
1662 xminTrack[3] =-1.; xmaxTrack[3]=1.; // 0.33 GeV cut
1666 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
1667 fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1668 fResHistoTPCvertex[0] = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1669 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
1670 fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1671 xminTrack[0] =-5; xmaxTrack[0]=5; //
1672 fResHistoTPCTOF[0] = new THnSparseS("TPCTOF#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1675 xminTrack[0] =-3.; xmaxTrack[0]=3.; //
1676 fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1677 fResHistoTPCvertex[1] = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1678 fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1679 xminTrack[0] =-5.; xmaxTrack[0]=5.; //
1680 fResHistoTPCTOF[1] = new THnSparseS("TPCTOF#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
1683 xminTrack[0] =-0.015; xmaxTrack[0]=0.015; //
1684 fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1685 fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1686 fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1687 fResHistoTPCTOF[2] = new THnSparseS("TPCTOF#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
1690 xminTrack[0] =-0.025; xmaxTrack[0]=0.025; //
1691 fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1692 fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1693 fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1694 fResHistoTPCTOF[3] = new THnSparseS("TPCTOF#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
1697 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
1698 fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1699 fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1700 fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1701 fResHistoTPCTOF[4] = new THnSparseS("TPCTOF#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
1703 for (Int_t ivar=0;ivar<4;ivar++){
1704 for (Int_t ivar2=0;ivar2<4;ivar2++){
1705 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1706 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1707 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1708 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1709 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
1710 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
1716 void AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
1718 // fill residual histograms pTPCIn-pITSOut
1719 // Histogram is filled only for primary tracks
1723 pTPCIn->GetXYZ(xyz);
1724 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1725 histoX[1]= pTPCIn->GetTgl();
1727 histoX[3]= pTPCIn->GetSnp();
1728 AliExternalTrackParam lits(*pITSOut);
1729 lits.Rotate(pTPCIn->GetAlpha());
1730 lits.PropagateTo(pTPCIn->GetX(),fMagF);
1732 for (Int_t ihisto=0; ihisto<5; ihisto++){
1733 histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
1734 fResHistoTPCITS[ihisto]->Fill(histoX);
1739 void AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
1741 // fill residual histograms pTPC - vertex
1742 // Histogram is filled only for primary tracks
1745 const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
1746 AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
1748 AliExternalTrackParam lits(*pTrack);
1749 if (TMath::Abs(pTrack->GetY())>3) return; // beam pipe
1750 pTPCvertex.Rotate(lits.GetAlpha());
1751 //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
1752 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
1753 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
1755 pTPCIn->GetXYZ(xyz);
1756 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1757 histoX[1]= pTPCIn->GetTgl();
1759 histoX[3]= pTPCIn->GetSnp();
1761 Float_t dca[2], cov[3];
1762 pTrack->GetImpactParametersTPC(dca,cov);
1763 for (Int_t ihisto=0; ihisto<5; ihisto++){
1764 histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
1765 // if (ihisto<2) histoX[0]=dca[ihisto];
1766 fResHistoTPCvertex[ihisto]->Fill(histoX);
1771 void AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
1773 // fill resuidual histogram TPCout-TRDin
1777 pTPCOut->GetXYZ(xyz);
1778 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1779 histoX[1]= pTPCOut->GetTgl();
1781 histoX[3]= pTPCOut->GetSnp();
1783 AliExternalTrackParam ltrd(*pTRDIn);
1784 ltrd.Rotate(pTPCOut->GetAlpha());
1785 // ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
1786 AliTracker::PropagateTrackToBxByBz(<rd,pTPCOut->GetX(),0.1,0.1,kFALSE);
1788 for (Int_t ihisto=0; ihisto<5; ihisto++){
1789 histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
1790 fResHistoTPCTRD[ihisto]->Fill(histoX);
1795 void AliTPCcalibTime::FillResHistoTPCTOF(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTOFIn ){
1797 // fill resuidual histogram TPCout-TOFin
1798 // track propagated to the TOF position
1802 AliExternalTrackParam ltpc(*pTPCOut);
1803 ltpc.Rotate(pTOFIn->GetAlpha());
1804 AliTracker::PropagateTrackToBxByBz(<pc,pTOFIn->GetX(),0.1,0.1,kFALSE);
1807 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
1808 histoX[1]= ltpc.GetTgl();
1810 histoX[3]= ltpc.GetSnp();
1812 for (Int_t ihisto=0; ihisto<2; ihisto++){
1813 histoX[0]=ltpc.GetParameter()[ihisto]-pTOFIn->GetParameter()[ihisto];
1814 fResHistoTPCTOF[ihisto]->Fill(histoX);