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);
27 #include "Riostream.h"
28 #include "TDatabasePDG.h"
29 #include "TGraphErrors.h"
31 #include "THnSparse.h"
34 #include "TTimeStamp.h"
40 #include "AliDCSSensor.h"
41 #include "AliDCSSensorArray.h"
42 #include "AliESDEvent.h"
43 #include "AliESDInputHandler.h"
44 #include "AliESDVertex.h"
45 #include "AliESDfriend.h"
47 #include "AliRelAlignerKalman.h"
48 #include "AliTPCCalROC.h"
49 #include "AliTPCParam.h"
50 #include "AliTPCTracklet.h"
51 #include "AliTPCcalibDB.h"
52 #include "AliTPCcalibLaser.h"
53 #include "AliTPCcalibTime.h"
54 #include "AliTPCclusterMI.h"
55 #include "AliTPCseed.h"
56 #include "AliTrackPointArray.h"
57 #include "AliTracker.h"
58 #include "AliKFVertex.h"
61 ClassImp(AliTPCcalibTime)
63 Double_t AliTPCcalibTime::fgResHistoMergeCut = 10000000.;
65 AliTPCcalibTime::AliTPCcalibTime()
67 fMemoryMode(1), // 0 -do not fill THnSparse with residuals 1- fill only important QA THn 2 - Fill all THnsparse for calibration
68 fLaser(0), // pointer to laser calibration
69 fDz(0), // current delta z
70 fCutMaxD(3), // maximal distance in rfi ditection
71 fCutMaxDz(25), // maximal distance in rfi ditection
72 fCutTheta(0.03), // maximal distan theta
73 fCutMinDir(-0.99), // direction vector products
75 fArrayLaserA(0), //laser fit parameters C
76 fArrayLaserC(0), //laser fit parameters A
77 fArrayDz(0), //NEW! Tmap of V drifts for different triggers
78 fAlignITSTPC(0), //alignemnt array ITS TPC match
79 fAlignTRDTPC(0), //alignemnt array TRD TPC match
80 fAlignTOFTPC(0), //alignemnt array TOF TPC match
81 fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
96 // default constructor
98 AliDebug(5,"Default Constructor");
99 for (Int_t i=0;i<3;i++) {
100 fHistVdriftLaserA[i]=0;
101 fHistVdriftLaserC[i]=0;
103 for (Int_t i=0;i<10;i++) {
104 fCosmiMatchingHisto[i]=0;
107 for (Int_t i=0;i<5;i++) {
109 fResHistoTPCITS[i]=0;
110 fResHistoTPCTRD[i]=0;
111 fResHistoTPCTOF[i]=0;
112 fResHistoTPCvertex[i]=0;
115 for (Int_t i=0;i<12;i++) {
118 for (Int_t i=0;i<5;i++) {
119 fTPCVertexCorrelation[i]=0;
121 static Int_t counter=0;
125 AliDebug(5,Form("Counter Constructor\t%d\t%d",counter,time));
131 AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeVdrift, Int_t memoryMode)
133 fMemoryMode(memoryMode), // 0 -do not fill THnSparse with residuals 1- fill only important QA THn 2 - Fill all THnsparse for calibration
134 fLaser(0), // pointer to laser calibration
135 fDz(0), // current delta z
136 fCutMaxD(5*0.5356), // maximal distance in rfi ditection
137 fCutMaxDz(40), // maximal distance in rfi ditection
138 fCutTheta(5*0.004644),// maximal distan theta
139 fCutMinDir(-0.99), // direction vector products
141 fArrayLaserA(new TObjArray(1000)), //laser fit parameters C
142 fArrayLaserC(new TObjArray(1000)), //laser fit parameters A
143 fArrayDz(0), //Tmap of V drifts for different triggers
144 fAlignITSTPC(0), //alignemnt array ITS TPC match
145 fAlignTRDTPC(0), //alignemnt array TRD TPC match
146 fAlignTOFTPC(0), //alignemnt array TOF TPC match
147 fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
162 // Non deafaul constructor - to be used in the Calibration setups
167 for (Int_t i=0;i<3;i++) {
168 fHistVdriftLaserA[i]=0;
169 fHistVdriftLaserC[i]=0;
172 for (Int_t i=0;i<5;i++) {
174 fResHistoTPCITS[i]=0;
175 fResHistoTPCTRD[i]=0;
176 fResHistoTPCTOF[i]=0;
177 fResHistoTPCvertex[i]=0;
181 AliDebug(5,"Non Default Constructor");
182 fTimeBins =(EndTime-StartTime)/deltaIntegrationTimeVdrift;
183 fTimeStart =StartTime; //(((TObjString*)(mapGRP->GetValue("fAliceStartTime")))->GetString()).Atoi();
184 fTimeEnd =EndTime; //(((TObjString*)(mapGRP->GetValue("fAliceStopTime")))->GetString()).Atoi();
195 Int_t binsVdriftLaser[4] = {fTimeBins , fPtBins , fVdriftBins*20, fRunBins };
196 Double_t xminVdriftLaser[4] = {fTimeStart, fPtStart, fVdriftStart , fRunStart};
197 Double_t xmaxVdriftLaser[4] = {fTimeEnd , fPtEnd , fVdriftEnd , fRunEnd };
198 TString axisTitle[4]={
204 TString histoName[3]={
211 for (Int_t i=0;i<3;i++) {
212 fHistVdriftLaserA[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
213 fHistVdriftLaserC[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
214 fHistVdriftLaserA[i]->SetName(histoName[i]);
215 fHistVdriftLaserC[i]->SetName(histoName[i]);
216 for (Int_t iaxis=0; iaxis<4;iaxis++){
217 fHistVdriftLaserA[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
218 fHistVdriftLaserC[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
221 fBinsVdrift[0] = fTimeBins;
222 fBinsVdrift[1] = fPtBins;
223 fBinsVdrift[2] = fVdriftBins;
224 fBinsVdrift[3] = fRunBins;
225 fXminVdrift[0] = fTimeStart;
226 fXminVdrift[1] = fPtStart;
227 fXminVdrift[2] = fVdriftStart;
228 fXminVdrift[3] = fRunStart;
229 fXmaxVdrift[0] = fTimeEnd;
230 fXmaxVdrift[1] = fPtEnd;
231 fXmaxVdrift[2] = fVdriftEnd;
232 fXmaxVdrift[3] = fRunEnd;
234 fArrayDz=new TObjArray();
235 fAlignITSTPC = new TObjArray; //alignemnt array ITS TPC match
236 fAlignTRDTPC = new TObjArray; //alignemnt array ITS TPC match
237 fAlignTOFTPC = new TObjArray; //alignemnt array ITS TPC match
238 fAlignITSTPC->SetOwner(kTRUE);
239 fAlignTRDTPC->SetOwner(kTRUE);
240 fAlignTOFTPC->SetOwner(kTRUE);
243 fCosmiMatchingHisto[0]=new TH1F("Cosmics matching","p0-all" ,100,-10*0.5356 ,10*0.5356 );
244 fCosmiMatchingHisto[1]=new TH1F("Cosmics matching","p1-all" ,100,-10*4.541 ,10*4.541 );
245 fCosmiMatchingHisto[2]=new TH1F("Cosmics matching","p2-all" ,100,-10*0.01134 ,10*0.01134 );
246 fCosmiMatchingHisto[3]=new TH1F("Cosmics matching","p3-all" ,100,-10*0.004644,10*0.004644);
247 fCosmiMatchingHisto[4]=new TH1F("Cosmics matching","p4-all" ,100,-10*0.03773 ,10*0.03773 );
248 fCosmiMatchingHisto[5]=new TH1F("Cosmics matching","p0-isPair",100,-10*0.5356 ,10*0.5356 );
249 fCosmiMatchingHisto[6]=new TH1F("Cosmics matching","p1-isPair",100,-10*4.541 ,10*4.541 );
250 fCosmiMatchingHisto[7]=new TH1F("Cosmics matching","p2-isPair",100,-10*0.01134 ,10*0.01134 );
251 fCosmiMatchingHisto[8]=new TH1F("Cosmics matching","p3-isPair",100,-10*0.004644,10*0.004644);
252 fCosmiMatchingHisto[9]=new TH1F("Cosmics matching","p4-isPair",100,-10*0.03773 ,10*0.03773 );
253 for (Int_t i=0;i<12;i++) {
256 for (Int_t i=0;i<5;i++) {
257 fTPCVertexCorrelation[i]=0;
259 BookDistortionMaps();
263 AliTPCcalibTime::~AliTPCcalibTime(){
265 // Virtual Destructor
267 static Int_t counter=0;
271 AliDebug(5,Form("Counter Destructor\t%s\t%d\t%d",GetName(),counter,time));
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<10;i++){
291 if(fCosmiMatchingHisto[i]){
292 delete fCosmiMatchingHisto[i];
293 fCosmiMatchingHisto[i]=NULL;
297 for (Int_t i=0;i<5;i++) {
298 delete fResHistoTPCCE[i];
299 delete fResHistoTPCITS[i];
300 delete fResHistoTPCTRD[i];
301 delete fResHistoTPCTOF[i];
302 delete fResHistoTPCvertex[i];
304 fResHistoTPCITS[i]=0;
305 fResHistoTPCTRD[i]=0;
306 fResHistoTPCTOF[i]=0;
307 fResHistoTPCvertex[i]=0;
310 for (Int_t i=0;i<12;i++) if (fTPCVertex[i]) delete fTPCVertex[i];
311 for (Int_t i=0;i<5;i++) if (fTPCVertexCorrelation[i]) delete fTPCVertexCorrelation[i];
314 fAlignITSTPC->SetOwner(kTRUE);
315 fAlignTRDTPC->SetOwner(kTRUE);
316 fAlignTOFTPC->SetOwner(kTRUE);
318 fAlignITSTPC->Delete();
319 fAlignTRDTPC->Delete();
320 fAlignTOFTPC->Delete();
327 fArrayLaserA->SetOwner();
328 fArrayLaserA->Delete();
333 fArrayLaserC->SetOwner();
334 fArrayLaserC->Delete();
340 // Bool_t AliTPCcalibTime::IsLaser(const AliESDEvent *const /*event*/) const{
342 // // Indicator is laser event not yet implemented - to be done using trigger info or event specie
344 // return kTRUE; //More accurate creteria to be added
346 // Bool_t AliTPCcalibTime::IsCosmics(const AliESDEvent *const /*event*/){
348 // // Indicator is cosmic event not yet implemented - to be done using trigger info or event specie
351 // return kTRUE; //More accurate creteria to be added
353 // Bool_t AliTPCcalibTime::IsBeam(const AliESDEvent *const /*event*/) const{
355 // // Indicator is physic event not yet implemented - to be done using trigger info or event specie
358 // return kTRUE; //More accurate creteria to be added
360 void AliTPCcalibTime::ResetCurrent(){
364 fDz=0; //Reset current dz
369 void AliTPCcalibTime::Process(AliESDEvent *event){
371 // main function to make calibration
374 if (event->GetNumberOfTracks()<2) return;
375 AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
379 if (ESDfriend->TestSkipBit()) return;
382 //if(IsLaser (event))
383 ProcessLaser (event);
384 //if(IsCosmics(event))
385 ProcessCosmic(event);
390 void AliTPCcalibTime::ProcessLaser(AliESDEvent *event){
392 // Fit drift velocity using laser
395 const Int_t kMinTracks = 40; // minimal number of laser tracks
396 const Int_t kMinTracksSide = 20; // minimal number of tracks per side
397 const Float_t kMaxDeltaZ = 30.; // maximal trigger delay
398 const Float_t kMaxDeltaV = 0.05; // maximal deltaV
399 const Float_t kMaxRMS = 0.1; // maximal RMS of tracks
402 TCut cutRMS("sqrt(laserA.fElements[4])<0.1&&sqrt(laserC.fElements[4])<0.1");
403 TCut cutZ("abs(laserA.fElements[0]-laserC.fElements[0])<3");
404 TCut cutV("abs(laserA.fElements[1]-laserC.fElements[1])<0.01");
405 TCut cutY("abs(laserA.fElements[2]-laserC.fElements[2])<2");
406 TCut cutAll = cutRMS+cutZ+cutV+cutY;
408 if (event->GetNumberOfTracks()<kMinTracks) return;
410 if(!fLaser) fLaser = new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
411 fLaser->Process(event);
412 if (fLaser->GetNtracks()<kMinTracks) return; // small amount of tracks cut
413 if (fLaser->fFitAside->GetNrows()==0 && fLaser->fFitCside->GetNrows()==0) return; // no fit neither a or C side
415 // debug streamer - activate stream level
416 // Use it for tuning of the cuts
418 // cuts to be applied
420 Int_t isReject[2]={0,0};
423 if (TMath::Abs((*fLaser->fFitAside)[3]) < kMinTracksSide) isReject[0]|=1;
424 if (TMath::Abs((*fLaser->fFitCside)[3]) < kMinTracksSide) isReject[1]|=1;
425 // unreasonable z offset
426 if (TMath::Abs((*fLaser->fFitAside)[0])>kMaxDeltaZ) isReject[0]|=2;
427 if (TMath::Abs((*fLaser->fFitCside)[0])>kMaxDeltaZ) isReject[1]|=2;
428 // unreasonable drift velocity
429 if (TMath::Abs((*fLaser->fFitAside)[1]-1)>kMaxDeltaV) isReject[0]|=4;
430 if (TMath::Abs((*fLaser->fFitCside)[1]-1)>kMaxDeltaV) isReject[1]|=4;
432 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitAside)[4]))>kMaxRMS ) isReject[0]|=8;
433 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitCside)[4]))>kMaxRMS ) isReject[1]|=8;
439 printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
441 TTreeSRedirector *cstream = GetDebugStreamer();
443 TTimeStamp tstamp(fTime);
444 (*cstream)<<"laserInfo"<<
445 "run="<<fRun<< // run number
446 "event="<<fEvent<< // event number
447 "time="<<fTime<< // time stamp of event
448 "trigger="<<fTrigger<< // trigger
449 "mag="<<fMagF<< // magnetic field
451 "rejectA="<<isReject[0]<<
452 "rejectC="<<isReject[1]<<
453 "laserA.="<<fLaser->fFitAside<<
454 "laserC.="<<fLaser->fFitCside<<
455 "laserAC.="<<fLaser->fFitACside<<
456 "trigger="<<event->GetFiredTriggerClasses()<<
463 TVectorD vdriftA(5), vdriftC(5),vdriftAC(6);
464 vdriftA=*(fLaser->fFitAside);
465 vdriftC=*(fLaser->fFitCside);
466 vdriftAC=*(fLaser->fFitACside);
467 Int_t npointsA=0, npointsC=0;
468 Float_t chi2A=0, chi2C=0;
469 npointsA= TMath::Nint(vdriftA[3]);
471 npointsC= TMath::Nint(vdriftC[3]);
474 if (npointsA>kMinTracksSide || npointsC>kMinTracksSide){
475 TVectorD *fitA = new TVectorD(6);
476 TVectorD *fitC = new TVectorD(6);
477 for (Int_t ipar=0; ipar<5; ipar++){
478 (*fitA)[ipar]=vdriftA[ipar];
479 (*fitC)[ipar]=vdriftC[ipar];
483 fArrayLaserA->AddLast(fitA);
484 fArrayLaserC->AddLast(fitC);
488 TTimeStamp tstamp(fTime);
489 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
490 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
491 Double_t driftA=0, driftC=0;
492 if (vdriftA[1]>1.-kMaxDeltaV) driftA = 1./vdriftA[1]-1.;
493 if (vdriftC[1]>1.-kMaxDeltaV) driftC = 1./vdriftC[1]-1.;
495 Double_t vecDriftLaserA[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftA,event->GetRunNumber()};
496 Double_t vecDriftLaserC[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftC,event->GetRunNumber()};
497 // Double_t vecDrift[4] ={fTime,(ptrelative0+ptrelative1)/2.0,1./((*(fLaser->fFitACside))[1])-1,event->GetRunNumber()};
499 for (Int_t icalib=0;icalib<3;icalib++){
500 if (icalib==0){ //z0 shift
501 vecDriftLaserA[2]=vdriftA[0]/250.;
502 vecDriftLaserC[2]=vdriftC[0]/250.;
504 if (icalib==1){ //vdrel shift
505 vecDriftLaserA[2]=driftA;
506 vecDriftLaserC[2]=driftC;
508 if (icalib==2){ //gy shift - full gy - full drift
509 vecDriftLaserA[2]=vdriftA[2]/250.;
510 vecDriftLaserC[2]=vdriftC[2]/250.;
512 //if (isReject[0]==0) fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
513 //if (isReject[1]==0) fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
514 fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
515 fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
519 void AliTPCcalibTime::ProcessCosmic(const AliESDEvent *const event){
521 // process Cosmic event - track matching A side C side
524 Printf("ERROR: ESD not available");
527 if (event->GetTimeStamp() == 0 ) {
528 Printf("no time stamp!");
535 // Track0 is choosen in upper TPC part
536 // Track1 is choosen in lower TPC part
538 const Int_t kMinClustersCross =30;
539 const Int_t kMinClusters =80;
540 Int_t ntracks=event->GetNumberOfTracks();
541 if (ntracks==0) return;
542 if (ntracks > fCutTracks) return;
544 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
545 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
547 TObjArray tpcSeeds(ntracks);
548 Double_t vtxx[3]={0,0,0};
549 Double_t svtxx[3]={0.000001,0.000001,100.};
550 AliESDVertex vtx(vtxx,svtxx);
554 TArrayI clusterSideA(ntracks);
555 TArrayI clusterSideC(ntracks);
556 for (Int_t i=0;i<ntracks;++i) {
559 AliESDtrack *track = event->GetTrack(i);
561 const AliExternalTrackParam * trackIn = track->GetInnerParam();
562 const AliExternalTrackParam * trackOut = track->GetOuterParam();
563 if (!trackIn) continue;
564 if (!trackOut) continue;
566 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
567 if (!friendTrack) continue;
568 if (friendTrack) ProcessSame(track,friendTrack,event);
569 if (friendTrack) ProcessAlignITS(track,friendTrack,event,esdFriend);
570 if (friendTrack) ProcessAlignTRD(track,friendTrack);
571 if (friendTrack) ProcessAlignTOF(track,friendTrack);
572 TObject *calibObject;
573 AliTPCseed *seed = 0;
574 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
576 tpcSeeds.AddAt(seed,i);
578 for (Int_t irow=159;irow>0;irow--) {
579 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
581 if ((cl->GetDetector()%36)<18) nA++;
582 if ((cl->GetDetector()%36)>=18) nC++;
588 if (ntracks<2) return;
593 for (Int_t i=0;i<ntracks;++i) {
594 AliESDtrack *track0 = event->GetTrack(i);
595 // track0 - choosen upper part
596 if (!track0) continue;
597 if (!track0->GetOuterParam()) continue;
598 if (track0->GetOuterParam()->GetAlpha()<0) continue;
600 track0->GetDirection(d1);
601 for (Int_t j=0;j<ntracks;++j) {
603 AliESDtrack *track1 = event->GetTrack(j);
605 if (!track1) continue;
606 if (!track1->GetOuterParam()) continue;
607 if (track0->GetTPCNcls()+ track1->GetTPCNcls()< kMinClusters) continue;
608 Int_t nAC = TMath::Max( TMath::Min(clusterSideA[i], clusterSideC[j]),
609 TMath::Min(clusterSideC[i], clusterSideA[j]));
610 if (nAC<kMinClustersCross) continue;
611 Int_t nA0=clusterSideA[i];
612 Int_t nC0=clusterSideC[i];
613 Int_t nA1=clusterSideA[j];
614 Int_t nC1=clusterSideC[j];
615 // if (track1->GetOuterParam()->GetAlpha()>0) continue;
618 track1->GetDirection(d2);
620 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
621 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
622 if (! seed0) continue;
623 if (! seed1) continue;
624 Float_t dir = (d1[0]*d2[0] + d1[1]*d2[1] + d1[2]*d2[2]);
625 Float_t dist0 = track0->GetLinearD(0,0);
626 Float_t dist1 = track1->GetLinearD(0,0);
628 // conservative cuts - convergence to be guarantied
629 // applying before track propagation
630 if (TMath::Abs(TMath::Abs(dist0)-TMath::Abs(dist1))>fCutMaxD) continue; // distance to the 0,0
631 if (TMath::Abs(dir)<TMath::Abs(fCutMinDir)) continue; // direction vector product
632 Float_t bz = AliTracker::GetBz();
633 Float_t dvertex0[2]; //distance to 0,0
634 Float_t dvertex1[2]; //distance to 0,0
635 track0->GetDZ(0,0,0,bz,dvertex0);
636 track1->GetDZ(0,0,0,bz,dvertex1);
637 if (TMath::Abs(dvertex0[1])>250) continue;
638 if (TMath::Abs(dvertex1[1])>250) continue;
642 Float_t dmax = TMath::Max(TMath::Abs(dist0),TMath::Abs(dist1));
643 AliExternalTrackParam param0(*track0);
644 AliExternalTrackParam param1(*track1);
646 // Propagate using Magnetic field and correct fo material budget
648 AliTracker::PropagateTrackTo(¶m0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
649 AliTracker::PropagateTrackTo(¶m1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
651 // Propagate rest to the 0,0 DCA - z should be ignored
654 param0.PropagateToDCA(&vtx,bz,1000);
656 param1.PropagateToDCA(&vtx,bz,1000);
657 param0.GetDZ(0,0,0,bz,dvertex0);
658 param1.GetDZ(0,0,0,bz,dvertex1);
663 Bool_t isPair = IsPair(¶m0,¶m1);
664 Bool_t isCross = IsCross(track0, track1);
665 Bool_t isSame = IsSame(track0, track1);
667 THnSparse* hist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
668 TString shortName=hist->ClassName();
669 shortName+="_MEAN_VDRIFT_COSMICS_";
673 if((isSame) || (isCross && isPair)){
674 if (track0->GetTPCNcls()+ track1->GetTPCNcls()> 80) {
675 fDz = param0.GetZ() - param1.GetZ();
676 Double_t sign=(nA0>nA1)? 1:-1;
678 TTimeStamp tstamp(fTime);
679 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
680 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
681 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
682 THnSparse* curHist=NULL;
686 name+=event->GetFiredTriggerClasses();
688 curHist=(THnSparseF*)fArrayDz->FindObject(name);
690 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
691 fArrayDz->AddLast(curHist);
693 // curHist=(THnSparseF*)(fMapDz->GetValue(event->GetFiredTriggerClasses()));
695 // curHist=new THnSparseF(event->GetFiredTriggerClasses(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
696 // fMapDz->Add(new TObjString(event->GetFiredTriggerClasses()),curHist);
698 curHist->Fill(vecDrift);
703 curHist=(THnSparseF*)fArrayDz->FindObject(name);
705 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
706 fArrayDz->AddLast(curHist);
708 // curHist=(THnSparseF*)(fMapDz->GetValue("all"));
710 // curHist=new THnSparseF("all","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
711 // fMapDz->Add(new TObjString("all"),curHist);
713 curHist->Fill(vecDrift);
716 TTreeSRedirector *cstream = GetDebugStreamer();
719 (*cstream)<<"trackInfo"<<
730 "isCross="<<isCross<<
738 } // end 2nd order loop
739 } // end 1st order loop
742 TTreeSRedirector *cstream = GetDebugStreamer();
744 (*cstream)<<"timeInfo"<<
745 "run="<<fRun<< // run number
746 "event="<<fEvent<< // event number
747 "time="<<fTime<< // time stamp of event
748 "trigger="<<fTrigger<< // trigger
749 "mag="<<fMagF<< // magnetic field
750 // Environment values
752 // accumulated values
754 "fDz="<<fDz<< //! current delta z
755 "trigger="<<event->GetFiredTriggerClasses()<<
759 if (GetDebugLevel()>20) printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
762 void AliTPCcalibTime::ProcessBeam(const AliESDEvent *const event){
764 // Process beam data - calculates vartex
765 // from A side and C side
766 // Histogram the differences
768 const Int_t kMinClusters =80;
769 const Int_t kMinTracks =2; // minimal number of tracks to define the vertex
770 const Int_t kMinTracksVertex=30; // minimal number of tracks to define the cumulative vertex
771 const Double_t kMaxTgl =1.2; // maximal Tgl (z angle)
772 const Double_t kMinPt =0.2; // minimal pt
773 const Double_t kMaxD0 =5.; // cut on distance to the primary vertex first guess
774 const Double_t kMaxZ0 =20;
775 const Double_t kMaxD =2.5; // cut on distance to the primary vertex
776 const Double_t kMaxZ =4; // maximal z distance between tracks form the same side
777 const Double_t kMaxChi2 =15; // maximal chi2 of the TPCvertex
778 const Double_t kCumulCovarXY=0.003; //increase the error of cumul vertex 30 microns profile
779 const Double_t kCumulCovarZ=250.; //increase the error of cumul vertex
780 const Double_t kMaxDvertex = 1.0; // cut to accept the vertex;
783 const Int_t kBuffSize=100;
784 static Double_t deltaZ[kBuffSize]={0};
785 static Int_t counterZ=0;
786 static AliKFVertex cumulVertexA, cumulVertexC, cumulVertexAC; // cumulative vertex
787 AliKFVertex vertexA, vertexC;
789 Float_t dca0[2]={0,0};
790 Double_t dcaVertex[2]={0,0};
791 Int_t ntracks=event->GetNumberOfTracks();
792 if (ntracks==0) return;
793 if (ntracks > fCutTracks) return;
795 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
797 // Divide tracks to A and C side tracks - using the cluster indexes
798 TObjArray tracksA(ntracks);
799 TObjArray tracksC(ntracks);
801 AliESDVertex *vertexSPD = (AliESDVertex *)event->GetPrimaryVertexSPD();
802 AliESDVertex *vertex = (AliESDVertex *)event->GetPrimaryVertex();
803 AliESDVertex *vertexTracks = (AliESDVertex *)event->GetPrimaryVertexTracks();
804 Double_t vertexZA[10000], vertexZC[10000];
809 for (Int_t itrack=0;itrack<ntracks;itrack++) {
810 AliESDtrack *track = event->GetTrack(itrack);
811 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(itrack);
812 if (!friendTrack) continue;
813 if (TMath::Abs(track->GetTgl())>kMaxTgl) continue;
814 if (TMath::Abs(track->Pt())<kMinPt) continue;
815 const AliExternalTrackParam * trackIn = track->GetInnerParam();
816 TObject *calibObject=0;
817 AliTPCseed *seed = 0;
819 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
821 for (Int_t irow=159;irow>0;irow--) {
822 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
824 if ((cl->GetDetector()%36)<18) nA++;
825 if ((cl->GetDetector()%36)>=18) nC++;
827 if ((nA>kMinClusters || nC>kMinClusters) && (nA*nC==0) ){
828 track->GetImpactParameters(dca0[0],dca0[1]);
829 if (TMath::Abs(dca0[0])>kMaxD0) continue;
830 if (TMath::Abs(dca0[1])>kMaxZ0) continue;
831 AliExternalTrackParam pTPCvertex(*trackIn);
832 if (!AliTracker::PropagateTrackToBxByBz(&pTPCvertex,4.+4.*TMath::Abs(dca0[0]),0.1,2,kTRUE)) continue;
833 pTPCvertex.PropagateToDCA(vertex,AliTracker::GetBz(), kMaxD, dcaVertex,0);
834 if (TMath::Abs(dcaVertex[0])>kMaxD) continue;
835 if (nA>kMinClusters &&nC==0) { tracksA.AddLast(pTPCvertex.Clone()); vertexZA[ntracksA++] = pTPCvertex.GetZ();}
836 if (nC>kMinClusters &&nA==0) {tracksC.AddLast(pTPCvertex.Clone()); vertexZC[ntracksC++] = pTPCvertex.GetZ();}
840 Double_t medianZA=TMath::Median(ntracksA, vertexZA); // tracks median
841 Double_t medianZC=TMath::Median(ntracksC, vertexZC); // tracks median
843 ntracksA= tracksA.GetEntriesFast();
844 ntracksC= tracksC.GetEntriesFast();
845 if (ntracksA>kMinTracks && ntracksC>kMinTracks){
846 deltaZ[counterZ%kBuffSize]=medianZA-medianZC;
848 Double_t medianDelta=(counterZ>=kBuffSize)? TMath::Median(kBuffSize, deltaZ): TMath::Median(counterZ, deltaZ);
849 if (TMath::Abs(medianDelta-(medianZA-medianZC))>kMaxZ) flags+=16;
850 // increse the error of cumulative vertex at the beginning of event
851 cumulVertexA.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
852 cumulVertexA.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
853 cumulVertexA.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
854 cumulVertexC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
855 cumulVertexC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
856 cumulVertexC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
857 cumulVertexAC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
858 cumulVertexAC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
859 cumulVertexAC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
861 for (Int_t iA=0; iA<ntracksA; iA++){
862 if (flags!=0) continue;
863 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksA.At(iA);
864 if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
865 AliKFParticle part(*aliTrack,211);
868 for (Int_t iC=0; iC<ntracksC; iC++){
869 if (flags!=0) continue;
870 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksC.At(iC);
871 if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
872 AliKFParticle part(*aliTrack,211);
876 if (vertexA.GetNDF()<kMinTracks) flags+=32;
877 if (vertexC.GetNDF()<kMinTracks) flags+=32;
878 if (TMath::Abs(vertexA.Z()-medianZA)>kMaxZ) flags+=1; //apply cuts
879 if (TMath::Abs(vertexC.Z()-medianZC)>kMaxZ) flags+=2;
880 if (TMath::Abs(vertexA.GetChi2()/vertexA.GetNDF()+vertexC.GetChi2()/vertexC.GetNDF())> kMaxChi2) flags+=4;
883 for (Int_t iA=0; iA<ntracksA; iA++){
884 if (flags!=0) continue;
885 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksA.At(iA);
886 if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
887 AliKFParticle part(*aliTrack,211);
891 for (Int_t iC=0; iC<ntracksC; iC++){
892 if (flags!=0) continue;
893 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksC.At(iC);
894 if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
895 AliKFParticle part(*aliTrack,211);
900 if (TMath::Abs(cumulVertexA.X()-vertexA.X())>kMaxDvertex) flags+=64;
901 if (TMath::Abs(cumulVertexA.Y()-vertexA.Y())>kMaxDvertex) flags+=64;
902 if (TMath::Abs(cumulVertexA.Z()-vertexA.Z())>kMaxDvertex) flags+=64;
904 if (TMath::Abs(cumulVertexC.X()-vertexC.X())>kMaxDvertex) flags+=64;
905 if (TMath::Abs(cumulVertexC.Y()-vertexC.Y())>kMaxDvertex) flags+=64;
906 if (TMath::Abs(cumulVertexC.Z()-vertexC.Z())>kMaxDvertex) flags+=64;
909 if ( flags==0 && cumulVertexC.GetNDF()>kMinTracksVertex&&cumulVertexA.GetNDF()>kMinTracksVertex){
910 Double_t cont[2]={0,fTime};
912 cont[0]= cumulVertexA.X();
913 fTPCVertex[0]->Fill(cont);
914 cont[0]= cumulVertexC.X();
915 fTPCVertex[1]->Fill(cont);
916 cont[0]= 0.5*(cumulVertexA.X()-cumulVertexC.X());
917 fTPCVertex[2]->Fill(cont);
918 cont[0]= 0.5*(cumulVertexA.X()+cumulVertexC.X())-vertexSPD->GetX();
919 fTPCVertex[3]->Fill(cont);
921 cont[0]= cumulVertexA.Y();
922 fTPCVertex[4]->Fill(cont);
923 cont[0]= cumulVertexC.Y();
924 fTPCVertex[5]->Fill(cont);
925 cont[0]= 0.5*(cumulVertexA.Y()-cumulVertexC.Y());
926 fTPCVertex[6]->Fill(cont);
927 cont[0]= 0.5*(cumulVertexA.Y()+cumulVertexC.Y())-vertexSPD->GetY();
928 fTPCVertex[7]->Fill(cont);
931 cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z());
932 fTPCVertex[8]->Fill(cont);
933 cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
934 fTPCVertex[9]->Fill(cont);
935 cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
936 fTPCVertex[10]->Fill(cont);
937 cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z())-vertexSPD->GetZ();
938 fTPCVertex[11]->Fill(cont);
940 Double_t correl[2]={0,0};
942 correl[0]=cumulVertexC.Z();
943 correl[1]=cumulVertexA.Z();
944 fTPCVertexCorrelation[0]->Fill(correl); // fill A side :TPC
945 correl[0]=cumulVertexA.Z();
946 correl[1]=cumulVertexC.Z();
947 fTPCVertexCorrelation[1]->Fill(correl); // fill C side :TPC
949 correl[0]=vertexSPD->GetZ();
950 correl[1]=cumulVertexA.Z()-correl[0];
951 fTPCVertexCorrelation[2]->Fill(correl); // fill A side :ITS
952 correl[1]=cumulVertexC.Z()-correl[0];
953 fTPCVertexCorrelation[3]->Fill(correl); // fill C side :ITS
954 correl[1]=0.5*(cumulVertexA.Z()+cumulVertexC.Z())-correl[0];
955 fTPCVertexCorrelation[4]->Fill(correl); // fill C side :ITS
958 TTreeSRedirector *cstream = GetDebugStreamer();
961 TCut cutChi2= "sqrt(vA.fChi2/vA.fNDF+vC.fChi2/vC.fNDF)<10"; // chi2 Cut e.g 10
962 TCut cutXY= "sqrt((vA.fP[0]-vC.fP[0])^2+(vA.fP[0]-vC.fP[1])^2)<5"; // vertex Cut
963 TCut cutZ= "abs(vA.fP[2]-mZA)<3&&abs(vC.fP[2]-mZC)<5"; // vertex Cut
964 tree->Draw("sqrt(vA.fChi2/vA.fNDF)","sqrt(vA.fChi2/vA.fNDF)<100","")
969 (*cstream)<<"vertexTPC"<<
970 "flags="<<flags<< // rejection flags
971 "vSPD.="<<vertexSPD<< // SPD vertex
972 "vT.="<<vertexTracks<< // track vertex
973 "v.="<<vertex<< // esd vertex
974 "mZA="<<medianZA<< // median Z position at vertex A side
975 "mZC="<<medianZC<< // median Z position at vertex C side
976 "mDelta="<<medianDelta<< // median delta A side -C side
977 "counter="<<counterZ<< // counter Z
979 "vA.="<<&vertexA<< // vertex A side
980 "vC.="<<&vertexC<< // vertex C side
981 "cvA.="<<&cumulVertexA<< // cumulative vertex A side
982 "cvC.="<<&cumulVertexC<< // cumulative vertex C side
983 "cvAC.="<<&cumulVertexAC<< // cumulative vertex A+C side
984 "nA="<<ntracksA<< // contributors
985 "nC="<<ntracksC<< // contributors
993 void AliTPCcalibTime::Analyze(){
995 // Special macro to analyze result of calibration and extract calibration entries
996 // Not yet ported to the Analyze function yet
1000 THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name) const
1003 // Get histogram for given trigger mask
1005 TIterator* iterator = fArrayDz->MakeIterator();
1007 TString newName=name;
1009 THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1010 THnSparse* addHist=NULL;
1011 while((addHist=(THnSparseF*)iterator->Next())){
1012 // if(!addHist) continue;
1013 TString histName=addHist->GetName();
1014 if(!histName.Contains(newName)) continue;
1016 newHist->Add(addHist);
1021 TObjArray* AliTPCcalibTime::GetHistoDrift() const
1024 // return array of histograms
1029 TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
1031 // Make a drift velocity (delta Z) graph
1033 THnSparse* histoDrift=GetHistoDrift(name);
1034 TGraphErrors* graphDrift=NULL;
1036 graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
1037 TString end=histoDrift->GetName();
1038 Int_t pos=end.Index("_");
1039 end=end(pos,end.Capacity()-pos);
1040 TString graphName=graphDrift->ClassName();
1042 graphName.ToUpper();
1043 graphDrift->SetName(graphName);
1048 TObjArray* AliTPCcalibTime::GetGraphDrift(){
1050 // make a array of drift graphs
1052 TObjArray* arrayGraphDrift=new TObjArray();
1053 TIterator* iterator=fArrayDz->MakeIterator();
1055 THnSparse* addHist=NULL;
1056 while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
1057 return arrayGraphDrift;
1060 AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
1062 // Make a fit AliSplinefit of drift velocity
1064 TGraph* graphDrift=GetGraphDrift(name);
1065 AliSplineFit* fitDrift=NULL;
1066 if(graphDrift && graphDrift->GetN()){
1067 fitDrift=new AliSplineFit();
1068 fitDrift->SetGraph(graphDrift);
1069 fitDrift->SetMinPoints(graphDrift->GetN()+1);
1070 fitDrift->InitKnots(graphDrift,2,0,0.001);
1071 fitDrift->SplineFit(0);
1072 TString end=graphDrift->GetName();
1073 Int_t pos=end.Index("_");
1074 end=end(pos,end.Capacity()-pos);
1075 TString fitName=fitDrift->ClassName();
1078 //fitDrift->SetName(fitName);
1086 Long64_t AliTPCcalibTime::Merge(TCollection *const li) {
1088 // Object specific merging procedure
1090 TIterator* iter = li->MakeIterator();
1091 AliTPCcalibTime* cal = 0;
1093 while ((cal = (AliTPCcalibTime*)iter->Next())) {
1094 if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
1095 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
1098 for (Int_t imeas=0; imeas<3; imeas++){
1099 if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
1100 fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
1101 fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
1105 if (fTPCVertexCorrelation[0] && cal->fTPCVertexCorrelation[0]){
1106 for (Int_t imeas=0; imeas<5; imeas++){
1107 if (fTPCVertexCorrelation[imeas] && cal->fTPCVertexCorrelation[imeas]) fTPCVertexCorrelation[imeas]->Add(cal->fTPCVertexCorrelation[imeas]);
1111 if (fTPCVertex[0] && cal->fTPCVertex[0])
1112 for (Int_t imeas=0; imeas<12; imeas++){
1113 if (fTPCVertex[imeas] && cal->fTPCVertex[imeas]) fTPCVertex[imeas]->Add(cal->fTPCVertex[imeas]);
1116 if (fMemoryMode>0) for (Int_t imeas=0; imeas<5; imeas++){
1118 if ( GetResHistoTPCCE(imeas) && cal->GetResHistoTPCCE(imeas)){
1119 if ((cal->GetResHistoTPCCE(imeas)->GetEntries()+GetResHistoTPCCE(imeas)->GetEntries()) < fgResHistoMergeCut)
1120 fResHistoTPCCE[imeas]->Add(cal->fResHistoTPCCE[imeas]);
1122 fResHistoTPCCE[imeas]=(THnSparse*)cal->fResHistoTPCCE[imeas]->Clone();
1126 if ((fMemoryMode>0) &&cal->GetResHistoTPCITS(imeas) && cal->GetResHistoTPCITS(imeas)){
1127 if (fMemoryMode>1 || (imeas%2)==1)
1129 if (fResHistoTPCITS[imeas]->GetEntries()+(cal->fResHistoTPCITS[imeas])->GetEntries() < fgResHistoMergeCut)
1130 fResHistoTPCITS[imeas]->Add(cal->fResHistoTPCITS[imeas]);
1134 if (fResHistoTPCvertex[imeas]->GetEntries()+(cal->fResHistoTPCvertex[imeas])->GetEntries() < fgResHistoMergeCut)
1135 fResHistoTPCvertex[imeas]->Add(cal->fResHistoTPCvertex[imeas]);
1139 if ((fMemoryMode>1) && cal->fResHistoTPCTRD[imeas]){
1140 if (fResHistoTPCTRD[imeas] && cal->fResHistoTPCTRD[imeas])
1142 if (fResHistoTPCTRD[imeas]->GetEntries()+(cal->fResHistoTPCTRD[imeas])->GetEntries() < fgResHistoMergeCut)
1143 fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
1146 fResHistoTPCTRD[imeas]=(THnSparse*)cal->fResHistoTPCTRD[imeas]->Clone();
1149 if ((fMemoryMode>1) && cal->fResHistoTPCTOF[imeas]){
1150 if (fResHistoTPCTOF[imeas] && cal->fResHistoTPCTOF[imeas])
1152 if (fResHistoTPCTOF[imeas]->GetEntries()+(cal->fResHistoTPCTOF[imeas])->GetEntries() < fgResHistoMergeCut)
1153 fResHistoTPCTOF[imeas]->Add(cal->fResHistoTPCTOF[imeas]);
1156 fResHistoTPCTOF[imeas]=(THnSparse*)cal->fResHistoTPCTOF[imeas]->Clone();
1159 if (cal->fArrayLaserA){
1160 fArrayLaserA->Expand(fArrayLaserA->GetEntriesFast()+cal->fArrayLaserA->GetEntriesFast());
1161 fArrayLaserC->Expand(fArrayLaserC->GetEntriesFast()+cal->fArrayLaserC->GetEntriesFast());
1162 for (Int_t ical=0; ical<cal->fArrayLaserA->GetEntriesFast(); ical++){
1163 if (cal->fArrayLaserA->UncheckedAt(ical)) fArrayLaserA->AddLast(cal->fArrayLaserA->UncheckedAt(ical)->Clone());
1164 if (cal->fArrayLaserC->UncheckedAt(ical)) fArrayLaserC->AddLast(cal->fArrayLaserC->UncheckedAt(ical)->Clone());
1169 // TObjArray* addArray=cal->GetHistoDrift();
1170 // if(!addArray) return 0;
1171 // TIterator* iterator = addArray->MakeIterator();
1172 // iterator->Reset();
1173 // THnSparse* addHist=NULL;
1174 // if ((fMemoryMode>1)) while((addHist=(THnSparseF*)iterator->Next())){
1175 // // if(!addHist) continue;
1176 // addHist->Print();
1177 // THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
1179 // localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1180 // fArrayDz->AddLast(localHist);
1182 // localHist->Add(addHist);
1185 for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
1189 const Int_t kMinUpdates=10;
1190 const Float_t kMaxOut=0.1;
1191 for (Int_t itype=0; itype<3; itype++){
1196 if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
1197 if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
1198 if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
1199 if (!arr1) continue;
1200 if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
1201 if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
1202 arr0->Expand(arr1->GetEntriesFast());
1204 for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
1205 AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
1206 AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
1207 if (!kalman1) continue;
1208 if (kalman1->GetNUpdates()<kMinUpdates) continue;
1209 if (kalman1->GetNOutliers()>(kalman1->GetNUpdates()*kMaxOut)) continue;
1210 if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
1211 kalman0->SetRejectOutliers(kFALSE);
1212 kalman0->Merge(kalman1);
1221 Bool_t AliTPCcalibTime::IsPair(const AliExternalTrackParam *tr0, const AliExternalTrackParam *tr1){
1223 // 0. Same direction - OPOSITE - cutDir +cutT
1224 TCut cutDir("cutDir","dir<-0.99")
1226 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
1229 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
1231 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
1234 const Double_t *p0 = tr0->GetParameter();
1235 const Double_t *p1 = tr1->GetParameter();
1236 fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
1237 fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
1238 fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1239 fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
1240 fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
1242 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
1243 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
1244 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
1245 Double_t d0[3], d1[3];
1246 tr0->GetDirection(d0);
1247 tr1->GetDirection(d1);
1248 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
1250 fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
1251 fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
1252 fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1253 fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
1254 fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
1258 Bool_t AliTPCcalibTime::IsCross(const AliESDtrack *const tr0, const AliESDtrack *const tr1){
1260 // check if the cosmic pair of tracks crossed A/C side
1262 Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
1263 if (result==kFALSE) return result;
1268 Bool_t AliTPCcalibTime::IsSame(const AliESDtrack *const tr0, const AliESDtrack *const tr1){
1270 // track crossing the CE
1271 // 0. minimal number of clusters
1272 // 1. Same sector +-1
1273 // 2. Inner and outer track param on opposite side
1274 // 3. Outer and inner track parameter close each to other
1276 Bool_t result=kTRUE;
1278 // inner and outer on opposite sides in z
1280 const Int_t knclCut0 = 30;
1281 const Double_t kalphaCut = 0.4;
1283 // 0. minimal number of clusters
1285 if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
1286 if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
1288 // 1. alpha cut - sector+-1
1290 if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
1294 if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
1295 if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
1296 if (result==kFALSE){
1301 const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
1302 const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
1303 const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1304 const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1306 if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD) result&=kFALSE;
1307 if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1308 if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1309 if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1310 if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD) result&=kFALSE;
1311 if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1312 if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1313 if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1315 result=kTRUE; // just to put break point here
1321 void AliTPCcalibTime::ProcessSame(const AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
1323 // Process TPC tracks crossing CE
1325 // 0. Select only track crossing the CE
1326 // 1. Cut on the track length
1327 // 2. Refit the the track on A and C side separatelly
1328 // 3. Fill time histograms
1329 const Int_t kMinNcl=100;
1330 const Int_t kMinNclS=25; // minimul number of clusters on the sides
1331 const Double_t pimass=TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1332 const Double_t kMaxDy=1; // maximal distance in y
1333 const Double_t kMaxDsnp=0.05; // maximal distance in snp
1334 const Double_t kMaxDtheta=0.05; // maximal distance in theta
1336 if (!friendTrack->GetTPCOut()) return;
1338 // 0. Select only track crossing the CE
1340 if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1342 // 1. cut on track length
1344 if (track->GetTPCNcls()<kMinNcl) return;
1346 // 2. Refit track sepparatel on A and C side
1348 TObject *calibObject;
1349 AliTPCseed *seed = 0;
1350 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1351 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1355 AliExternalTrackParam trackIn(*track->GetInnerParam());
1356 AliExternalTrackParam trackOut(*track->GetOuterParam());
1357 Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1358 Double_t xyz[3]={0,0.,0.0};
1360 Int_t nclIn=0,nclOut=0;
1361 trackIn.ResetCovariance(1000.);
1362 trackOut.ResetCovariance(1000.);
1367 for (Int_t irow=0;irow<159;irow++) {
1368 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1370 if (cl->GetX()<80) continue;
1372 if (cl->GetDetector()%36<18) sideIn=1;
1373 if (cl->GetDetector()%36>=18) sideIn=-1;
1375 if (sideIn== -1 && (cl->GetDetector()%36)<18) break;
1376 if (sideIn== 1 &&(cl->GetDetector()%36)>=18) break;
1377 Int_t sector = cl->GetDetector();
1378 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1379 if (TMath::Abs(dalpha)>0.01){
1380 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1382 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1383 trackIn.GetXYZ(xyz);
1384 bz = AliTracker::GetBz(xyz);
1385 AliTracker::PropagateTrackToBxByBz(&trackIn,r[0],pimass,1.,kFALSE);
1386 if (!trackIn.PropagateTo(r[0],bz)) break;
1388 trackIn.Update(&r[1],cov);
1394 for (Int_t irow=159;irow>0;irow--) {
1395 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1397 if (cl->GetX()<80) continue;
1399 if (cl->GetDetector()%36<18) sideOut=1;
1400 if (cl->GetDetector()%36>=18) sideOut=-1;
1401 if (sideIn==sideOut) break;
1403 if (sideOut== -1 && (cl->GetDetector()%36)<18) break;
1404 if (sideOut== 1 &&(cl->GetDetector()%36)>=18) break;
1406 Int_t sector = cl->GetDetector();
1407 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1408 if (TMath::Abs(dalpha)>0.01){
1409 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1411 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1412 trackOut.GetXYZ(xyz);
1413 bz = AliTracker::GetBz(xyz);
1414 AliTracker::PropagateTrackToBxByBz(&trackOut,r[0],pimass,1.,kFALSE);
1415 if (!trackOut.PropagateTo(r[0],bz)) break;
1417 trackOut.Update(&r[1],cov);
1419 trackOut.Rotate(trackIn.GetAlpha());
1420 Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1421 trackIn.PropagateTo(meanX,bz);
1422 trackOut.PropagateTo(meanX,bz);
1423 if (TMath::Abs(trackIn.GetY()-trackOut.GetY())>kMaxDy) return;
1424 if (TMath::Abs(trackIn.GetSnp()-trackOut.GetSnp())>kMaxDsnp) return;
1425 if (TMath::Abs(trackIn.GetTgl()-trackOut.GetTgl())>kMaxDtheta) return;
1426 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1427 FillResHistoTPCCE(&trackIn,&trackOut);
1429 TTreeSRedirector *cstream = GetDebugStreamer();
1432 trackIn.GetXYZ(gxyz.GetMatrixArray());
1433 TTimeStamp tstamp(fTime);
1434 (*cstream)<<"tpctpc"<<
1435 "run="<<fRun<< // run number
1436 "event="<<fEvent<< // event number
1437 "time="<<fTime<< // time stamp of event
1438 "trigger="<<fTrigger<< // trigger
1439 "mag="<<fMagF<< // magnetic field
1441 "sideIn="<<sideIn<< // side at inner part
1442 "sideOut="<<sideOut<< // side at puter part
1443 "xyz.="<<&gxyz<< // global position
1444 "tIn.="<<&trackIn<< // refitterd track in
1445 "tOut.="<<&trackOut<< // refitter track out
1446 "nclIn="<<nclIn<< //
1447 "nclOut="<<nclOut<< //
1451 // 3. Fill time histograms
1452 // Debug stremaer expression
1453 // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1454 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1455 fDz = trackOut.GetZ()-trackIn.GetZ();
1456 if (trackOut.GetTgl()<0) fDz*=-1.;
1457 TTimeStamp tstamp(fTime);
1458 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1459 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1460 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1462 // fill histograms per trigger class and itegrated
1464 THnSparse* curHist=NULL;
1465 for (Int_t itype=0; itype<2; itype++){
1466 TString name="MEAN_VDRIFT_CROSS_";
1468 name+=event->GetFiredTriggerClasses();
1473 curHist=(THnSparseF*)fArrayDz->FindObject(name);
1475 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1476 fArrayDz->AddLast(curHist);
1478 curHist->Fill(vecDrift);
1484 void AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
1486 // Process track - Update TPC-ITS alignment
1488 // 0. Apply standartd cuts
1489 // 1. Recalucluate the current statistic median/RMS
1490 // 2. Apply median+-rms cut
1491 // 3. Update kalman filter
1493 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1494 const Int_t kMinITS = 3; // minimal number of ITS cluster
1495 const Double_t kMinZ = 10; // maximal dz distance
1496 const Double_t kMaxDy = 2.; // maximal dy distance
1497 const Double_t kMaxAngle= 0.07; // maximal angular distance
1498 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1499 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1500 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
1501 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1502 const Double_t kOutCut = 3.0; // outlyer cut in AliRelAlgnmentKalman
1503 const Double_t kMinPt = 0.3; // minimal pt
1504 const Double_t kMax1Pt=0.5; //maximal 1/pt distance
1505 const Int_t kN=50; // deepnes of history
1506 static Int_t kglast=0;
1507 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1509 // 0. Apply standard cuts
1511 Int_t dummycl[1000];
1512 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1513 if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
1514 if (!track->GetInnerParam()) return;
1515 if (!track->GetOuterParam()) return;
1516 if (track->GetInnerParam()->Pt()<kMinPt) return;
1517 // exclude crossing track
1518 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1519 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ/3.) return;
1520 if (track->GetInnerParam()->GetX()>90) return;
1522 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
1524 AliExternalTrackParam pITS; // ITS standalone if possible
1525 AliExternalTrackParam pITS2; //TPC-ITS track
1526 if (friendTrack->GetITSOut()){
1527 pITS2=(*(friendTrack->GetITSOut())); //TPC-ITS track - snapshot ITS out
1528 pITS2.Rotate(pTPC.GetAlpha());
1529 AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1532 AliESDfriendTrack *itsfriendTrack=0;
1534 // try to find standalone ITS track corresponing to the TPC if possible
1536 Bool_t hasAlone=kFALSE;
1537 Int_t ntracks=event->GetNumberOfTracks();
1538 for (Int_t i=0; i<ntracks; i++){
1539 AliESDtrack * trackITS = event->GetTrack(i);
1540 if (!trackITS) continue;
1541 if (trackITS->GetITSclusters(dummycl)<kMinITS) continue; // minimal amount of clusters
1542 itsfriendTrack = esdFriend->GetTrack(i);
1543 if (!itsfriendTrack) continue;
1544 if (!itsfriendTrack->GetITSOut()) continue;
1546 if (TMath::Abs(pTPC.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1547 if (TMath::Abs(pTPC.GetSigned1Pt()-itsfriendTrack->GetITSOut()->GetSigned1Pt())> kMax1Pt) continue;
1548 pITS=(*(itsfriendTrack->GetITSOut()));
1550 pITS.Rotate(pTPC.GetAlpha());
1551 AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
1552 if (TMath::Abs(pTPC.GetY()-pITS.GetY())> kMaxDy) continue;
1553 if (TMath::Abs(pTPC.GetSnp()-pITS.GetSnp())> kMaxAngle) continue;
1557 if (track->GetITSclusters(dummycl)<kMinITS) return;
1558 pITS=pITS2; // use combined track if it has ITS
1561 if (TMath::Abs(pITS.GetY()-pTPC.GetY()) >kMaxDy) return;
1562 if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1563 if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1565 // 1. Update median and RMS info
1567 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1568 TVectorD vecDeltaN(5);
1569 Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1571 for (Int_t i=0;i<4;i++){
1572 vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1573 kgdP[i][kglast%kN]=vecDelta[i];
1576 Int_t entries=(kglast<kN)?kglast:kN;
1577 for (Int_t i=0;i<4;i++){
1578 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1579 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1582 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1583 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1587 // 2. Apply median+-rms cut
1589 if (kglast<3) return; //median and RMS to be defined
1590 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1592 // 3. Update alignment
1594 Int_t htime = (fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time bins number
1595 if (fAlignITSTPC->GetEntriesFast()<htime){
1596 fAlignITSTPC->Expand(htime*2+20);
1598 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1600 // make Alignment object if doesn't exist
1601 align=new AliRelAlignerKalman();
1602 align->SetRunNumber(fRun);
1603 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1604 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1605 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1606 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1607 // align->SetRejectOutliers(kFALSE);
1608 align->SetRejectOutliers(kTRUE);
1609 align->SetRejectOutliersSigma2Median(kTRUE);
1611 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1612 align->SetMagField(fMagF);
1613 fAlignITSTPC->AddAt(align,htime);
1615 align->AddTrackParams(&pITS,&pTPC);
1616 Double_t averageTime = fTime;
1617 if (align->GetTimeStamp()>0&&align->GetNUpdates()>0){
1618 averageTime=((Double_t(align->GetTimeStamp())*Double_t(align->GetNUpdates())+Double_t(fTime)))/(Double_t(align->GetNUpdates())+1.);
1620 align->SetTimeStamp(Int_t(averageTime));
1622 align->SetRunNumber(fRun );
1623 Float_t dca[2],cov[3];
1624 track->GetImpactParameters(dca,cov);
1625 if (TMath::Abs(dca[0])<kMaxDy){
1626 FillResHistoTPCITS(&pTPC,&pITS);
1627 FillResHistoTPC(track);
1630 Int_t nupdates=align->GetNUpdates();
1631 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates+1));
1632 // align->SetRejectOutliers(kFALSE);
1633 align->SetRejectOutliers(kTRUE);
1634 align->SetRejectOutliersSigma2Median(kTRUE);
1636 TTreeSRedirector *cstream = GetDebugStreamer();
1637 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1638 TVectorD gpTPC(3), gdTPC(3);
1639 TVectorD gpITS(3), gdITS(3);
1640 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1641 pTPC.GetDirection(gdTPC.GetMatrixArray());
1642 pITS.GetXYZ(gpITS.GetMatrixArray());
1643 pITS.GetDirection(gdITS.GetMatrixArray());
1644 (*cstream)<<"itstpc"<<
1645 "run="<<fRun<< // run number
1646 "event="<<fEvent<< // event number
1647 "time="<<fTime<< // time stamp of event
1648 "trigger="<<fTrigger<< // trigger
1649 "mag="<<fMagF<< // magnetic field
1651 "hasAlone="<<hasAlone<< // has ITS standalone ?
1652 "track.="<<track<< // track info
1653 "nmed="<<kglast<< // number of entries to define median and RMS
1654 "vMed.="<<&vecMedian<< // median of deltas
1655 "vRMS.="<<&vecRMS<< // rms of deltas
1656 "vDelta.="<<&vecDelta<< // delta in respect to median
1657 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1658 "a.="<<align<< // current alignment
1659 "pITS.="<<&pITS<< // track param ITS
1660 "pITS2.="<<&pITS2<< // track param ITS+TPC
1661 "pTPC.="<<&pTPC<< // track param TPC
1662 "gpTPC.="<<&gpTPC<< // global position TPC
1663 "gdTPC.="<<&gdTPC<< // global direction TPC
1664 "gpITS.="<<&gpITS<< // global position ITS
1665 "gdITS.="<<&gdITS<< // global position ITS
1673 void AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack){
1675 // Process track - Update TPC-TRD alignment
1677 // 0. Apply standartd cuts
1678 // 1. Recalucluate the current statistic median/RMS
1679 // 2. Apply median+-rms cut
1680 // 3. Update kalman filter
1682 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1683 const Int_t kMinTRD = 60; // minimal number of TRD cluster
1684 // const Double_t kMinZ = 20; // maximal dz distance
1685 const Double_t kMaxDy = 5.; // maximal dy distance
1686 const Double_t kMaxAngle= 0.1; // maximal angular distance
1687 const Double_t kSigmaCut= 10; // maximal sigma distance to median
1688 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1689 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
1690 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1691 const Double_t kOutCut = 3.0; // outlyer cut in AliRelAlgnmentKalman
1692 const Double_t kRefX = 330; // reference X
1693 const Int_t kN=50; // deepnes of history
1694 static Int_t kglast=0;
1695 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1697 // 0. Apply standard cuts
1699 Int_t dummycl[1000];
1700 if (track->GetTRDclusters(dummycl)<kMinTRD) return; // minimal amount of clusters
1701 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1702 // if (!friendTrack->GetTRDIn()) return;
1703 // if (!track->IsOn(AliESDtrack::kTRDrefit)) return;
1704 if (!track->IsOn(AliESDtrack::kTRDout)) return;
1705 if (!track->GetInnerParam()) return;
1706 if (!friendTrack->GetTPCOut()) return;
1707 // exclude crossing track
1708 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1710 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1711 AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
1712 AliExternalTrackParam *pTRDtrack = 0;
1713 TObject *calibObject=0;
1714 for (Int_t l=0;(calibObject=((AliESDfriendTrack*)friendTrack)->GetCalibObject(l));++l) {
1715 if ((dynamic_cast< AliTPCseed*>(calibObject))) continue;
1716 if ((pTRDtrack=dynamic_cast< AliExternalTrackParam*>(calibObject))) break;
1718 if (!pTRDtrack) return;
1719 // AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
1720 AliExternalTrackParam pTRD(*(pTRDtrack));
1721 pTRD.Rotate(pTPC.GetAlpha());
1722 // pTRD.PropagateTo(pTPC.GetX(),fMagF);
1723 AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1725 ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.; // increas sys errors
1726 ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1728 if (TMath::Abs(pTRD.GetY()-pTPC.GetY()) >kMaxDy) return;
1729 if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1730 // if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1732 // 1. Update median and RMS info
1734 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1735 TVectorD vecDeltaN(5);
1736 Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1738 for (Int_t i=0;i<4;i++){
1739 vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1740 kgdP[i][kglast%kN]=vecDelta[i];
1743 Int_t entries=(kglast<kN)?kglast:kN;
1744 for (Int_t i=0;i<4;i++){
1745 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1747 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1750 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1751 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1755 // 2. Apply median+-rms cut
1757 if (kglast<3) return; //median and RMS to be defined
1758 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1760 // 3. Update alignment
1762 //Int_t htime = fTime/3600; //time in hours
1763 Int_t htime = (Int_t)(fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time in half hour
1764 if (fAlignTRDTPC->GetEntriesFast()<htime){
1765 fAlignTRDTPC->Expand(htime*2+20);
1767 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1769 // make Alignment object if doesn't exist
1770 align=new AliRelAlignerKalman();
1771 align->SetRunNumber(fRun);
1772 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1773 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1774 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1775 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1776 // align->SetRejectOutliers(kFALSE);
1777 align->SetRejectOutliers(kTRUE);
1778 align->SetRejectOutliersSigma2Median(kTRUE);
1780 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1781 align->SetMagField(fMagF);
1782 fAlignTRDTPC->AddAt(align,htime);
1784 align->AddTrackParams(&pTRD,&pTPC);
1785 //align->SetTimeStamp(fTime);
1786 Double_t averageTime = fTime;
1787 if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1788 averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1789 //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1791 align->SetTimeStamp((Int_t)averageTime);
1793 //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1795 align->SetRunNumber(fRun );
1796 Float_t dca[2],cov[3];
1797 track->GetImpactParameters(dca,cov);
1798 if (TMath::Abs(dca[0])<kMaxDy){
1799 FillResHistoTPCTRD(&pTPC,&pTRD); //only primaries
1802 Int_t nupdates=align->GetNUpdates();
1803 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates+1));
1804 // align->SetRejectOutliers(kFALSE);
1805 align->SetRejectOutliers(kTRUE);
1806 align->SetRejectOutliersSigma2Median(kTRUE);
1808 TTreeSRedirector *cstream = GetDebugStreamer();
1809 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1810 TVectorD gpTPC(3), gdTPC(3);
1811 TVectorD gpTRD(3), gdTRD(3);
1812 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1813 pTPC.GetDirection(gdTPC.GetMatrixArray());
1814 pTRD.GetXYZ(gpTRD.GetMatrixArray());
1815 pTRD.GetDirection(gdTRD.GetMatrixArray());
1816 (*cstream)<<"trdtpc"<<
1817 "run="<<fRun<< // run number
1818 "event="<<fEvent<< // event number
1819 "time="<<fTime<< // time stamp of event
1820 "trigger="<<fTrigger<< // trigger
1821 "mag="<<fMagF<< // magnetic field
1823 "nmed="<<kglast<< // number of entries to define median and RMS
1824 "vMed.="<<&vecMedian<< // median of deltas
1825 "vRMS.="<<&vecRMS<< // rms of deltas
1826 "vDelta.="<<&vecDelta<< // delta in respect to median
1827 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1828 "t.="<<track<< // ful track - find proper cuts
1829 "a.="<<align<< // current alignment
1830 "pTRD.="<<&pTRD<< // track param TRD
1831 "pTPC.="<<&pTPC<< // track param TPC
1832 "gpTPC.="<<&gpTPC<< // global position TPC
1833 "gdTPC.="<<&gdTPC<< // global direction TPC
1834 "gpTRD.="<<&gpTRD<< // global position TRD
1835 "gdTRD.="<<&gdTRD<< // global position TRD
1841 void AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack){
1844 // Process track - Update TPC-TOF alignment
1846 // -1. Make a TOF "track"
1847 // 0. Apply standartd cuts
1848 // 1. Recalucluate the current statistic median/RMS
1849 // 2. Apply median+-rms cut
1850 // 3. Update kalman filter
1852 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1853 // const Double_t kMinZ = 10; // maximal dz distance
1854 const Double_t kMaxDy = 5.; // maximal dy distance
1855 const Double_t kMaxAngle= 0.05; // maximal angular distance
1856 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1857 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1858 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
1859 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1861 const Double_t kOutCut = 3.0; // outlyer cut in AliRelAlgnmentKalman
1862 const Int_t kN=50; // deepnes of history
1863 static Int_t kglast=0;
1864 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1866 // -1. Make a TOF track-
1867 // Clusters are not in friends - use alingment points
1869 if (track->GetTOFsignal()<=0) return;
1870 if (!friendTrack->GetTPCOut()) return;
1871 if (!track->GetInnerParam()) return;
1872 if (!friendTrack->GetTPCOut()) return;
1873 const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1874 if (!points) return;
1875 AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
1876 AliExternalTrackParam pTOF(pTPC);
1877 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1878 Int_t npoints = points->GetNPoints();
1879 AliTrackPoint point;
1882 for (Int_t ipoint=0;ipoint<npoints;ipoint++){
1883 points->GetPoint(point,ipoint);
1886 Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1887 if (r<350) continue;
1888 if (r>400) continue;
1889 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1890 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);
1891 AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1892 pTPC.PropagateTo(lpoint.GetX(),fMagF);
1894 ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1895 ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
1896 ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1897 ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
1898 ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1899 ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
1902 if (naccept==0) return; // no tof match clusters
1904 // 0. Apply standard cuts
1906 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1907 // exclude crossing track
1908 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1910 if (TMath::Abs(pTOF.GetY()-pTPC.GetY()) >kMaxDy) return;
1911 if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1912 if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1914 // 1. Update median and RMS info
1916 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1917 TVectorD vecDeltaN(5);
1918 Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1920 for (Int_t i=0;i<4;i++){
1921 vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1922 kgdP[i][kglast%kN]=vecDelta[i];
1925 Int_t entries=(kglast<kN)?kglast:kN;
1927 for (Int_t i=0;i<4;i++){
1928 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1929 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1932 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1933 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1934 if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
1938 // 2. Apply median+-rms cut
1940 if (kglast<10) return; //median and RMS to be defined
1943 // 3. Update alignment
1945 //Int_t htime = fTime/3600; //time in hours
1946 Int_t htime = (Int_t)(fTime-fTimeKalmanBin)/fTimeKalmanBin; //time bin
1947 if (fAlignTOFTPC->GetEntriesFast()<htime){
1948 fAlignTOFTPC->Expand(htime*2+20);
1950 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1952 // make Alignment object if doesn't exist
1953 align=new AliRelAlignerKalman();
1954 align->SetRunNumber(fRun);
1955 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1956 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1957 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1958 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1959 // align->SetRejectOutliers(kFALSE);
1960 align->SetRejectOutliers(kTRUE);
1961 align->SetRejectOutliersSigma2Median(kTRUE);
1963 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1964 align->SetMagField(fMagF);
1965 fAlignTOFTPC->AddAt(align,htime);
1967 align->AddTrackParams(&pTOF,&pTPC);
1968 Float_t dca[2],cov[3];
1969 track->GetImpactParameters(dca,cov);
1970 if (TMath::Abs(dca[0])<kMaxDy){
1971 FillResHistoTPCTOF(&pTPC,&pTOF);
1973 //align->SetTimeStamp(fTime);
1974 Double_t averageTime = fTime;
1975 if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1976 averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1977 //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1979 align->SetTimeStamp((Int_t)averageTime);
1981 //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1983 align->SetRunNumber(fRun );
1985 Int_t nupdates=align->GetNUpdates();
1986 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates+1));
1987 // align->SetRejectOutliers(kFALSE);
1988 align->SetRejectOutliers(kTRUE);
1989 align->SetRejectOutliersSigma2Median(kTRUE);
1991 TTreeSRedirector *cstream = GetDebugStreamer();
1992 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1993 TVectorD gpTPC(3), gdTPC(3);
1994 TVectorD gpTOF(3), gdTOF(3);
1995 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1996 pTPC.GetDirection(gdTPC.GetMatrixArray());
1997 pTOF.GetXYZ(gpTOF.GetMatrixArray());
1998 pTOF.GetDirection(gdTOF.GetMatrixArray());
1999 (*cstream)<<"toftpc"<<
2000 "run="<<fRun<< // run number
2001 "event="<<fEvent<< // event number
2002 "time="<<fTime<< // time stamp of event
2003 "trigger="<<fTrigger<< // trigger
2004 "mag="<<fMagF<< // magnetic field
2006 "nmed="<<kglast<< // number of entries to define median and RMS
2007 "vMed.="<<&vecMedian<< // median of deltas
2008 "vRMS.="<<&vecRMS<< // rms of deltas
2009 "vDelta.="<<&vecDelta<< // delta in respect to median
2010 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
2011 "t.="<<track<< // ful track - find proper cuts
2012 "a.="<<align<< // current alignment
2013 "pTOF.="<<&pTOF<< // track param TOF
2014 "pTPC.="<<&pTPC<< // track param TPC
2015 "gpTPC.="<<&gpTPC<< // global position TPC
2016 "gdTPC.="<<&gdTPC<< // global direction TPC
2017 "gpTOF.="<<&gpTOF<< // global position TOF
2018 "gdTOF.="<<&gdTOF<< // global position TOF
2024 void AliTPCcalibTime::BookDistortionMaps(){
2026 // Book ndimensional histograms of distortions/residuals
2027 // Only primary tracks are selected for analysis
2030 Double_t xminTrack[5], xmaxTrack[5];
2032 TString axisName[5];
2033 TString axisTitle[5];
2036 axisName[0] ="#Delta";
2037 axisTitle[0] ="#Delta";
2040 xminTrack[1] =-1.1; xmaxTrack[1]=1.1;
2041 axisName[1] ="tanTheta";
2042 axisTitle[1] ="tan(#Theta)";
2045 xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi();
2047 axisTitle[2] ="#phi";
2050 xminTrack[3] =-1.; xmaxTrack[3]=1.; // 0.33 GeV cut
2052 axisTitle[3] ="snp";
2055 xminTrack[4] =120.; xmaxTrack[4]=215.; // crossing radius for CE only
2057 axisTitle[4] ="r(cm)";
2060 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
2061 fResHistoTPCCE[0] = new THnSparseS("TPCCE#Delta_{Y} (cm)","#Delta_{Y} (cm)", 5, binsTrack,xminTrack, xmaxTrack);
2062 fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2063 fResHistoTPCvertex[0] = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2064 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
2065 fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2066 xminTrack[0] =-5; xmaxTrack[0]=5; //
2067 fResHistoTPCTOF[0] = new THnSparseS("TPCTOF#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2070 xminTrack[0] =-6.; xmaxTrack[0]=6.; //
2071 fResHistoTPCCE[1] = new THnSparseS("TPCCE#Delta_{Z} (cm)","#Delta_{Z} (cm)", 5, binsTrack,xminTrack, xmaxTrack);
2072 fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2073 fResHistoTPCvertex[1] = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2074 fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2075 xminTrack[0] =-5.; xmaxTrack[0]=5.; //
2076 fResHistoTPCTOF[1] = new THnSparseS("TPCTOF#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2079 xminTrack[0] =-0.015; xmaxTrack[0]=0.015; //
2080 fResHistoTPCCE[2] = new THnSparseS("TPCCE#Delta_{#phi}","#Delta_{#phi}", 5, binsTrack,xminTrack, xmaxTrack);
2081 fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2082 fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2083 fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2084 fResHistoTPCTOF[2] = new THnSparseS("TPCTOF#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2087 xminTrack[0] =-0.05; xmaxTrack[0]=0.05; //
2088 fResHistoTPCCE[3] = new THnSparseS("TPCCE#Delta_{#theta}","#Delta_{#theta}", 5, binsTrack,xminTrack, xmaxTrack);
2089 fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2090 fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2091 fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2092 fResHistoTPCTOF[3] = new THnSparseS("TPCTOF#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2095 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
2096 fResHistoTPCCE[4] = new THnSparseS("TPCCE#Delta_{1/pt}","#Delta_{1/pt}", 5, binsTrack,xminTrack, xmaxTrack);
2097 fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
2098 fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
2099 fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
2100 fResHistoTPCTOF[4] = new THnSparseS("TPCTOF#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
2102 for (Int_t ivar=0;ivar<4;ivar++){
2103 for (Int_t ivar2=0;ivar2<5;ivar2++){
2104 fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2105 fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2107 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2108 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2109 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2110 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2111 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2112 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2117 // Book vertex: time histograms
2119 Int_t binsVertex[2]={500, fTimeBins};
2120 Double_t aminVertex[2]={-5,fTimeStart};
2121 Double_t amaxVertex[2]={5, fTimeEnd};
2122 const char* hnames[12]={"TPCXAside", "TPCXCside","TPCXACdiff","TPCXAPCdiff",
2123 "TPCYAside", "TPCYCside","TPCYACdiff","TPCYAPCdiff",
2124 "TPCZAPCside", "TPCZAMCside","TPCZACdiff","TPCZAPCdiff"};
2125 const char* anames[12]={"x (cm) - A side ", "x (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{x} (cm) - TPC-Common",
2126 "y (cm) - A side ", "y (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{y} (cm) - TPC-Common",
2127 "z (cm)", "#Delta_{Z} (cm) A-C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{Z} (cm) TPC-common"};
2128 for (Int_t ihis=0; ihis<12; ihis++) {
2129 if (ihis>=8) aminVertex[0]=-20.;
2130 if (ihis>=8) amaxVertex[0]=20.;
2131 fTPCVertex[ihis]=new THnSparseF(hnames[ihis],hnames[ihis],2,binsVertex,aminVertex,amaxVertex);
2132 fTPCVertex[ihis]->GetAxis(1)->SetTitle("Time");
2133 fTPCVertex[ihis]->GetAxis(0)->SetTitle(anames[ihis]);
2136 Int_t binsVertexC[2]={40, 300};
2137 Double_t aminVertexC[2]={-20,-30};
2138 Double_t amaxVertexC[2]={20,30};
2139 const char* hnamesC[5]={"TPCA_TPC","TPCC_TPC","TPCA_ITS","TPCC_ITS","TPC_ITS"};
2140 for (Int_t ihis=0; ihis<5; ihis++) {
2141 fTPCVertexCorrelation[ihis]=new THnSparseF(hnamesC[ihis],hnamesC[ihis],2,binsVertexC,aminVertexC,amaxVertexC);
2142 fTPCVertexCorrelation[ihis]->GetAxis(1)->SetTitle("z (cm)");
2143 fTPCVertexCorrelation[ihis]->GetAxis(0)->SetTitle("z (cm)");
2148 void AliTPCcalibTime::FillResHistoTPCCE(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pTPCOut ){
2150 // fill residual histograms pTPCOut-pTPCin - trac crossing CE
2153 if (fMemoryMode<2) return;
2156 pTPCIn->GetXYZ(xyz);
2157 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2158 histoX[1]= pTPCIn->GetTgl();
2160 histoX[3]= pTPCIn->GetSnp();
2161 histoX[4]= pTPCIn->GetX();
2162 AliExternalTrackParam lout(*pTPCOut);
2163 lout.Rotate(pTPCIn->GetAlpha());
2164 lout.PropagateTo(pTPCIn->GetX(),fMagF);
2166 for (Int_t ihisto=0; ihisto<5; ihisto++){
2167 histoX[0]=lout.GetParameter()[ihisto]-pTPCIn->GetParameter()[ihisto];
2168 fResHistoTPCCE[ihisto]->Fill(histoX);
2171 void AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
2173 // fill residual histograms pTPCIn-pITSOut
2174 // Histogram is filled only for primary tracks
2178 pTPCIn->GetXYZ(xyz);
2179 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2180 histoX[1]= pTPCIn->GetTgl();
2182 histoX[3]= pTPCIn->GetSnp();
2183 AliExternalTrackParam lits(*pITSOut);
2184 lits.Rotate(pTPCIn->GetAlpha());
2185 lits.PropagateTo(pTPCIn->GetX(),fMagF);
2187 for (Int_t ihisto=0; ihisto<5; ihisto++){
2188 histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2189 fResHistoTPCITS[ihisto]->Fill(histoX);
2194 void AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
2196 // fill residual histograms pTPC - vertex
2197 // Histogram is filled only for primary tracks
2199 if (fMemoryMode<2) return;
2201 const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
2202 AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
2204 if (!(pTrack->GetConstrainedParam())) return;
2205 AliExternalTrackParam lits(*(pTrack->GetConstrainedParam()));
2206 if (TMath::Abs(pTrack->GetY())>3) return; // beam pipe
2207 pTPCvertex.Rotate(lits.GetAlpha());
2208 //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
2209 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
2210 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
2212 pTPCIn->GetXYZ(xyz);
2213 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2214 histoX[1]= pTPCIn->GetTgl();
2216 histoX[3]= pTPCIn->GetSnp();
2218 Float_t dca[2], cov[3];
2219 pTrack->GetImpactParametersTPC(dca,cov);
2220 for (Int_t ihisto=0; ihisto<5; ihisto++){
2221 histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2222 // if (ihisto<2) histoX[0]=dca[ihisto];
2223 fResHistoTPCvertex[ihisto]->Fill(histoX);
2228 void AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
2230 // fill resuidual histogram TPCout-TRDin
2232 if (fMemoryMode<2) return;
2235 pTPCOut->GetXYZ(xyz);
2236 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2237 histoX[1]= pTPCOut->GetTgl();
2239 histoX[3]= pTPCOut->GetSnp();
2241 AliExternalTrackParam ltrd(*pTRDIn);
2242 ltrd.Rotate(pTPCOut->GetAlpha());
2243 // ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
2244 AliTracker::PropagateTrackToBxByBz(<rd,pTPCOut->GetX(),0.1,0.1,kFALSE);
2246 for (Int_t ihisto=0; ihisto<5; ihisto++){
2247 histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
2248 fResHistoTPCTRD[ihisto]->Fill(histoX);
2253 void AliTPCcalibTime::FillResHistoTPCTOF(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTOFIn ){
2255 // fill resuidual histogram TPCout-TOFin
2256 // track propagated to the TOF position
2257 if (fMemoryMode<2) return;
2261 AliExternalTrackParam ltpc(*pTPCOut);
2262 ltpc.Rotate(pTOFIn->GetAlpha());
2263 AliTracker::PropagateTrackToBxByBz(<pc,pTOFIn->GetX(),0.1,0.1,kFALSE);
2266 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2267 histoX[1]= ltpc.GetTgl();
2269 histoX[3]= ltpc.GetSnp();
2271 for (Int_t ihisto=0; ihisto<2; ihisto++){
2272 histoX[0]=ltpc.GetParameter()[ihisto]-pTOFIn->GetParameter()[ihisto];
2273 fResHistoTPCTOF[ihisto]->Fill(histoX);