]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliVertexerTracks.cxx
Bugfix. Default argument was taken (J.Belikov)
[u/mrichter/AliRoot.git] / STEER / AliVertexerTracks.cxx
CommitLineData
2d57349e 1/**************************************************************************
2 * Copyright(c) 1998-2003, 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 vertexer from ESD tracks
18//
19// Origin: AliITSVertexerTracks
20// A.Dainese, Padova,
21// andrea.dainese@pd.infn.it
22// M.Masera, Torino,
23// massimo.masera@to.infn.it
24// Moved to STEER and adapted to ESD tracks:
25// F.Prino, Torino, prino@to.infn.it
26//-----------------------------------------------------------------
27
28//---- Root headers --------
29#include <TTree.h>
30//---- AliRoot headers -----
31#include "AliStrLine.h"
32#include "AliVertexerTracks.h"
c5e3e5d1 33#include "AliESD.h"
2d57349e 34#include "AliESDtrack.h"
35
36ClassImp(AliVertexerTracks)
37
38
39//----------------------------------------------------------------------------
40AliVertexerTracks::AliVertexerTracks():
41 TObject(),fVert()
42{
43//
44// Default constructor
45//
46 SetVtxStart();
47 SetMinTracks();
48 fDCAcut=0;
49 fAlgo=1;
50}
51//-----------------------------------------------------------------------------
52AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
53 TObject(),fVert()
54{
55//
56// Standard constructor
57//
58 SetVtxStart(xStart,yStart);
59 SetMinTracks();
60 fDCAcut=0;
61 fAlgo=1;
62}
63//-----------------------------------------------------------------------------
64AliVertexerTracks::~AliVertexerTracks() {
65 // Default Destructor
66 // The objects poited by the following pointers are not owned
67 // by this class and are not deleted
68}
69//----------------------------------------------------------------------------
70Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree) {
71//
72// Propagate tracks to initial vertex position and store them in a TObjArray
73//
74 Double_t alpha,xlStart;
75 Int_t nTrks = 0;
76
77 Int_t nEntries = (Int_t)trkTree.GetEntries();
78 if(!fTrkArray.IsEmpty()) fTrkArray.Clear();
79 fTrkArray.Expand(nEntries);
80
81 for(Int_t i=0; i<nEntries; i++) {
82 AliESDtrack *track = new AliESDtrack;
83 trkTree.SetBranchAddress("tracks",&track);
84 trkTree.GetEvent(i);
85
86 // propagate track to vtxSeed
87 alpha = track->GetAlpha();
88 xlStart = fNominalPos[0]*TMath::Cos(alpha)+fNominalPos[1]*TMath::Sin(alpha);
89 track->PropagateTo(xlStart,GetField()); // to vtxSeed
90 fTrkArray.AddLast(track);
91 nTrks++;
92 }
93 return nTrks;
94}
95//----------------------------------------------------------------------------
96AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) {
97//
98// Return vertex from tracks in trkTree
99//
100
101 // get tracks and propagate them to initial vertex position
102 Int_t nTrks = PrepareTracks(*trkTree);
103 if(nTrks < fMinTracks) {
104 printf("TooFewTracks\n");
105 Double_t vtx[3]={0,0,0};
106 fVert.SetXYZ(vtx);
107 fVert.SetDispersion(999);
108 fVert.SetNContributors(-5);
109 return &fVert;
110 }
111
112 // Set initial vertex position from ESD
113 if(fAlgo==1) StrLinVertexFinderMinDist(1);
114 if(fAlgo==2) StrLinVertexFinderMinDist(0);
115 if(fAlgo==3) HelixVertexFinder();
116 if(fAlgo==4) VertexFinder(1);
117 if(fAlgo==5) VertexFinder(0);
118 return &fVert;
119}
c5e3e5d1 120
121//----------------------------------------------------------------------------
122AliESDVertex* AliVertexerTracks::FindVertex(const AliESD *event) {
123//
124// This is a simple wrapping (by Jouri.Belikov@cern.ch) over the original
125// code by the authors of this class.
126//
127 Int_t nt=event->GetNumberOfTracks(), nacc=0;
128 while (nt--) {
129 AliESDtrack *t=event->GetTrack(nt);
130 if ((t->GetStatus()&AliESDtrack::kITSrefit)==0) continue;
131 fTrkArray.AddLast(t);
132 nacc++;
133 }
134
135 // get tracks and propagate them to initial vertex position
136 if(nacc < fMinTracks) {
137 printf("TooFewTracks\n");
138 Double_t vtx[3]={0,0,0};
139 fVert.SetXYZ(vtx);
140 fVert.SetDispersion(999);
141 fVert.SetNContributors(-5);
142 } else
143 switch (fAlgo) {
144 case 1: StrLinVertexFinderMinDist(1); break;
145 case 2: StrLinVertexFinderMinDist(0); break;
146 case 3: HelixVertexFinder(); break;
147 case 4: VertexFinder(1); break;
148 case 5: VertexFinder(0); break;
149 default: printf("Wrong algorithm\n"); break;
150 }
151
152 fTrkArray.Clear();
153 return &fVert;
154}
155
156
2d57349e 157//----------------------------------------------------------------------------
158AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) {
159//
160// Return vertex from array of tracks
161//
162
163 // get tracks and propagate them to initial vertex position
164 Int_t nTrks = trkArray->GetEntriesFast();
165 if(nTrks < fMinTracks) {
166 printf("TooFewTracks\n");
167 Double_t vtx[3]={0,0,0};
168 fVert.SetXYZ(vtx);
169 fVert.SetDispersion(999);
170 fVert.SetNContributors(-5);
171 return &fVert;
172 }
173 TTree *trkTree = new TTree("TreeT","tracks");
174 AliESDtrack *esdTrack = 0;
175 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
176 for(Int_t i=0; i<nTrks; i++){
177 esdTrack = (AliESDtrack*)trkArray->At(i);
178 trkTree->Fill();
179 }
180
2e5c8fc7 181 AliVertex *vtx = VertexForSelectedTracks(trkTree);
182 delete trkTree;
183 return vtx;
2d57349e 184}
185
186//---------------------------------------------------------------------------
187void AliVertexerTracks::VertexFinder(Int_t optUseWeights) {
188
189 // Get estimate of vertex position in (x,y) from tracks DCA
190
191 Double_t initPos[3];
192 initPos[2] = 0.;
193 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
194 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
195 Double_t aver[3]={0.,0.,0.};
196 Double_t aversq[3]={0.,0.,0.};
197 Double_t sigmasq[3]={0.,0.,0.};
198 Double_t sigma=0;
199 Int_t ncombi = 0;
200 AliESDtrack *track1;
201 AliESDtrack *track2;
202 Double_t alpha,mindist;
203 Double_t field=GetField();
204
205 for(Int_t i=0; i<nacc; i++){
206 track1 = (AliESDtrack*)fTrkArray.At(i);
207 alpha=track1->GetAlpha();
208 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
08df6187 209
210 Double_t pos[3]; track1->GetXYZAt(mindist,field,pos);
211 Double_t dir[3]; track1->GetPxPyPzAt(mindist,field,dir);
212 AliStrLine *line1 = new AliStrLine(pos,dir);
2d57349e 213
214 for(Int_t j=i+1; j<nacc; j++){
215 track2 = (AliESDtrack*)fTrkArray.At(j);
216 alpha=track2->GetAlpha();
217 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
08df6187 218
219 Double_t pos[3]; track2->GetXYZAt(mindist,field,pos);
220 Double_t dir[3]; track2->GetPxPyPzAt(mindist,field,dir);
221 AliStrLine *line2 = new AliStrLine(pos,dir);
222
2d57349e 223 Double_t distCA=line2->GetDCA(line1);
224 if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
225 Double_t pnt1[3],pnt2[3],crosspoint[3];
226
227 if(optUseWeights<=0){
228 Int_t retcode = line2->Cross(line1,crosspoint);
229 if(retcode>=0){
230 ncombi++;
231 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
232 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
233 }
234 }
235 if(optUseWeights>0){
236 Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
237 if(retcode>=0){
238 Double_t alpha, cs, sn;
239 alpha=track1->GetAlpha();
240 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
241 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
242 alpha=track2->GetAlpha();
243 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
244 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
245 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
246 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
247 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
248 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
249 crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0];
250 crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1];
251 crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
252
253 ncombi++;
254 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
255 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
256 }
257 }
258 }
259 delete line2;
260 }
261 delete line1;
262 }
263 if(ncombi>0){
264 for(Int_t jj=0;jj<3;jj++){
265 initPos[jj] = aver[jj]/ncombi;
266 aversq[jj]/=ncombi;
267 sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
268 sigma+=sigmasq[jj];
269 }
270 sigma=TMath::Sqrt(TMath::Abs(sigma));
271 }
272 else {
273 Warning("VertexFinder","Finder did not succed");
274 sigma=999;
275 }
276 fVert.SetXYZ(initPos);
277 fVert.SetDispersion(sigma);
278 fVert.SetNContributors(ncombi);
279}
280//---------------------------------------------------------------------------
281void AliVertexerTracks::HelixVertexFinder() {
282
283 // Get estimate of vertex position in (x,y) from tracks DCA
284
285
286 Double_t initPos[3];
287 initPos[2] = 0.;
288 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
289 Double_t field=GetField();
290
291 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
292
293 Double_t aver[3]={0.,0.,0.};
294 Double_t averquad[3]={0.,0.,0.};
295 Double_t sigmaquad[3]={0.,0.,0.};
296 Double_t sigma=0;
297 Int_t ncombi = 0;
298 AliESDtrack *track1;
299 AliESDtrack *track2;
300 Double_t distCA;
301 Double_t x, par[5];
302 Double_t alpha, cs, sn;
303 Double_t crosspoint[3];
304 for(Int_t i=0; i<nacc; i++){
305 track1 = (AliESDtrack*)fTrkArray.At(i);
306
307
308 for(Int_t j=i+1; j<nacc; j++){
309 track2 = (AliESDtrack*)fTrkArray.At(j);
310
311 distCA=track2->PropagateToDCA(track1,field);
312
313 if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
314 track1->GetExternalParameters(x,par);
315 alpha=track1->GetAlpha();
316 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
317 Double_t x1=x*cs - par[0]*sn;
318 Double_t y1=x*sn + par[0]*cs;
319 Double_t z1=par[1];
320 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
321 track2->GetExternalParameters(x,par);
322 alpha=track2->GetAlpha();
323 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
324 Double_t x2=x*cs - par[0]*sn;
325 Double_t y2=x*sn + par[0]*cs;
326 Double_t z2=par[1];
327 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
328 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
329 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
330 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
331 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
332 crosspoint[0]=wx1*x1 + wx2*x2;
333 crosspoint[1]=wy1*y1 + wy2*y2;
334 crosspoint[2]=wz1*z1 + wz2*z2;
335
336 ncombi++;
337 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
338 for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
339 }
340 }
341
342 }
343 if(ncombi>0){
344 for(Int_t jj=0;jj<3;jj++){
345 initPos[jj] = aver[jj]/ncombi;
346 averquad[jj]/=ncombi;
347 sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
348 sigma+=sigmaquad[jj];
349 }
350 sigma=TMath::Sqrt(TMath::Abs(sigma));
351 }
352 else {
353 Warning("HelixVertexFinder","Finder did not succed");
354 sigma=999;
355 }
356 fVert.SetXYZ(initPos);
357 fVert.SetDispersion(sigma);
358 fVert.SetNContributors(ncombi);
359}
360//---------------------------------------------------------------------------
361void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights){
362
363 // Calculate the point at minimum distance to prepared tracks
364
365 Double_t initPos[3];
366 initPos[2] = 0.;
367 Double_t sigma=0;
368 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
369 const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
370 Double_t field=GetField();
371
372 AliESDtrack *track1;
373 Double_t (*vectP0)[3]=new Double_t [knacc][3];
374 Double_t (*vectP1)[3]=new Double_t [knacc][3];
375
376 Double_t sum[3][3];
377 Double_t dsum[3]={0,0,0};
378 for(Int_t i=0;i<3;i++)
379 for(Int_t j=0;j<3;j++)sum[i][j]=0;
380 for(Int_t i=0; i<knacc; i++){
381 track1 = (AliESDtrack*)fTrkArray.At(i);
382 Double_t alpha=track1->GetAlpha();
383 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
08df6187 384
385 Double_t pos[3]; track1->GetXYZAt(mindist,field,pos);
7f101aa0 386 Double_t dir[3]; track1->GetPxPyPzAt(mindist,field,dir);
08df6187 387 AliStrLine *line1 = new AliStrLine(pos,dir);
2d57349e 388
389 Double_t p0[3],cd[3];
390 line1->GetP0(p0);
391 line1->GetCd(cd);
392 Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
393 vectP0[i][0]=p0[0];
394 vectP0[i][1]=p0[1];
395 vectP0[i][2]=p0[2];
396 vectP1[i][0]=p1[0];
397 vectP1[i][1]=p1[1];
398 vectP1[i][2]=p1[2];
399
400 Double_t matr[3][3];
401 Double_t dknow[3];
402 if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
403 if(optUseWeights==1){
404 Double_t sigmasq[3];
405 sigmasq[0]=track1->GetSigmaY2();
406 sigmasq[1]=track1->GetSigmaY2();
407 sigmasq[2]=track1->GetSigmaZ2();
408 GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
409 }
410
411 for(Int_t iii=0;iii<3;iii++){
412 dsum[iii]+=dknow[iii];
413 for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
414 }
415 delete line1;
416 }
417
418 Double_t vett[3][3];
419 Double_t det=GetDeterminant3X3(sum);
420
421 if(det!=0){
422 for(Int_t zz=0;zz<3;zz++){
423 for(Int_t ww=0;ww<3;ww++){
424 for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
425 }
426 for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
427 initPos[zz]=GetDeterminant3X3(vett)/det;
428 }
429
430
431 for(Int_t i=0; i<knacc; i++){
432 Double_t p0[3]={0,0,0},p1[3]={0,0,0};
433 for(Int_t ii=0;ii<3;ii++){
434 p0[ii]=vectP0[i][ii];
435 p1[ii]=vectP1[i][ii];
436 }
437 sigma+=GetStrLinMinDist(p0,p1,initPos);
438 }
439
440 sigma=TMath::Sqrt(sigma);
441 }else{
442 Warning("StrLinVertexFinderMinDist","Finder did not succed");
443 sigma=999;
444 }
445 delete vectP0;
446 delete vectP1;
447 fVert.SetXYZ(initPos);
448 fVert.SetDispersion(sigma);
449 fVert.SetNContributors(knacc);
450}
451//_______________________________________________________________________
452Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3]){
453 //
454 Double_t det=matr[0][0]*matr[1][1]*matr[2][2]-matr[0][0]*matr[1][2]*matr[2][1]-matr[0][1]*matr[1][0]*matr[2][2]+matr[0][1]*matr[1][2]*matr[2][0]+matr[0][2]*matr[1][0]*matr[2][1]-matr[0][2]*matr[1][1]*matr[2][0];
455 return det;
456}
457//____________________________________________________________________________
458void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t m[][3],Double_t *d){
459
460 //
461 Double_t x12=p0[0]-p1[0];
462 Double_t y12=p0[1]-p1[1];
463 Double_t z12=p0[2]-p1[2];
464 Double_t kk=x12*x12+y12*y12+z12*z12;
465 m[0][0]=2-2/kk*x12*x12;
466 m[0][1]=-2/kk*x12*y12;
467 m[0][2]=-2/kk*x12*z12;
468 m[1][0]=-2/kk*x12*y12;
469 m[1][1]=2-2/kk*y12*y12;
470 m[1][2]=-2/kk*y12*z12;
471 m[2][0]=-2/kk*x12*z12;
472 m[2][1]=-2*y12*z12;
473 m[2][2]=2-2/kk*z12*z12;
474 d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
475 d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
476 d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
477
478}
479//____________________________________________________________________________
480void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t m[][3],Double_t *d){
481 //
482 Double_t x12=p1[0]-p0[0];
483 Double_t y12=p1[1]-p0[1];
484 Double_t z12=p1[2]-p0[2];
485
486 Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
487
488 Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
489
490 Double_t cc[3];
491 cc[0]=-x12/sigmasq[0];
492 cc[1]=-y12/sigmasq[1];
493 cc[2]=-z12/sigmasq[2];
494
495 Double_t ww=(-p0[0]*x12*sigmasq[1]*sigmasq[2]-p0[1]*y12*sigmasq[0]*sigmasq[2]-p0[2]*z12*sigmasq[0]*sigmasq[1])/den;
496
497 Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
498
499 Double_t aa[3];
500 aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
501 aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
502 aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
503
504 m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
505 m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
506 m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
507
508 m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
509 m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
510 m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
511
512 m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
513 m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
514 m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
515
516 d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
517 d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
518 d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
519
520}
521//_____________________________________________________________________________
522Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0){
523 //
524 Double_t x12=p0[0]-p1[0];
525 Double_t y12=p0[1]-p1[1];
526 Double_t z12=p0[2]-p1[2];
527 Double_t x10=p0[0]-x0[0];
528 Double_t y10=p0[1]-x0[1];
529 Double_t z10=p0[2]-x0[2];
530 return ((x10*x10+y10*y10+z10*z10)*(x12*x12+y12*y12+z12*z12)-(x10*x12+y10*y12+z10*z12)*(x10*x12+y10*y12+z10*z12))/(x12*x12+y12*y12+z12*z12);
531}