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:
18 1. What do we calibrate.
19 2. How to interpret results
21 4. Analysis using debug streamers.
26 // To make cosmic scan the user interaction neccessary
29 gSystem->Load("libANALYSIS");
30 gSystem->Load("libTPCcalib");
31 TFile fcalib("CalibObjects.root");
32 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
33 AliTPCcalibCosmic * cosmic = ( AliTPCcalibCosmic *)array->FindObject("cosmicTPC");
41 #include "Riostream.h"
51 #include "THnSparse.h"
53 #include "AliTPCclusterMI.h"
54 #include "AliTPCseed.h"
55 #include "AliESDVertex.h"
56 #include "AliESDEvent.h"
57 #include "AliESDfriend.h"
58 #include "AliESDInputHandler.h"
59 #include "AliAnalysisManager.h"
61 #include "AliTracker.h"
63 #include "AliTPCCalROC.h"
67 #include "AliTPCcalibCosmic.h"
68 #include "TTreeStream.h"
69 #include "AliTPCTracklet.h"
70 #include "AliESDcosmic.h"
73 ClassImp(AliTPCcalibCosmic)
76 AliTPCcalibCosmic::AliTPCcalibCosmic()
85 fCutMaxD(5), // maximal distance in rfi ditection
86 fCutMaxDz(40), // maximal distance in z ditection
87 fCutTheta(0.03), // maximal distan theta
88 fCutMinDir(-0.99) // direction vector products
90 AliInfo("Default Constructor");
91 for (Int_t ihis=0; ihis<6;ihis++){
95 for (Int_t ihis=0; ihis<4;ihis++){
96 fHistodEdxMax[ihis] =0;
97 fHistodEdxTot[ihis] =0;
102 AliTPCcalibCosmic::AliTPCcalibCosmic(const Text_t *name, const Text_t *title)
111 fCutMaxD(5), // maximal distance in rfi ditection
112 fCutMaxDz(40), // maximal distance in z ditection
113 fCutTheta(0.03), // maximal distan theta
114 fCutMinDir(-0.99) // direction vector products
119 fHistNTracks = new TH1F("ntracks","Number of Tracks per Event; number of tracks per event; number of tracks",501,-0.5,500.5);
120 fClusters = new TH1F("signal","Number of Clusters per track; number of clusters per track n_{cl}; counts",160,0,160);
121 fModules = new TH2F("sector","Acorde hits; z (cm); x(cm)",1200,-650,650,600,-700,700);
122 fHistPt = new TH1F("Pt","Pt distribution; p_{T} (GeV); counts",2000,0,50);
123 fDeDx = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",500,0.01,100.,500,2.,1000);
125 fDeDxMIP = new TH1F("DeDxMIP","MIP region; TPC signal (a.u.);counts ",500,2.,1000);
127 AliInfo("Non Default Constructor");
131 AliTPCcalibCosmic::~AliTPCcalibCosmic(){
135 for (Int_t ihis=0; ihis<6;ihis++){
136 delete fHistoDelta[ihis];
137 delete fHistoPull[ihis];
139 for (Int_t ihis=0; ihis<4;ihis++){
140 delete fHistodEdxTot[ihis];
141 delete fHistodEdxMax[ihis];
144 delete fHistNTracks; // histogram showing number of ESD tracks per event
145 delete fClusters; // histogram showing the number of clusters per track
146 delete fModules; // 2d histogram of tracks which are propagated to the ACORDE scintillator array
147 delete fHistPt; // Pt histogram of reconstructed tracks
148 delete fDeDx; // dEdx spectrum showing the different particle types
149 delete fDeDxMIP; // TPC signal close to the MIP region of muons 0.4 < p < 0.45 GeV
153 void AliTPCcalibCosmic::Init(){
156 // Make performance histograms
159 // tracking performance bins
160 // 0 - delta of interest
161 // 1 - min (track0, track1) number of clusters
162 // 2 - R - vertex radius
164 // 4 - P2 - snp(phi) at inner wall of TPC
165 // 5 - P3 - tan(theta) at inner wall of TPC
166 // 6 - P4 - 1/pt mean
170 Double_t xminTrack[9], xmaxTrack[9];
175 axisName[0] ="#Delta";
178 xminTrack[1] =80; xmaxTrack[1]=160;
179 axisName[1] ="N_{cl}";
182 xminTrack[2] =0; xmaxTrack[2]=90; //
183 axisName[2] ="dca_{r} (cm)";
186 xminTrack[3] =-250; xmaxTrack[3]=250; //
187 axisName[3] ="z (cm)";
190 xminTrack[4] =-0.8; xmaxTrack[4]=0.8; //
191 axisName[4] ="sin(#phi)";
194 xminTrack[5] =-1; xmaxTrack[5]=1; //
195 axisName[5] ="tan(#theta)";
198 xminTrack[6] =0; xmaxTrack[6]=2; //
199 axisName[6] ="1/pt (1/GeV)";
202 xminTrack[7] =0.2; xmaxTrack[7]=50; //
203 axisName[7] ="pt (GeV)";
206 xminTrack[8] =0; xmaxTrack[8]=TMath::Pi(); //
207 axisName[8] ="alpha";
210 xminTrack[0] =-1; xmaxTrack[0]=1; //
211 fHistoDelta[0] = new THnSparseS("#Delta_{Y} (cm)","#Delta_{Y} (cm)", 9, binsTrack,xminTrack, xmaxTrack);
212 xminTrack[0] =-5; xmaxTrack[0]=5; //
213 fHistoPull[0] = new THnSparseS("#Delta_{Y} (unit)","#Delta_{Y} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
216 xminTrack[0] =-1; xmaxTrack[0]=1; //
217 fHistoDelta[1] = new THnSparseS("#Delta_{Z} (cm)","#Delta_{Z} (cm)", 9, binsTrack,xminTrack, xmaxTrack);
218 xminTrack[0] =-5; xmaxTrack[0]=5; //
219 fHistoPull[1] = new THnSparseS("#Delta_{Z} (unit)","#Delta_{Z} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
222 xminTrack[0] =-10; xmaxTrack[0]=10; //
223 fHistoDelta[2] = new THnSparseS("#Delta_{#phi} (mrad)","#Delta_{#phi} (mrad)", 9, binsTrack,xminTrack, xmaxTrack);
224 xminTrack[0] =-5; xmaxTrack[0]=5; //
225 fHistoPull[2] = new THnSparseS("#Delta_{#phi} (unit)","#Delta_{#phi} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
228 xminTrack[0] =-10; xmaxTrack[0]=10; //
229 fHistoDelta[3] = new THnSparseS("#Delta_{#theta} (mrad)","#Delta_{#theta} (mrad)", 9, binsTrack,xminTrack, xmaxTrack);
230 xminTrack[0] =-5; xmaxTrack[0]=5; //
231 fHistoPull[3] = new THnSparseS("#Delta_{#theta} (unit)","#Delta_{#theta} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
234 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
235 fHistoDelta[4] = new THnSparseS("#Delta_{1/pt} (1/GeV)","#Delta_{1/pt} (1/GeV)", 9, binsTrack,xminTrack, xmaxTrack);
236 xminTrack[0] =-5; xmaxTrack[0]=5; //
237 fHistoPull[4] = new THnSparseS("#Delta_{1/pt} (unit)","#Delta_{1/pt} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
241 xminTrack[0] =-0.5; xmaxTrack[0]=0.5; //
242 fHistoDelta[5] = new THnSparseS("#Delta_{pt}/p_{t}","#Delta_{pt}/p_{t}", 9, binsTrack,xminTrack, xmaxTrack);
243 xminTrack[0] =-5; xmaxTrack[0]=5; //
244 fHistoPull[5] = new THnSparseS("#Delta_{pt}/p_{t} (unit)","#Delta_{pt}/p_{t} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
247 for (Int_t idedx=0;idedx<4;idedx++){
248 xminTrack[0] =0.5; xmaxTrack[0]=1.5; //
250 xminTrack[1] =10; xmaxTrack[1]=160;
252 fHistodEdxMax[idedx] = new THnSparseS(Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
253 Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
254 9, binsTrack,xminTrack, xmaxTrack);
255 fHistodEdxTot[idedx] = new THnSparseS(Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
256 Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
257 9, binsTrack,xminTrack, xmaxTrack);
262 for (Int_t ivar=0;ivar<6;ivar++){
263 for (Int_t ivar2=0;ivar2<9;ivar2++){
264 fHistoDelta[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
265 fHistoDelta[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
266 fHistoPull[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
267 fHistoPull[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
268 BinLogX(fHistoDelta[ivar],7);
269 BinLogX(fHistoPull[ivar],7);
271 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
272 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
273 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
274 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
275 BinLogX(fHistodEdxMax[ivar],7);
276 BinLogX(fHistodEdxTot[ivar],7);
283 void AliTPCcalibCosmic::Add(const AliTPCcalibCosmic* cosmic){
287 for (Int_t ivar=0; ivar<6;ivar++){
288 if (fHistoDelta[ivar] && cosmic->fHistoDelta[ivar]){
289 fHistoDelta[ivar]->Add(cosmic->fHistoDelta[ivar]);
291 if (fHistoPull[ivar] && cosmic->fHistoPull[ivar]){
292 fHistoPull[ivar]->Add(cosmic->fHistoPull[ivar]);
295 for (Int_t ivar=0; ivar<4;ivar++){
296 if (fHistodEdxMax[ivar] && cosmic->fHistodEdxMax[ivar]){
297 fHistodEdxMax[ivar]->Add(cosmic->fHistodEdxMax[ivar]);
299 if (fHistodEdxTot[ivar] && cosmic->fHistodEdxTot[ivar]){
300 fHistodEdxTot[ivar]->Add(cosmic->fHistodEdxTot[ivar]);
308 void AliTPCcalibCosmic::Process(AliESDEvent *event) {
313 Printf("ERROR: ESD not available");
316 AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
318 Printf("ERROR: ESDfriend not available");
323 FindPairs(event); // nearly everything takes place in find pairs...
325 if (GetDebugLevel()>20) printf("Hallo world: Im here and processing an event\n");
326 Int_t ntracks=event->GetNumberOfTracks();
327 fHistNTracks->Fill(ntracks);
328 if (ntracks==0) return;
329 AliESDcosmic cosmicESD;
330 TTreeSRedirector * cstream = GetDebugStreamer();
331 cosmicESD.SetDebugStreamer(cstream);
332 cosmicESD.ProcessEvent(event);
333 if (cstream) cosmicESD.DumpToTree();
339 void AliTPCcalibCosmic::FillHistoPerformance(AliExternalTrackParam *par0, AliExternalTrackParam *par1, AliExternalTrackParam *inner0, AliExternalTrackParam *inner1, AliTPCseed *seed0, AliTPCseed *seed1){
343 Int_t kMinCldEdx =20;
344 Int_t ncl0 = seed0->GetNumberOfClusters();
345 Int_t ncl1 = seed1->GetNumberOfClusters();
347 const Double_t kpullCut = 10;
353 Double_t radius0 = TMath::Sqrt(xyz0[0]*xyz0[0]+xyz0[1]*xyz0[1]);
354 Double_t radius1 = TMath::Sqrt(xyz1[0]*xyz1[0]+xyz1[1]*xyz1[1]);
355 inner0->GetXYZ(xyz0);
356 Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
358 x[1] = TMath::Min(ncl0,ncl1);
359 x[2] = (radius0+radius1)*0.5;
360 x[3] = (inner0->GetZ()+inner1->GetZ())*0.5;
361 x[4] = (inner0->GetSnp()-inner1->GetSnp())*0.5;
362 x[5] = (inner0->GetTgl()-inner1->GetTgl())*0.5;
363 x[6] = (1/par0->Pt()+1/par1->Pt())*0.5;
364 x[7] = (par0->Pt()+par1->Pt())*0.5;
369 delta[0] = (par0->GetY()+par1->GetY());
370 delta[1] = (par0->GetZ()-par1->GetZ());
371 delta[2] = (par0->GetAlpha()-par1->GetAlpha()-TMath::Pi());
372 delta[3] = (par0->GetTgl()+par1->GetTgl());
373 delta[4] = (par0->GetParameter()[4]+par1->GetParameter()[4]);
374 delta[5] = (par0->Pt()-par1->Pt())/((par0->Pt()+par1->Pt())*0.5);
376 sigma[0] = TMath::Sqrt(par0->GetSigmaY2()+par1->GetSigmaY2());
377 sigma[1] = TMath::Sqrt(par0->GetSigmaZ2()+par1->GetSigmaZ2());
378 sigma[2] = TMath::Sqrt(par0->GetSigmaSnp2()+par1->GetSigmaSnp2());
379 sigma[3] = TMath::Sqrt(par0->GetSigmaTgl2()+par1->GetSigmaTgl2());
380 sigma[4] = TMath::Sqrt(par0->GetSigma1Pt2()+par1->GetSigma1Pt2());
381 sigma[5] = sigma[4]*((par0->Pt()+par1->Pt())*0.5);
384 for (Int_t ivar=0;ivar<6;ivar++){
385 if (sigma[ivar]==0) isOK=kFALSE;
386 x[0]= delta[ivar]/sigma[ivar];
387 if (TMath::Abs(x[0])>kpullCut) isOK = kFALSE;
391 if (isOK) for (Int_t ivar=0;ivar<6;ivar++){
392 x[0]= delta[ivar]/TMath::Sqrt(2);
393 if (ivar==2 || ivar ==3) x[0]*=1000;
394 fHistoDelta[ivar]->Fill(x);
396 x[0]= delta[ivar]/sigma[ivar];
397 fHistoPull[ivar]->Fill(x);
402 // Fill dedx performance
404 for (Int_t ipad=0; ipad<4;ipad++){
410 if (ipad==0) row1=63;
411 if (ipad==1) {row0=63; row1=63+64;}
412 if (ipad==2) {row0=128;}
413 Int_t nclUp = TMath::Nint(seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
414 Int_t nclDown = TMath::Nint(seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
415 Int_t minCl = TMath::Min(nclUp,nclDown);
416 if (minCl<kMinCldEdx) continue;
419 Float_t dEdxTotUp = seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
420 Float_t dEdxTotDown = seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
421 Float_t dEdxMaxUp = seed0->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
422 Float_t dEdxMaxDown = seed1->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
424 if (dEdxTotDown<=0) continue;
425 if (dEdxMaxDown<=0) continue;
426 x[0]=dEdxTotUp/dEdxTotDown;
427 fHistodEdxTot[ipad]->Fill(x);
428 x[0]=dEdxMaxUp/dEdxMaxDown;
429 fHistodEdxMax[ipad]->Fill(x);
438 void AliTPCcalibCosmic::Analyze() {
440 fMIPvalue = CalculateMIPvalue(fDeDxMIP);
448 void AliTPCcalibCosmic::FindPairs(AliESDEvent *event) {
452 // Track0 is choosen in upper TPC part
453 // Track1 is choosen in lower TPC part
455 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
456 AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
457 Int_t ntracks=event->GetNumberOfTracks();
458 TObjArray tpcSeeds(ntracks);
459 if (ntracks==0) return;
460 Double_t vtxx[3]={0,0,0};
461 Double_t svtxx[3]={0.000001,0.000001,100.};
462 AliESDVertex vtx(vtxx,svtxx);
466 for (Int_t i=0;i<ntracks;++i) {
467 AliESDtrack *track = event->GetTrack(i);
468 fClusters->Fill(track->GetTPCNcls());
470 const AliExternalTrackParam * trackIn = track->GetInnerParam();
471 const AliExternalTrackParam * trackOut = track->GetOuterParam();
472 if (!trackIn) continue;
473 if (!trackOut) continue;
474 if (ntracks>4 && TMath::Abs(trackIn->GetTgl())<0.0015) continue; // filter laser
477 AliESDfriendTrack *friendTrack = ESDfriend->GetTrack(i);
478 TObject *calibObject;
479 AliTPCseed *seed = 0;
480 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
481 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
483 if (seed) tpcSeeds.AddAt(seed,i);
485 Double_t meanP = 0.5*(trackIn->GetP() + trackOut->GetP());
486 if (seed && track->GetTPCNcls() > 80 + 60/(1+TMath::Exp(-meanP+5))) {
487 fDeDx->Fill(meanP, seed->CookdEdxNorm(0.0,0.45,0,0,159));
489 if (meanP > 0.4 && meanP < 0.45) fDeDxMIP->Fill(seed->CookdEdxNorm(0.0,0.45,0,0,159));
491 if (GetDebugLevel()>0&&meanP>0.2&&seed->CookdEdxNorm(0.0,0.45,0,0,159)>300) {
492 TFile *curfile = AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile();
493 if (curfile) printf(">>> p+ in file: %s \t event: %i \t Number of ESD tracks: %i \n", curfile->GetName(), (int)event->GetEventNumberInFile(), (int)ntracks);
494 if (track->GetOuterParam()->GetAlpha()<0) cout << " Polartiy: " << track->GetSign() << endl;
501 if (ntracks<2) return;
505 for (Int_t i=0;i<ntracks;++i) {
506 AliESDtrack *track0 = event->GetTrack(i);
507 // track0 - choosen upper part
508 if (!track0) continue;
509 if (!track0->GetOuterParam()) continue;
510 if (track0->GetOuterParam()->GetAlpha()<0) continue;
512 track0->GetDirection(dir0);
513 for (Int_t j=0;j<ntracks;++j) {
515 AliESDtrack *track1 = event->GetTrack(j);
517 if (!track1) continue;
518 if (!track1->GetOuterParam()) continue;
519 if (track1->GetOuterParam()->GetAlpha()>0) continue;
522 track1->GetDirection(dir1);
524 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
525 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
526 if (! seed0) continue;
527 if (! seed1) continue;
528 Float_t dedx0 = seed0->CookdEdxNorm(0.05,0.55,0,0,159);
529 Float_t dedx1 = seed1->CookdEdxNorm(0.05,0.55,0,0,159);
531 Float_t dedx0I = seed0->CookdEdxNorm(0.05,0.55,0,0,63);
532 Float_t dedx1I = seed1->CookdEdxNorm(0.05,0.55,0,0,63);
534 Float_t dedx0O = seed0->CookdEdxNorm(0.05,0.55,0,64,159);
535 Float_t dedx1O = seed1->CookdEdxNorm(0.05,0.55,0,64,159);
537 Float_t dir = (dir0[0]*dir1[0] + dir0[1]*dir1[1] + dir0[2]*dir1[2]);
538 Float_t d0 = track0->GetLinearD(0,0);
539 Float_t d1 = track1->GetLinearD(0,0);
541 // conservative cuts - convergence to be guarantied
542 // applying before track propagation
543 if (TMath::Abs(d0+d1)>fCutMaxD) continue; // distance to the 0,0
544 if (dir>fCutMinDir) continue; // direction vector product
545 Float_t bz = AliTracker::GetBz();
546 Float_t dvertex0[2]; //distance to 0,0
547 Float_t dvertex1[2]; //distance to 0,0
548 track0->GetDZ(0,0,0,bz,dvertex0);
549 track1->GetDZ(0,0,0,bz,dvertex1);
550 if (TMath::Abs(dvertex0[1])>250) continue;
551 if (TMath::Abs(dvertex1[1])>250) continue;
555 Float_t dmax = TMath::Max(TMath::Abs(d0),TMath::Abs(d1));
556 AliExternalTrackParam param0(*track0);
557 AliExternalTrackParam param1(*track1);
559 // Propagate using Magnetic field and correct fo material budget
561 AliTracker::PropagateTrackTo(¶m0,dmax+1,0.0005,3,kTRUE);
562 AliTracker::PropagateTrackTo(¶m1,dmax+1,0.0005,3,kTRUE);
564 // Propagate rest to the 0,0 DCA - z should be ignored
566 Bool_t b0 = param0.PropagateToDCA(&vtx,bz,1000);
567 Bool_t b1 = param1.PropagateToDCA(&vtx,bz,1000);
569 param0.GetDZ(0,0,0,bz,dvertex0);
570 param1.GetDZ(0,0,0,bz,dvertex1);
571 if (TMath::Abs(param0.GetZ()-param1.GetZ())>fCutMaxDz) continue;
573 Double_t xyz0[3];//,pxyz0[3];
574 Double_t xyz1[3];//,pxyz1[3];
577 Bool_t isPair = IsPair(¶m0,¶m1);
579 if (isPair) FillAcordeHist(track0);
581 // combined track params
583 AliExternalTrackParam *par0U=MakeCombinedTrack(¶m0,¶m1);
584 AliExternalTrackParam *par1U=MakeCombinedTrack(¶m1,¶m0);
589 TTreeSRedirector * cstream = GetDebugStreamer();
590 //printf("My stream=%p\n",(void*)cstream);
591 AliExternalTrackParam *ip0 = (AliExternalTrackParam *)track0->GetInnerParam();
592 AliExternalTrackParam *ip1 = (AliExternalTrackParam *)track1->GetInnerParam();
593 AliExternalTrackParam *op0 = (AliExternalTrackParam *)track0->GetOuterParam();
594 AliExternalTrackParam *op1 = (AliExternalTrackParam *)track1->GetOuterParam();
595 Bool_t isCrossI = ip0->GetZ()*ip1->GetZ()<0;
596 Bool_t isCrossO = op0->GetZ()*op1->GetZ()<0;
597 Double_t alpha0 = TMath::ATan2(dir0[1],dir0[0]);
598 Double_t alpha1 = TMath::ATan2(dir1[1],dir1[0]);
602 FillHistoPerformance(¶m0, ¶m1, ip0, ip1, seed0, seed1);
605 (*cstream) << "Track0" <<
606 "run="<<fRun<< // run number
607 "event="<<fEvent<< // event number
608 "time="<<fTime<< // time stamp of event
609 "trigger="<<fTrigger<< // trigger
610 "triggerClass="<<&fTriggerClass<< // trigger
611 "mag="<<fMagF<< // magnetic field
612 "dir="<<dir<< // direction
613 "OK="<<isPair<< // will be accepted
614 "b0="<<b0<< // propagate status
615 "b1="<<b1<< // propagate status
616 "crossI="<<isCrossI<< // cross inner
617 "crossO="<<isCrossO<< // cross outer
619 "Orig0.=" << track0 << // original track 0
620 "Orig1.=" << track1 << // original track 1
621 "Tr0.="<<¶m0<< // track propagated to the DCA 0,0
622 "Tr1.="<<¶m1<< // track propagated to the DCA 0,0
623 "Ip0.="<<ip0<< // inner param - upper
624 "Ip1.="<<ip1<< // inner param - lower
625 "Op0.="<<op0<< // outer param - upper
626 "Op1.="<<op1<< // outer param - lower
627 "Up0.="<<par0U<< // combined track 0
628 "Up1.="<<par1U<< // combined track 1
630 "v00="<<dvertex0[0]<< // distance using kalman
631 "v01="<<dvertex0[1]<< //
632 "v10="<<dvertex1[0]<< //
633 "v11="<<dvertex1[1]<< //
634 "d0="<<d0<< // linear distance to 0,0
635 "d1="<<d1<< // linear distance to 0,0
639 "x00="<<xyz0[0]<< // global position close to vertex
643 "x10="<<xyz1[0]<< // global position close to vertex
649 "dir00="<<dir0[0]<< // direction upper
653 "dir10="<<dir1[0]<< // direction lower
658 "Seed0.=" << seed0 << // original seed 0
659 "Seed1.=" << seed1 << // original seed 1
661 "dedx0="<<dedx0<< // dedx0 - all
662 "dedx1="<<dedx1<< // dedx1 - all
664 "dedx0I="<<dedx0I<< // dedx0 - inner ROC
665 "dedx1I="<<dedx1I<< // dedx1 - inner ROC
667 "dedx0O="<<dedx0O<< // dedx0 - outer ROC
668 "dedx1O="<<dedx1O<< // dedx1 - outer ROC
681 void AliTPCcalibCosmic::FillAcordeHist(AliESDtrack *upperTrack) {
683 // Pt cut to select straight tracks which can be easily propagated to ACORDE which is outside the magnetic field
684 if (upperTrack->Pt() < 10 || upperTrack->GetTPCNcls() < 80) return;
686 const Double_t AcordePlane = 850.; // distance of the central Acorde detectors to the beam line at y =0
687 const Double_t roof = 210.5; // distance from x =0 to end of magnet roof
690 upperTrack->GetXYZ(r);
692 upperTrack->GetDirection(d);
694 z = r[2] + (d[2]/d[1])*(AcordePlane - r[1]);
695 x = r[0] + (d[0]/d[1])*(AcordePlane - r[1]);
698 x = r[0] + (d[0]/(d[0]+d[1]))*(AcordePlane+roof-r[0]-r[1]);
699 z = r[2] + (d[2]/(d[0]+d[1]))*(AcordePlane+roof-r[0]-r[1]);
702 x = r[0] + (d[0]/(d[1]-d[0]))*(AcordePlane+roof+r[0]-r[1]);
703 z = r[2] + (d[2]/(d[1]-d[0]))*(AcordePlane+roof+r[0]-r[1]);
706 fModules->Fill(z, x);
712 Long64_t AliTPCcalibCosmic::Merge(TCollection *li) {
714 TIterator* iter = li->MakeIterator();
715 AliTPCcalibCosmic* cal = 0;
717 while ((cal = (AliTPCcalibCosmic*)iter->Next())) {
718 if (!cal->InheritsFrom(AliTPCcalibCosmic::Class())) {
719 //Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
723 fHistNTracks->Add(cal->GetHistNTracks());
724 fClusters->Add(cal-> GetHistClusters());
725 fModules->Add(cal->GetHistAcorde());
726 fHistPt->Add(cal->GetHistPt());
727 fDeDx->Add(cal->GetHistDeDx());
728 fDeDxMIP->Add(cal->GetHistMIP());
736 Bool_t AliTPCcalibCosmic::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
740 // 0. Same direction - OPOSITE - cutDir +cutT
741 TCut cutDir("cutDir","dir<-0.99")
743 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
746 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
750 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
753 const Double_t *p0 = tr0->GetParameter();
754 const Double_t *p1 = tr1->GetParameter();
755 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
756 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
757 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
759 Double_t d0[3], d1[3];
760 tr0->GetDirection(d0);
761 tr1->GetDirection(d1);
762 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
769 Double_t AliTPCcalibCosmic::CalculateMIPvalue(TH1F * hist) {
771 TF1 * funcDoubleGaus = new TF1("funcDoubleGaus", "gaus(0)+gaus(3)",0,1000);
772 funcDoubleGaus->SetParameters(hist->GetEntries()*0.75,hist->GetMean()/1.3,hist->GetMean()*0.10,
773 hist->GetEntries()*0.25,hist->GetMean()*1.3,hist->GetMean()*0.10);
774 hist->Fit(funcDoubleGaus);
775 Double_t MIPvalue = TMath::Min(funcDoubleGaus->GetParameter(1),funcDoubleGaus->GetParameter(4));
777 delete funcDoubleGaus;
786 void AliTPCcalibCosmic::CalculateBetheParams(TH2F */*hist*/, Double_t * /*initialParam*/) {
788 // Not implemented yet
795 void AliTPCcalibCosmic::BinLogX(THnSparse *h, Int_t axisDim) {
797 // Method for the correct logarithmic binning of histograms
799 TAxis *axis = h->GetAxis(axisDim);
800 int bins = axis->GetNbins();
802 Double_t from = axis->GetXmin();
803 Double_t to = axis->GetXmax();
804 Double_t *new_bins = new Double_t[bins + 1];
807 Double_t factor = pow(to/from, 1./bins);
809 for (int i = 1; i <= bins; i++) {
810 new_bins[i] = factor * new_bins[i-1];
812 axis->Set(bins, new_bins);
818 void AliTPCcalibCosmic::BinLogX(TH1 *h) {
820 // Method for the correct logarithmic binning of histograms
822 TAxis *axis = h->GetXaxis();
823 int bins = axis->GetNbins();
825 Double_t from = axis->GetXmin();
826 Double_t to = axis->GetXmax();
827 Double_t *new_bins = new Double_t[bins + 1];
830 Double_t factor = pow(to/from, 1./bins);
832 for (int i = 1; i <= bins; i++) {
833 new_bins[i] = factor * new_bins[i-1];
835 axis->Set(bins, new_bins);
841 AliExternalTrackParam *AliTPCcalibCosmic::MakeTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
845 AliExternalTrackParam *par1R= new AliExternalTrackParam(*track1);
846 par1R->Rotate(track0->GetAlpha());
847 par1R->PropagateTo(track0->GetX(),AliTracker::GetBz());
850 Double_t * param = (Double_t*)par1R->GetParameter();
851 Double_t * covar = (Double_t*)par1R->GetCovariance();
859 covar[6] *=-1.; covar[7] *=-1.; covar[8] *=-1.;
860 //covar[10]*=-1.; covar[11]*=-1.; covar[12]*=-1.;
865 AliExternalTrackParam *AliTPCcalibCosmic::MakeCombinedTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
867 // Make combined track
870 AliExternalTrackParam * par1T = MakeTrack(track0,track1);
871 AliExternalTrackParam * par0U = new AliExternalTrackParam(*track0);
873 UpdateTrack(*par0U,*par1T);
879 void AliTPCcalibCosmic::UpdateTrack(AliExternalTrackParam &track1, const AliExternalTrackParam &track2){
881 // Update track 1 with track 2
885 TMatrixD vecXk(5,1); // X vector
886 TMatrixD covXk(5,5); // X covariance
887 TMatrixD matHk(5,5); // vector to mesurement
888 TMatrixD measR(5,5); // measurement error
889 TMatrixD vecZk(5,1); // measurement
891 TMatrixD vecYk(5,1); // Innovation or measurement residual
892 TMatrixD matHkT(5,5);
893 TMatrixD matSk(5,5); // Innovation (or residual) covariance
894 TMatrixD matKk(5,5); // Optimal Kalman gain
895 TMatrixD mat1(5,5); // update covariance matrix
896 TMatrixD covXk2(5,5); //
897 TMatrixD covOut(5,5);
899 Double_t *param1=(Double_t*) track1.GetParameter();
900 Double_t *covar1=(Double_t*) track1.GetCovariance();
901 Double_t *param2=(Double_t*) track2.GetParameter();
902 Double_t *covar2=(Double_t*) track2.GetCovariance();
904 // copy data to the matrix
905 for (Int_t ipar=0; ipar<5; ipar++){
906 for (Int_t jpar=0; jpar<5; jpar++){
907 covXk(ipar,jpar) = covar1[track1.GetIndex(ipar, jpar)];
908 measR(ipar,jpar) = covar2[track2.GetIndex(ipar, jpar)];
912 vecXk(ipar,0) = param1[ipar];
913 vecZk(ipar,0) = param2[ipar];
922 vecYk = vecZk-matHk*vecXk; // Innovation or measurement residual
923 matHkT=matHk.T(); matHk.T();
924 matSk = (matHk*(covXk*matHkT))+measR; // Innovation (or residual) covariance
926 matKk = (covXk*matHkT)*matSk; // Optimal Kalman gain
927 vecXk += matKk*vecYk; // updated vector
928 covXk2 = (mat1-(matKk*matHk));
929 covOut = covXk2*covXk;
933 // copy from matrix to parameters
946 for (Int_t ipar=0; ipar<5; ipar++){
947 param1[ipar]= vecXk(ipar,0) ;
948 for (Int_t jpar=0; jpar<5; jpar++){
949 covar1[track1.GetIndex(ipar, jpar)]=covOut(ipar,jpar);