]>
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> |
7f6ddf58 | 22 | #include <TFile.h> |
116cbefd | 23 | #include <TPDGCode.h> |
7f6ddf58 | 24 | #include <TObjArray.h> |
116cbefd | 25 | #include <TTree.h> |
7f6ddf58 | 26 | |
116cbefd | 27 | #include "AliITStrackV2.h" |
7f6ddf58 | 28 | #include "AliV0vertex.h" |
29 | #include "AliV0vertexer.h" | |
7f6ddf58 | 30 | |
31 | ClassImp(AliV0vertexer) | |
32 | ||
33 | Int_t AliV0vertexer::Tracks2V0vertices(const TFile *inp, TFile *out) { | |
34 | //-------------------------------------------------------------------- | |
35 | //This function reconstructs V0 vertices | |
36 | //-------------------------------------------------------------------- | |
37 | TFile *in=(TFile*)inp; | |
38 | TDirectory *savedir=gDirectory; | |
39 | ||
40 | if (!in->IsOpen()) { | |
41 | cerr<<"AliV0vertexer::Tracks2V0vertices(): "; | |
42 | cerr<<"file with ITS tracks has not been open !\n"; | |
43 | return 1; | |
44 | } | |
45 | ||
46 | if (!out->IsOpen()) { | |
47 | cerr<<"AliV0vertexer::Tracks2V0vertices(): "; | |
48 | cerr<<"file for V0 vertices has not been open !\n"; | |
49 | return 2; | |
50 | } | |
51 | ||
52 | in->cd(); | |
53 | ||
5d390eda | 54 | Char_t name[100]; |
55 | sprintf(name,"TreeT_ITS_%d",fEventN); | |
56 | TTree *trkTree=(TTree*)in->Get(name); | |
7f6ddf58 | 57 | TBranch *branch=trkTree->GetBranch("tracks"); |
58 | Int_t nentr=(Int_t)trkTree->GetEntries(); | |
59 | ||
60 | TObjArray negtrks(nentr/2); | |
61 | TObjArray postrks(nentr/2); | |
62 | ||
63 | Int_t nneg=0, npos=0, nvtx=0; | |
64 | ||
65 | Int_t i; | |
66 | for (i=0; i<nentr; i++) { | |
67 | AliITStrackV2 *iotrack=new AliITStrackV2; | |
68 | branch->SetAddress(&iotrack); | |
69 | trkTree->GetEvent(i); | |
4ab260d6 | 70 | |
71 | iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.); | |
72 | ||
7f6ddf58 | 73 | if (iotrack->Get1Pt() > 0.) {nneg++; negtrks.AddLast(iotrack);} |
74 | else {npos++; postrks.AddLast(iotrack);} | |
75 | } | |
76 | ||
77 | ||
78 | out->cd(); | |
5d390eda | 79 | sprintf(name,"TreeV%d",fEventN); |
80 | TTree vtxTree(name,"Tree with V0 vertices"); | |
7f6ddf58 | 81 | AliV0vertex *ioVertex=0; |
82 | vtxTree.Branch("vertices","AliV0vertex",&ioVertex,32000,0); | |
83 | ||
84 | ||
85 | for (i=0; i<nneg; i++) { | |
86 | if (i%10==0) cerr<<nneg-i<<'\r'; | |
87 | AliITStrackV2 *ntrk=(AliITStrackV2 *)negtrks.UncheckedAt(i); | |
a9a2d814 | 88 | |
5d390eda | 89 | if (TMath::Abs(ntrk->GetD(fX,fY))<fDPmin) continue; |
90 | if (TMath::Abs(ntrk->GetD(fX,fY))>fRmax) continue; | |
4ab260d6 | 91 | |
7f6ddf58 | 92 | for (Int_t k=0; k<npos; k++) { |
93 | AliITStrackV2 *ptrk=(AliITStrackV2 *)postrks.UncheckedAt(k); | |
94 | ||
5d390eda | 95 | if (TMath::Abs(ptrk->GetD(fX,fY))<fDPmin) continue; |
96 | if (TMath::Abs(ptrk->GetD(fX,fY))>fRmax) continue; | |
4ab260d6 | 97 | |
5d390eda | 98 | if (TMath::Abs(ntrk->GetD(fX,fY))<fDNmin) |
99 | if (TMath::Abs(ptrk->GetD(fX,fY))<fDNmin) continue; | |
4ab260d6 | 100 | |
7f6ddf58 | 101 | |
102 | AliITStrackV2 nt(*ntrk), pt(*ptrk), *pnt=&nt, *ppt=&pt; | |
103 | ||
104 | Double_t dca=PropagateToDCA(pnt,ppt); | |
105 | if (dca > fDCAmax) continue; | |
106 | ||
107 | AliV0vertex vertex(*pnt,*ppt); | |
108 | if (vertex.GetChi2() > fChi2max) continue; | |
109 | ||
110 | /* Think of something better here ! | |
111 | nt.PropagateToVertex(); if (TMath::Abs(nt.GetZ())<0.04) continue; | |
112 | pt.PropagateToVertex(); if (TMath::Abs(pt.GetZ())<0.04) continue; | |
113 | */ | |
114 | ||
115 | Double_t x,y,z; vertex.GetXYZ(x,y,z); | |
116 | Double_t r2=x*x + y*y; | |
117 | if (r2 > fRmax*fRmax) continue; | |
118 | if (r2 < fRmin*fRmin) continue; | |
119 | ||
120 | Double_t px,py,pz; vertex.GetPxPyPz(px,py,pz); | |
121 | Double_t p2=px*px+py*py+pz*pz; | |
5d390eda | 122 | Double_t cost=((x-fX)*px + (y-fY)*py + (z-fZ)*pz)/ |
123 | TMath::Sqrt(p2*((x-fX)*(x-fX) + (y-fY)*(y-fY) + (z-fZ)*(z-fZ))); | |
4ab260d6 | 124 | |
5d390eda | 125 | //if (cost < (5*fCPAmax-0.9-TMath::Sqrt(r2)*(fCPAmax-1))/4.1) continue; |
126 | if (cost < fCPAmax) continue; | |
7f6ddf58 | 127 | |
128 | //vertex.ChangeMassHypothesis(); //default is Lambda0 | |
129 | ||
130 | ioVertex=&vertex; vtxTree.Fill(); | |
131 | ||
132 | nvtx++; | |
133 | } | |
134 | } | |
135 | ||
136 | cerr<<"Number of reconstructed V0 vertices: "<<nvtx<<endl; | |
137 | ||
138 | vtxTree.Write(); | |
139 | ||
140 | negtrks.Delete(); | |
141 | postrks.Delete(); | |
142 | ||
143 | savedir->cd(); | |
144 | ||
145 | delete trkTree; | |
146 | ||
147 | return 0; | |
148 | } | |
149 | ||
7f6ddf58 | 150 | Double_t AliV0vertexer::PropagateToDCA(AliITStrackV2 *n, AliITStrackV2 *p) { |
151 | //-------------------------------------------------------------------- | |
152 | // This function returns the DCA between two tracks | |
153 | // The tracks will be moved to the point of DCA ! | |
154 | //-------------------------------------------------------------------- | |
49a7a79a | 155 | return n->PropagateToDCA(p); |
7f6ddf58 | 156 | } |
157 | ||
158 | ||
159 | ||
160 | ||
161 | ||
162 | ||
163 | ||
164 | ||
165 | ||
166 | ||
167 | ||
168 | ||
169 | ||
170 | ||
171 | ||
172 |