]>
Commit | Line | Data |
---|---|---|
7f6ddf58 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
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 | **************************************************************************/ | |
15 | ||
16 | //------------------------------------------------------------------------- | |
17 | // Implementation of the V0 vertexer class | |
18 | // | |
19 | // Origin: Iouri Belikov, IReS, Strasbourg, Jouri.Belikov@cern.ch | |
20 | //------------------------------------------------------------------------- | |
116cbefd | 21 | #include <Riostream.h> |
116cbefd | 22 | #include <TPDGCode.h> |
7f6ddf58 | 23 | #include <TObjArray.h> |
116cbefd | 24 | #include <TTree.h> |
7f6ddf58 | 25 | |
116cbefd | 26 | #include "AliITStrackV2.h" |
7f6ddf58 | 27 | #include "AliV0vertex.h" |
28 | #include "AliV0vertexer.h" | |
7f6ddf58 | 29 | |
30 | ClassImp(AliV0vertexer) | |
31 | ||
566abf75 | 32 | Int_t AliV0vertexer::Tracks2V0vertices(TTree *tTree, TTree *vTree) { |
7f6ddf58 | 33 | //-------------------------------------------------------------------- |
34 | //This function reconstructs V0 vertices | |
35 | //-------------------------------------------------------------------- | |
566abf75 | 36 | TBranch *branch=tTree->GetBranch("tracks"); |
37 | if (!branch) { | |
38 | Error("Tracks2V0vertices","Can't get the branch !"); | |
39 | return 1; | |
7f6ddf58 | 40 | } |
566abf75 | 41 | Int_t nentr=(Int_t)tTree->GetEntries(); |
7f6ddf58 | 42 | |
43 | TObjArray negtrks(nentr/2); | |
44 | TObjArray postrks(nentr/2); | |
45 | ||
46 | Int_t nneg=0, npos=0, nvtx=0; | |
47 | ||
48 | Int_t i; | |
49 | for (i=0; i<nentr; i++) { | |
50 | AliITStrackV2 *iotrack=new AliITStrackV2; | |
51 | branch->SetAddress(&iotrack); | |
566abf75 | 52 | tTree->GetEvent(i); |
4ab260d6 | 53 | |
54 | iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.); | |
55 | ||
7f6ddf58 | 56 | if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);} |
57 | else {npos++; postrks.AddLast(iotrack);} | |
58 | } | |
59 | ||
60 | ||
7f6ddf58 | 61 | AliV0vertex *ioVertex=0; |
566abf75 | 62 | branch=vTree->GetBranch("vertices"); |
63 | if (!branch) vTree->Branch("vertices","AliV0vertex",&ioVertex,32000,3); | |
64 | else branch->SetAddress(&ioVertex); | |
7f6ddf58 | 65 | |
66 | ||
67 | for (i=0; i<nneg; i++) { | |
68 | if (i%10==0) cerr<<nneg-i<<'\r'; | |
69 | AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i); | |
a9a2d814 | 70 | |
5d390eda | 71 | if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue; |
72 | if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue; | |
4ab260d6 | 73 | |
7f6ddf58 | 74 | for (Int_t k=0; k<npos; k++) { |
75 | AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k); | |
76 | ||
5d390eda | 77 | if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue; |
78 | if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue; | |
4ab260d6 | 79 | |
5d390eda | 80 | if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin) |
81 | if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue; | |
4ab260d6 | 82 | |
7f6ddf58 | 83 | |
84 | AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt; | |
85 | ||
86 | Double_t dca=PropagateToDCA(pnt,ppt); | |
87 | if (dca > fDCAmax) continue; | |
88 | ||
89 | AliV0vertex vertex(*pnt,*ppt); | |
90 | if (vertex.GetChi2() > fChi2max) continue; | |
91 | ||
92 | /* Think of something better here ! | |
93 | nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue; | |
94 | pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue; | |
95 | */ | |
96 | ||
97 | Double_t x,y,z; vertex.GetXYZ(x,y,z); | |
98 | Double_t r2=x*x + y*y; | |
99 | if (r2 > fRmax*fRmax) continue; | |
100 | if (r2 < fRmin*fRmin) continue; | |
101 | ||
102 | Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz); | |
103 | Double_t p2=px*px+py*py+pz*pz; | |
5d390eda | 104 | Double_t cost=((x-fX)*px + (y-fY)*py + (z-fZ)*pz)/ |
105 | TMath::Sqrt(p2*((x-fX)*(x-fX) + (y-fY)*(y-fY) + (z-fZ)*(z-fZ))); | |
4ab260d6 | 106 | |
5d390eda | 107 | //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue; |
108 | if (cost < fCPAmax) continue; | |
7f6ddf58 | 109 | |
110 | //vertex.ChangeMassHypothesis(); //default is Lambda0 | |
111 | ||
566abf75 | 112 | ioVertex=&vertex; vTree->Fill(); |
7f6ddf58 | 113 | |
114 | nvtx++; | |
115 | } | |
116 | } | |
117 | ||
118 | cerr<<"Number of reconstructed V0 vertices: "<<nvtx<<endl; | |
119 | ||
7f6ddf58 | 120 | negtrks.Delete(); |
121 | postrks.Delete(); | |
122 | ||
7f6ddf58 | 123 | return 0; |
124 | } | |
125 | ||
7f6ddf58 | 126 | Double_t AliV0vertexer::PropagateToDCA(AliITStrackV2 *n, AliITStrackV2 *p) { |
127 | //-------------------------------------------------------------------- | |
128 | // This function returns the DCA between two tracks | |
129 | // The tracks will be moved to the point of DCA ! | |
130 | //-------------------------------------------------------------------- | |
49a7a79a | 131 | return n->PropagateToDCA(p); |
7f6ddf58 | 132 | } |
133 | ||
134 | ||
135 | ||
136 | ||
137 | ||
138 | ||
139 | ||
140 | ||
141 | ||
142 | ||
143 | ||
144 | ||
145 | ||
146 | ||
147 | ||
148 |