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");
38 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
39 gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
40 AliXRDPROOFtoolkit tool;
41 TChain * chainCosmic = tool.MakeChainRandom("cosmicF.txt","Track0",0,10000);
43 TCut cutITSN="min(Orig0.fITSncls,Orig1.fITSncls)>2";
44 TCut cutTPCN="min(Orig0.fTPCncls,Orig1.fTPCncls)>120";
46 chainCosmic->Draw(">>listITS",cutITSN+cutTPCN,"entryList");
47 TEntryList *elistITS = (TEntryList*)gDirectory->Get("listITS");
48 chainCosmic->SetEntryList(elistITS);
54 #include "Riostream.h"
64 #include "THnSparse.h"
65 #include "TDatabasePDG.h"
67 #include "AliTPCclusterMI.h"
68 #include "AliTPCseed.h"
69 #include "AliESDVertex.h"
70 #include "AliESDEvent.h"
71 #include "AliESDfriend.h"
72 #include "AliESDInputHandler.h"
73 #include "AliAnalysisManager.h"
75 #include "AliTracker.h"
77 #include "AliTPCCalROC.h"
81 #include "AliTPCcalibCosmic.h"
82 #include "TTreeStream.h"
83 #include "AliTPCTracklet.h"
84 //#include "AliESDcosmic.h"
87 ClassImp(AliTPCcalibCosmic)
90 AliTPCcalibCosmic::AliTPCcalibCosmic()
99 fCutMaxD(5), // maximal distance in rfi ditection
100 fCutMaxDz(40), // maximal distance in z ditection
101 fCutTheta(0.03), // maximal distan theta
102 fCutMinDir(-0.99) // direction vector products
104 AliInfo("Default Constructor");
105 for (Int_t ihis=0; ihis<6;ihis++){
109 for (Int_t ihis=0; ihis<4;ihis++){
110 fHistodEdxMax[ihis] =0;
111 fHistodEdxTot[ihis] =0;
116 AliTPCcalibCosmic::AliTPCcalibCosmic(const Text_t *name, const Text_t *title)
125 fCutMaxD(5), // maximal distance in rfi ditection
126 fCutMaxDz(40), // maximal distance in z ditection
127 fCutTheta(0.03), // maximal distan theta
128 fCutMinDir(-0.99) // direction vector products
133 fHistNTracks = new TH1F("ntracks","Number of Tracks per Event; number of tracks per event; number of tracks",501,-0.5,500.5);
134 fClusters = new TH1F("signal","Number of Clusters per track; number of clusters per track n_{cl}; counts",160,0,160);
135 fModules = new TH2F("sector","Acorde hits; z (cm); x(cm)",1200,-650,650,600,-700,700);
136 fHistPt = new TH1F("Pt","Pt distribution; p_{T} (GeV); counts",2000,0,50);
137 fDeDx = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",500,0.01,100.,500,2.,1000);
139 fDeDxMIP = new TH1F("DeDxMIP","MIP region; TPC signal (a.u.);counts ",500,2.,1000);
141 AliInfo("Non Default Constructor");
145 AliTPCcalibCosmic::~AliTPCcalibCosmic(){
149 for (Int_t ihis=0; ihis<6;ihis++){
150 delete fHistoDelta[ihis];
151 delete fHistoPull[ihis];
153 for (Int_t ihis=0; ihis<4;ihis++){
154 delete fHistodEdxTot[ihis];
155 delete fHistodEdxMax[ihis];
158 delete fHistNTracks; // histogram showing number of ESD tracks per event
159 delete fClusters; // histogram showing the number of clusters per track
160 delete fModules; // 2d histogram of tracks which are propagated to the ACORDE scintillator array
161 delete fHistPt; // Pt histogram of reconstructed tracks
162 delete fDeDx; // dEdx spectrum showing the different particle types
163 delete fDeDxMIP; // TPC signal close to the MIP region of muons 0.4 < p < 0.45 GeV
167 void AliTPCcalibCosmic::Init(){
170 // Make performance histograms
173 // tracking performance bins
174 // 0 - delta of interest
175 // 1 - min (track0, track1) number of clusters
176 // 2 - R - vertex radius
178 // 4 - P2 - snp(phi) at inner wall of TPC
179 // 5 - P3 - tan(theta) at inner wall of TPC
180 // 6 - P4 - 1/pt mean
184 Double_t xminTrack[9], xmaxTrack[9];
189 axisName[0] ="#Delta";
192 xminTrack[1] =80; xmaxTrack[1]=160;
193 axisName[1] ="N_{cl}";
196 xminTrack[2] =0; xmaxTrack[2]=90; //
197 axisName[2] ="dca_{r} (cm)";
200 xminTrack[3] =-250; xmaxTrack[3]=250; //
201 axisName[3] ="z (cm)";
204 xminTrack[4] =-0.8; xmaxTrack[4]=0.8; //
205 axisName[4] ="sin(#phi)";
208 xminTrack[5] =-1; xmaxTrack[5]=1; //
209 axisName[5] ="tan(#theta)";
212 xminTrack[6] =-2; xmaxTrack[6]=2; //
213 axisName[6] ="1/pt (1/GeV)";
216 xminTrack[7] =0.2; xmaxTrack[7]=50; //
217 axisName[7] ="pt (GeV)";
220 xminTrack[8] =0; xmaxTrack[8]=TMath::Pi(); //
221 axisName[8] ="alpha";
224 xminTrack[0] =-1; xmaxTrack[0]=1; //
225 fHistoDelta[0] = new THnSparseS("#Delta_{Y} (cm)","#Delta_{Y} (cm)", 9, binsTrack,xminTrack, xmaxTrack);
226 xminTrack[0] =-5; xmaxTrack[0]=5; //
227 fHistoPull[0] = new THnSparseS("#Delta_{Y} (unit)","#Delta_{Y} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
230 xminTrack[0] =-1; xmaxTrack[0]=1; //
231 fHistoDelta[1] = new THnSparseS("#Delta_{Z} (cm)","#Delta_{Z} (cm)", 9, binsTrack,xminTrack, xmaxTrack);
232 xminTrack[0] =-5; xmaxTrack[0]=5; //
233 fHistoPull[1] = new THnSparseS("#Delta_{Z} (unit)","#Delta_{Z} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
236 xminTrack[0] =-10; xmaxTrack[0]=10; //
237 fHistoDelta[2] = new THnSparseS("#Delta_{#phi} (mrad)","#Delta_{#phi} (mrad)", 9, binsTrack,xminTrack, xmaxTrack);
238 xminTrack[0] =-5; xmaxTrack[0]=5; //
239 fHistoPull[2] = new THnSparseS("#Delta_{#phi} (unit)","#Delta_{#phi} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
242 xminTrack[0] =-10; xmaxTrack[0]=10; //
243 fHistoDelta[3] = new THnSparseS("#Delta_{#theta} (mrad)","#Delta_{#theta} (mrad)", 9, binsTrack,xminTrack, xmaxTrack);
244 xminTrack[0] =-5; xmaxTrack[0]=5; //
245 fHistoPull[3] = new THnSparseS("#Delta_{#theta} (unit)","#Delta_{#theta} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
248 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
249 fHistoDelta[4] = new THnSparseS("#Delta_{1/pt} (1/GeV)","#Delta_{1/pt} (1/GeV)", 9, binsTrack,xminTrack, xmaxTrack);
250 xminTrack[0] =-5; xmaxTrack[0]=5; //
251 fHistoPull[4] = new THnSparseS("#Delta_{1/pt} (unit)","#Delta_{1/pt} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
255 xminTrack[0] =-0.5; xmaxTrack[0]=0.5; //
256 fHistoDelta[5] = new THnSparseS("#Delta_{pt}/p_{t}","#Delta_{pt}/p_{t}", 9, binsTrack,xminTrack, xmaxTrack);
257 xminTrack[0] =-5; xmaxTrack[0]=5; //
258 fHistoPull[5] = new THnSparseS("#Delta_{pt}/p_{t} (unit)","#Delta_{pt}/p_{t} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
261 for (Int_t idedx=0;idedx<4;idedx++){
262 xminTrack[0] =0.5; xmaxTrack[0]=1.5; //
264 xminTrack[1] =10; xmaxTrack[1]=160;
266 fHistodEdxMax[idedx] = new THnSparseS(Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
267 Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
268 9, binsTrack,xminTrack, xmaxTrack);
269 fHistodEdxTot[idedx] = new THnSparseS(Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
270 Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
271 9, binsTrack,xminTrack, xmaxTrack);
276 for (Int_t ivar=0;ivar<6;ivar++){
277 for (Int_t ivar2=0;ivar2<9;ivar2++){
278 fHistoDelta[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
279 fHistoDelta[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
280 fHistoPull[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
281 fHistoPull[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
282 BinLogX(fHistoDelta[ivar],7);
283 BinLogX(fHistoPull[ivar],7);
285 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
286 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
287 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
288 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
289 BinLogX(fHistodEdxMax[ivar],7);
290 BinLogX(fHistodEdxTot[ivar],7);
297 void AliTPCcalibCosmic::Add(const AliTPCcalibCosmic* cosmic){
301 for (Int_t ivar=0; ivar<6;ivar++){
302 if (fHistoDelta[ivar] && cosmic->fHistoDelta[ivar]){
303 fHistoDelta[ivar]->Add(cosmic->fHistoDelta[ivar]);
305 if (fHistoPull[ivar] && cosmic->fHistoPull[ivar]){
306 fHistoPull[ivar]->Add(cosmic->fHistoPull[ivar]);
309 for (Int_t ivar=0; ivar<4;ivar++){
310 if (fHistodEdxMax[ivar] && cosmic->fHistodEdxMax[ivar]){
311 fHistodEdxMax[ivar]->Add(cosmic->fHistodEdxMax[ivar]);
313 if (fHistodEdxTot[ivar] && cosmic->fHistodEdxTot[ivar]){
314 fHistodEdxTot[ivar]->Add(cosmic->fHistodEdxTot[ivar]);
322 void AliTPCcalibCosmic::Process(AliESDEvent *event) {
327 Printf("ERROR: ESD not available");
330 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
332 Printf("ERROR: esdFriend not available");
337 FindPairs(event); // nearly everything takes place in find pairs...
339 if (GetDebugLevel()>20) printf("Hallo world: Im here and processing an event\n");
340 Int_t ntracks=event->GetNumberOfTracks();
341 fHistNTracks->Fill(ntracks);
342 if (ntracks==0) return;
343 // AliESDcosmic cosmicESD;
344 // TTreeSRedirector * cstream = GetDebugStreamer();
345 // cosmicESD.SetDebugStreamer(cstream);
346 // cosmicESD.ProcessEvent(event);
347 // if (cstream) cosmicESD.DumpToTree();
353 void AliTPCcalibCosmic::FillHistoPerformance(const AliExternalTrackParam *par0, const AliExternalTrackParam *par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam */*inner1*/, AliTPCseed *seed0, AliTPCseed *seed1, const AliExternalTrackParam *param0Combined ){
355 // par0,par1 - parameter of tracks at DCA 0
356 // inner0,inner1 - parameter of tracks at the TPC entrance
357 // seed0, seed1 - detailed track information
358 // param0Combined - Use combined track parameters for binning
360 Int_t kMinCldEdx =20;
361 Int_t ncl0 = seed0->GetNumberOfClusters();
362 Int_t ncl1 = seed1->GetNumberOfClusters();
364 const Double_t kpullCut = 10;
370 Double_t radius0 = TMath::Sqrt(xyz0[0]*xyz0[0]+xyz0[1]*xyz0[1]);
371 Double_t radius1 = TMath::Sqrt(xyz1[0]*xyz1[0]+xyz1[1]*xyz1[1]);
372 inner0->GetXYZ(xyz0);
373 Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
375 x[1] = TMath::Min(ncl0,ncl1);
376 x[2] = (radius0+radius1)*0.5;
377 x[3] = param0Combined->GetZ();
378 x[4] = inner0->GetSnp();
379 x[5] = param0Combined->GetTgl();
380 x[6] = param0Combined->GetSigned1Pt();
381 x[7] = param0Combined->Pt();
386 delta[0] = (par0->GetY()+par1->GetY());
387 delta[1] = (par0->GetZ()-par1->GetZ());
388 delta[2] = (par0->GetAlpha()-par1->GetAlpha()-TMath::Pi());
389 delta[3] = (par0->GetTgl()+par1->GetTgl());
390 delta[4] = (par0->GetParameter()[4]+par1->GetParameter()[4]);
391 delta[5] = (par0->Pt()-par1->Pt())/((par0->Pt()+par1->Pt())*0.5);
393 sigma[0] = TMath::Sqrt(par0->GetSigmaY2()+par1->GetSigmaY2());
394 sigma[1] = TMath::Sqrt(par0->GetSigmaZ2()+par1->GetSigmaZ2());
395 sigma[2] = TMath::Sqrt(par0->GetSigmaSnp2()+par1->GetSigmaSnp2());
396 sigma[3] = TMath::Sqrt(par0->GetSigmaTgl2()+par1->GetSigmaTgl2());
397 sigma[4] = TMath::Sqrt(par0->GetSigma1Pt2()+par1->GetSigma1Pt2());
398 sigma[5] = sigma[4]*((par0->Pt()+par1->Pt())*0.5);
401 for (Int_t ivar=0;ivar<6;ivar++){
402 if (sigma[ivar]==0) isOK=kFALSE;
403 x[0]= delta[ivar]/sigma[ivar];
404 if (TMath::Abs(x[0])>kpullCut) isOK = kFALSE;
408 if (isOK) for (Int_t ivar=0;ivar<6;ivar++){
409 x[0]= delta[ivar]/TMath::Sqrt(2);
410 if (ivar==2 || ivar ==3) x[0]*=1000;
411 fHistoDelta[ivar]->Fill(x);
413 x[0]= delta[ivar]/sigma[ivar];
414 fHistoPull[ivar]->Fill(x);
419 // Fill dedx performance
421 for (Int_t ipad=0; ipad<4;ipad++){
427 if (ipad==0) row1=63;
428 if (ipad==1) {row0=63; row1=63+64;}
429 if (ipad==2) {row0=128;}
430 Int_t nclUp = TMath::Nint(seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
431 Int_t nclDown = TMath::Nint(seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
432 Int_t minCl = TMath::Min(nclUp,nclDown);
433 if (minCl<kMinCldEdx) continue;
436 Float_t dEdxTotUp = seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
437 Float_t dEdxTotDown = seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
438 Float_t dEdxMaxUp = seed0->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
439 Float_t dEdxMaxDown = seed1->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
441 if (dEdxTotDown<=0) continue;
442 if (dEdxMaxDown<=0) continue;
443 x[0]=dEdxTotUp/dEdxTotDown;
444 fHistodEdxTot[ipad]->Fill(x);
445 x[0]=dEdxMaxUp/dEdxMaxDown;
446 fHistodEdxMax[ipad]->Fill(x);
453 void AliTPCcalibCosmic::MaterialBudgetDump(AliExternalTrackParam *const par0, AliExternalTrackParam *const par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam *inner1, AliTPCseed *const seed0, AliTPCseed *const seed1, AliExternalTrackParam *const param0Combined, AliExternalTrackParam *const param1Combined){
455 // matrial budget AOD dump
457 // par0,par1 - parameter of tracks at DCA 0
458 // inner0,inner1 - parameter of tracks at the TPC entrance
459 // seed0, seed1 - detailed track information
460 // param0Combined - Use combined track parameters for binning
462 Double_t p0In = inner0->GetP();
463 Double_t p1In = inner1->GetP();
464 Double_t p0V = par0->GetP();
465 Double_t p1V = par1->GetP();
467 Double_t pt0In = inner0->Pt();
468 Double_t pt1In = inner1->Pt();
469 Double_t pt0V = par0->Pt();
470 Double_t pt1V = par1->Pt();
471 Int_t ncl0 = seed0->GetNumberOfClusters();
472 Int_t ncl1 = seed1->GetNumberOfClusters();
473 Int_t nclmin=TMath::Min(ncl0,ncl1);
474 Double_t sign = (param0Combined->GetSigned1Pt()>0) ? 1:-1.;
476 TTreeSRedirector * pcstream = GetDebugStreamer();
478 (*pcstream)<<"material"<<
479 "run="<<fRun<< // run number
480 "event="<<fEvent<< // event number
481 "time="<<fTime<< // time stamp of event
482 "trigger="<<fTrigger<< // trigger
483 "triggerClass="<<&fTriggerClass<< // trigger
484 "mag="<<fMagF<< // magnetic field
485 "sign="<<sign<< // sign of the track
501 "up0.="<<param0Combined<<
502 "up1.="<<param1Combined<<
508 void AliTPCcalibCosmic::Analyze() {
510 fMIPvalue = CalculateMIPvalue(fDeDxMIP);
518 void AliTPCcalibCosmic::FindPairs(AliESDEvent *event) {
522 // Track0 is choosen in upper TPC part
523 // Track1 is choosen in lower TPC part
525 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
526 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
527 Int_t ntracks=event->GetNumberOfTracks();
528 TObjArray tpcSeeds(ntracks);
529 if (ntracks==0) return;
530 Double_t vtxx[3]={0,0,0};
531 Double_t svtxx[3]={0.000001,0.000001,100.};
532 AliESDVertex vtx(vtxx,svtxx);
536 for (Int_t i=0;i<ntracks;++i) {
537 AliESDtrack *track = event->GetTrack(i);
538 fClusters->Fill(track->GetTPCNcls());
540 const AliExternalTrackParam * trackIn = track->GetInnerParam();
541 const AliExternalTrackParam * trackOut = track->GetOuterParam();
542 if (!trackIn) continue;
543 if (!trackOut) continue;
544 if (ntracks>4 && TMath::Abs(trackIn->GetTgl())<0.0015) continue; // filter laser
547 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
548 if (!friendTrack) continue;
549 TObject *calibObject;
550 AliTPCseed *seed = 0;
551 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
552 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
554 if (seed) tpcSeeds.AddAt(seed,i);
556 Double_t meanP = 0.5*(trackIn->GetP() + trackOut->GetP());
557 if (seed && track->GetTPCNcls() > 80 + 60/(1+TMath::Exp(-meanP+5))) {
558 fDeDx->Fill(meanP, seed->CookdEdxNorm(0.0,0.45,0,0,159));
560 if (meanP > 0.4 && meanP < 0.45) fDeDxMIP->Fill(seed->CookdEdxNorm(0.0,0.45,0,0,159));
562 // if (GetDebugLevel()>0&&meanP>0.2&&seed->CookdEdxNorm(0.0,0.45,0,0,159)>300) {
563 // //TFile *curfile = AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile();
564 // //if (curfile) printf(">>> p+ in file: %s \t event: %i \t Number of ESD tracks: %i \n", curfile->GetName(), (int)event->GetEventNumberInFile(), (int)ntracks);
565 // // if (track->GetOuterParam()->GetAlpha()<0) cout << " Polartiy: " << track->GetSign() << endl;
572 if (ntracks<2) return;
576 for (Int_t i=0;i<ntracks;++i) {
577 AliESDtrack *track0 = event->GetTrack(i);
578 // track0 - choosen upper part
579 if (!track0) continue;
580 if (!track0->GetOuterParam()) continue;
581 if (track0->GetOuterParam()->GetAlpha()<0) continue;
583 track0->GetDirection(dir0);
584 for (Int_t j=0;j<ntracks;++j) {
586 AliESDtrack *track1 = event->GetTrack(j);
588 if (!track1) continue;
589 if (!track1->GetOuterParam()) continue;
590 if (track1->GetOuterParam()->GetAlpha()>0) continue;
593 track1->GetDirection(dir1);
595 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
596 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
597 if (! seed0) continue;
598 if (! seed1) continue;
599 Float_t dedx0 = seed0->CookdEdxNorm(0.05,0.55,0,0,159);
600 Float_t dedx1 = seed1->CookdEdxNorm(0.05,0.55,0,0,159);
602 Float_t dedx0I = seed0->CookdEdxNorm(0.05,0.55,0,0,63);
603 Float_t dedx1I = seed1->CookdEdxNorm(0.05,0.55,0,0,63);
605 Float_t dedx0O = seed0->CookdEdxNorm(0.05,0.55,0,64,159);
606 Float_t dedx1O = seed1->CookdEdxNorm(0.05,0.55,0,64,159);
608 Float_t dir = (dir0[0]*dir1[0] + dir0[1]*dir1[1] + dir0[2]*dir1[2]);
609 Float_t d0 = track0->GetLinearD(0,0);
610 Float_t d1 = track1->GetLinearD(0,0);
612 // conservative cuts - convergence to be guarantied
613 // applying before track propagation
614 if (TMath::Abs(d0+d1)>fCutMaxD) continue; // distance to the 0,0
615 if (dir>fCutMinDir) continue; // direction vector product
616 Float_t bz = AliTracker::GetBz();
617 Float_t dvertex0[2]; //distance to 0,0
618 Float_t dvertex1[2]; //distance to 0,0
619 track0->GetDZ(0,0,0,bz,dvertex0);
620 track1->GetDZ(0,0,0,bz,dvertex1);
621 if (TMath::Abs(dvertex0[1])>250) continue;
622 if (TMath::Abs(dvertex1[1])>250) continue;
626 Float_t dmax = TMath::Max(TMath::Abs(d0),TMath::Abs(d1));
627 AliExternalTrackParam param0(*track0);
628 AliExternalTrackParam param1(*track1);
630 // Propagate using Magnetic field and correct fo material budget
634 Double_t maxsnp=0.90;
635 AliTracker::PropagateTrackToBxByBz(¶m0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign0);
636 AliTracker::PropagateTrackToBxByBz(¶m1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign1);
638 // Propagate rest to the 0,0 DCA - z should be ignored
640 Bool_t b0 = param0.PropagateToDCA(&vtx,bz,1000);
641 Bool_t b1 = param1.PropagateToDCA(&vtx,bz,1000);
643 param0.GetDZ(0,0,0,bz,dvertex0);
644 param1.GetDZ(0,0,0,bz,dvertex1);
645 if (TMath::Abs(param0.GetZ()-param1.GetZ())>fCutMaxDz) continue;
647 Double_t xyz0[3];//,pxyz0[3];
648 Double_t xyz1[3];//,pxyz1[3];
651 Bool_t isPair = IsPair(¶m0,¶m1);
653 if (isPair) FillAcordeHist(track0);
655 // combined track params
657 AliExternalTrackParam *par0U=MakeCombinedTrack(¶m0,¶m1);
658 AliExternalTrackParam *par1U=MakeCombinedTrack(¶m1,¶m0);
663 TTreeSRedirector * cstream = GetDebugStreamer();
664 //printf("My stream=%p\n",(void*)cstream);
665 AliExternalTrackParam *ip0 = (AliExternalTrackParam *)track0->GetInnerParam();
666 AliExternalTrackParam *ip1 = (AliExternalTrackParam *)track1->GetInnerParam();
667 AliExternalTrackParam *op0 = (AliExternalTrackParam *)track0->GetOuterParam();
668 AliExternalTrackParam *op1 = (AliExternalTrackParam *)track1->GetOuterParam();
669 Bool_t isCrossI = ip0->GetZ()*ip1->GetZ()<0;
670 Bool_t isCrossO = op0->GetZ()*op1->GetZ()<0;
671 Double_t alpha0 = TMath::ATan2(dir0[1],dir0[0]);
672 Double_t alpha1 = TMath::ATan2(dir1[1],dir1[0]);
676 FillHistoPerformance(¶m0, ¶m1, ip0, ip1, seed0, seed1,par0U);
677 MaterialBudgetDump(¶m0, ¶m1, ip0, ip1, seed0, seed1,par0U,par1U);
679 (*cstream) << "Track0" <<
680 "run="<<fRun<< // run number
681 "event="<<fEvent<< // event number
682 "time="<<fTime<< // time stamp of event
683 "trigger="<<fTrigger<< // trigger
684 "triggerClass="<<&fTriggerClass<< // trigger
685 "mag="<<fMagF<< // magnetic field
686 "dir="<<dir<< // direction
687 "OK="<<isPair<< // will be accepted
688 "b0="<<b0<< // propagate status
689 "b1="<<b1<< // propagate status
690 "crossI="<<isCrossI<< // cross inner
691 "crossO="<<isCrossO<< // cross outer
693 "Orig0.=" << track0 << // original track 0
694 "Orig1.=" << track1 << // original track 1
695 "Tr0.="<<¶m0<< // track propagated to the DCA 0,0
696 "Tr1.="<<¶m1<< // track propagated to the DCA 0,0
697 "Ip0.="<<ip0<< // inner param - upper
698 "Ip1.="<<ip1<< // inner param - lower
699 "Op0.="<<op0<< // outer param - upper
700 "Op1.="<<op1<< // outer param - lower
701 "Up0.="<<par0U<< // combined track 0
702 "Up1.="<<par1U<< // combined track 1
704 "v00="<<dvertex0[0]<< // distance using kalman
705 "v01="<<dvertex0[1]<< //
706 "v10="<<dvertex1[0]<< //
707 "v11="<<dvertex1[1]<< //
708 "d0="<<d0<< // linear distance to 0,0
709 "d1="<<d1<< // linear distance to 0,0
713 "x00="<<xyz0[0]<< // global position close to vertex
717 "x10="<<xyz1[0]<< // global position close to vertex
723 "dir00="<<dir0[0]<< // direction upper
727 "dir10="<<dir1[0]<< // direction lower
732 "Seed0.=" << seed0 << // original seed 0
733 "Seed1.=" << seed1 << // original seed 1
735 "dedx0="<<dedx0<< // dedx0 - all
736 "dedx1="<<dedx1<< // dedx1 - all
738 "dedx0I="<<dedx0I<< // dedx0 - inner ROC
739 "dedx1I="<<dedx1I<< // dedx1 - inner ROC
741 "dedx0O="<<dedx0O<< // dedx0 - outer ROC
742 "dedx1O="<<dedx1O<< // dedx1 - outer ROC
755 void AliTPCcalibCosmic::FillAcordeHist(AliESDtrack *upperTrack) {
757 // Pt cut to select straight tracks which can be easily propagated to ACORDE which is outside the magnetic field
758 if (upperTrack->Pt() < 10 || upperTrack->GetTPCNcls() < 80) return;
760 const Double_t acordePlane = 850.; // distance of the central Acorde detectors to the beam line at y =0
761 const Double_t roof = 210.5; // distance from x =0 to end of magnet roof
764 upperTrack->GetXYZ(r);
766 upperTrack->GetDirection(d);
768 z = r[2] + (d[2]/d[1])*(acordePlane - r[1]);
769 x = r[0] + (d[0]/d[1])*(acordePlane - r[1]);
772 x = r[0] + (d[0]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
773 z = r[2] + (d[2]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
776 x = r[0] + (d[0]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);
777 z = r[2] + (d[2]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);
780 fModules->Fill(z, x);
786 Long64_t AliTPCcalibCosmic::Merge(TCollection *const li) {
788 TIterator* iter = li->MakeIterator();
789 AliTPCcalibCosmic* cal = 0;
791 while ((cal = (AliTPCcalibCosmic*)iter->Next())) {
792 if (!cal->InheritsFrom(AliTPCcalibCosmic::Class())) {
793 //Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
797 fHistNTracks->Add(cal->GetHistNTracks());
798 fClusters->Add(cal-> GetHistClusters());
799 fModules->Add(cal->GetHistAcorde());
800 fHistPt->Add(cal->GetHistPt());
801 fDeDx->Add(cal->GetHistDeDx());
802 fDeDxMIP->Add(cal->GetHistMIP());
810 Bool_t AliTPCcalibCosmic::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
814 // 0. Same direction - OPOSITE - cutDir +cutT
815 TCut cutDir("cutDir","dir<-0.99")
817 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
820 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
824 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
827 const Double_t *p0 = tr0->GetParameter();
828 const Double_t *p1 = tr1->GetParameter();
829 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
830 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
831 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
833 Double_t d0[3], d1[3];
834 tr0->GetDirection(d0);
835 tr1->GetDirection(d1);
836 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
843 Double_t AliTPCcalibCosmic::CalculateMIPvalue(TH1F * hist) {
845 TF1 * funcDoubleGaus = new TF1("funcDoubleGaus", "gaus(0)+gaus(3)",0,1000);
846 funcDoubleGaus->SetParameters(hist->GetEntries()*0.75,hist->GetMean()/1.3,hist->GetMean()*0.10,
847 hist->GetEntries()*0.25,hist->GetMean()*1.3,hist->GetMean()*0.10);
848 hist->Fit(funcDoubleGaus);
849 Double_t aMIPvalue = TMath::Min(funcDoubleGaus->GetParameter(1),funcDoubleGaus->GetParameter(4));
851 delete funcDoubleGaus;
860 void AliTPCcalibCosmic::CalculateBetheParams(TH2F */*hist*/, Double_t * /*initialParam*/) {
862 // Not implemented yet
869 void AliTPCcalibCosmic::BinLogX(THnSparse *const h, Int_t axisDim) {
871 // Method for the correct logarithmic binning of histograms
873 TAxis *axis = h->GetAxis(axisDim);
874 int bins = axis->GetNbins();
876 Double_t from = axis->GetXmin();
877 Double_t to = axis->GetXmax();
878 Double_t *newBins = new Double_t[bins + 1];
881 Double_t factor = pow(to/from, 1./bins);
883 for (int i = 1; i <= bins; i++) {
884 newBins[i] = factor * newBins[i-1];
886 axis->Set(bins, newBins);
892 void AliTPCcalibCosmic::BinLogX(TH1 *const h) {
894 // Method for the correct logarithmic binning of histograms
896 TAxis *axis = h->GetXaxis();
897 int bins = axis->GetNbins();
899 Double_t from = axis->GetXmin();
900 Double_t to = axis->GetXmax();
901 Double_t *newBins = new Double_t[bins + 1];
904 Double_t factor = pow(to/from, 1./bins);
906 for (int i = 1; i <= bins; i++) {
907 newBins[i] = factor * newBins[i-1];
909 axis->Set(bins, newBins);
915 AliExternalTrackParam *AliTPCcalibCosmic::MakeTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
919 AliExternalTrackParam *par1R= new AliExternalTrackParam(*track1);
920 par1R->Rotate(track0->GetAlpha());
921 par1R->PropagateTo(track0->GetX(),AliTracker::GetBz());
924 Double_t * param = (Double_t*)par1R->GetParameter();
925 Double_t * covar = (Double_t*)par1R->GetCovariance();
933 covar[6] *=-1.; covar[7] *=-1.; covar[8] *=-1.;
934 //covar[10]*=-1.; covar[11]*=-1.; covar[12]*=-1.;
939 AliExternalTrackParam *AliTPCcalibCosmic::MakeCombinedTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
941 // Make combined track
944 AliExternalTrackParam * par1T = MakeTrack(track0,track1);
945 AliExternalTrackParam * par0U = new AliExternalTrackParam(*track0);
947 UpdateTrack(*par0U,*par1T);
953 void AliTPCcalibCosmic::UpdateTrack(AliExternalTrackParam &track1, const AliExternalTrackParam &track2){
955 // Update track 1 with track 2
959 TMatrixD vecXk(5,1); // X vector
960 TMatrixD covXk(5,5); // X covariance
961 TMatrixD matHk(5,5); // vector to mesurement
962 TMatrixD measR(5,5); // measurement error
963 TMatrixD vecZk(5,1); // measurement
965 TMatrixD vecYk(5,1); // Innovation or measurement residual
966 TMatrixD matHkT(5,5);
967 TMatrixD matSk(5,5); // Innovation (or residual) covariance
968 TMatrixD matKk(5,5); // Optimal Kalman gain
969 TMatrixD mat1(5,5); // update covariance matrix
970 TMatrixD covXk2(5,5); //
971 TMatrixD covOut(5,5);
973 Double_t *param1=(Double_t*) track1.GetParameter();
974 Double_t *covar1=(Double_t*) track1.GetCovariance();
975 Double_t *param2=(Double_t*) track2.GetParameter();
976 Double_t *covar2=(Double_t*) track2.GetCovariance();
978 // copy data to the matrix
979 for (Int_t ipar=0; ipar<5; ipar++){
980 for (Int_t jpar=0; jpar<5; jpar++){
981 covXk(ipar,jpar) = covar1[track1.GetIndex(ipar, jpar)];
982 measR(ipar,jpar) = covar2[track2.GetIndex(ipar, jpar)];
986 vecXk(ipar,0) = param1[ipar];
987 vecZk(ipar,0) = param2[ipar];
996 vecYk = vecZk-matHk*vecXk; // Innovation or measurement residual
997 matHkT=matHk.T(); matHk.T();
998 matSk = (matHk*(covXk*matHkT))+measR; // Innovation (or residual) covariance
1000 matKk = (covXk*matHkT)*matSk; // Optimal Kalman gain
1001 vecXk += matKk*vecYk; // updated vector
1002 covXk2 = (mat1-(matKk*matHk));
1003 covOut = covXk2*covXk;
1007 // copy from matrix to parameters
1020 for (Int_t ipar=0; ipar<5; ipar++){
1021 param1[ipar]= vecXk(ipar,0) ;
1022 for (Int_t jpar=0; jpar<5; jpar++){
1023 covar1[track1.GetIndex(ipar, jpar)]=covOut(ipar,jpar);