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. *
16 **************************************************************************/
20 // This class performs a fast fit of helices going through the <=6 *
21 // points of the ITS, with the goal of studying tracking and *
22 // vertexing performances. *
23 // Generated kinematics is used to take into account different weights *
24 // associated to points in different layers (with different multiple *
25 // scattering-originated errors). *
27 // Based on the work by A. Strandlie, R. Fruhwirth *
29 // First implementation by N. Bustreo, R. Turrisi - July 2000 *
31 // Further modifications by A. Dainese, R. Turrisi *
33 // Contact: Rosario Turrisi, rosario.turrisi@pd.infn.it *
35 // ************************************************************************
38 // Modified November, 7th 2001 by Rosario Turrisi
39 // (rosario.turrisi@pd.infn.it)
41 // FitHelix returns different values. 0=ok, >0 =problem
42 // void FitLinear -> Int_t FitLinear to give feedback of errors to FitHelix
45 // Modified July, 30th 2001 by Rosario Turrisi
46 // (rosario.turrisi@pd.infn.it)
48 // Fit for z now in (z,s) plane.
49 // Returns parameters in order to write the helix equation
50 // and find the right phase/initial point.
52 // "PROPER WEIGHTS": (1+R^2)^2/(\sigma_x^2 + \sigma_y^2 + \sigma_MS^2)
57 #include <Riostream.h>
59 #include "AliITSRiemannFit.h"
61 #include "TClonesArray.h"
64 #include "Riostream.h"
67 #include "TGraphErrors.h"
69 #include "TParticle.h"
72 #include "AliITSRecPoint.h"
73 #include "AliITSgeom.h"
74 #include "AliITSmodule.h"
78 ClassImp(AliITSRiemannFit)
81 AliITSRiemannFit::AliITSRiemannFit() {
82 ///////////////////////////////////////////////////////////
83 // Default constructor.
84 // Set everything to zero.
85 ////////////////////////////////////////////////////////////
95 for(Int_t i=0;i<6;i++)fPLay[i] = 0;
99 //______________________________________________________________________
100 AliITSRiemannFit::AliITSRiemannFit(const AliITSRiemannFit &rf) : TObject(rf) {
102 // Copies are not allowed. The method is protected to avoid misuse.
103 Error("AliITSRiemannFit","Copy constructor not allowed\n");
106 //______________________________________________________________________
107 AliITSRiemannFit& AliITSRiemannFit::operator=(const AliITSRiemannFit& /* rf */){
108 // Assignment operator
109 // Assignment is not allowed. The method is protected to avoid misuse.
110 Error("= operator","Assignment operator not allowed\n");
114 //______________________________________________________________________
115 AliITSRiemannFit::~AliITSRiemannFit() {
116 ///////////////////////////////////////////////////////////
117 // Default destructor.
118 // if arrays exist delete them. Then set everything to zero.
119 ////////////////////////////////////////////////////////////
121 for(Int_t i=0;i<fSizeEvent;i++) delete[] fPointRecs[i];
123 } // end if fPointRecs!=0
132 for(Int_t i=0;i<6;i++)fPLay[i] = 0;
135 //----------------------------------------------------------------------
137 AliITSRiemannFit::AliITSRiemannFit(Int_t size,Int_t ntracks) {
138 ///////////////////////////////////////////////////////////
140 // Set fSizeEvent to size and fPrimaryTracks to ntracks.
142 ////////////////////////////////////////////////////////////
146 fPrimaryTracks = ntracks;
151 AliPointtl *first = new AliPointtl[fSizeEvent];
152 AliPointtl **pointRecs = new AliPointtl*[fSizeEvent];
153 for(Int_t i=0;i<6;i++)fPLay[i] = 0;
154 for(Int_t j=0;j<fSizeEvent;j++) // create an array of struct
155 pointRecs[j] = &(first[j]);
158 // ---------------------------------------------------------------------
159 AliITSRiemannFit::AliPointtl::AliPointtl(){
160 // default constructor
183 // ---------------------------------------------------------------------
185 void FillPoints(AliITSRiemannFit::AliPointtl **Points,Int_t &index,Float_t *xpoint,
187 TLorentzVector pE,TLorentzVector oT,Int_t *id,
188 Int_t track, Char_t *name,Int_t code,
190 ///////////////////////////////////////////////////////////////////////
191 // Fill the structure AliPointtl with the proper data
193 //////////////////////////////////////////////////////////////////////
194 Float_t pPI2 = 2.0*TMath::Pi();
202 phi = TMath::ATan2(y,x);
203 if(phi<0.0) phi += pPI2;
204 Points[i]->SetPhi(phi);
205 Points[i]->SetEta(-0.5*tan(0.5*TMath::ATan2(r,z)));
209 Points[i]->SetdX(error[0]);
210 Points[i]->SetdY(error[1]);
211 Points[i]->SetdZ(error[2]);
213 Points[i]->SetTrack(track);
214 Points[i]->SetLay(id[0]);
215 Points[i]->SetLad(id[1]);
216 Points[i]->SetDet(id[2]);
217 Points[i]->SetMomentum(&pE);
218 Points[i]->SetOrigin(&oT);
219 Points[i]->SetPt(sqrt(pE.X()*pE.X()+pE.Y()*pE.Y()));
220 Points[i]->SetCode(code);
221 Points[i]->SetName(name);
222 Points[i]->SetVertexPhi(phiorigin);
227 // -----------------------------------------------------------------------
229 void AliITSRiemannFit::InitPoints(Int_t ntracks,AliITS *ITS,
230 TTree *TR,Int_t nparticles){
231 //////////////////////////////////////////////////////////////////////
232 // Fill the class member fPointRecs with the reconstructed points
233 // Set All other members to the real values
235 /////////////////////////////////////////////////////////////////////
236 printf("\n ************* Starting Init Points *************\n");
238 AliITSgeom *gm = (AliITSgeom*)ITS->GetITSgeom();
239 //get pointer to modules array
240 TObjArray *iTSmodules = ITS->GetModules();
241 Int_t nmodules=iTSmodules->GetEntriesFast();
242 printf("nmodules = %d \n",nmodules);
243 // Get the points from points file
244 AliITSmodule *itsModule;
247 AliITSRecPoint *recp;
248 nent=TR->GetEntries();
249 TClonesArray *iTSrec = ITS->RecPoints();
252 for (mod=0; mod<nmodules; mod++) {
253 itsModule=(AliITSmodule*)iTSmodules->At(mod);
254 ITS->ResetRecPoints();
256 Int_t nrecp = iTSrec->GetEntries();
262 fPrimaryTracks = ntracks;
263 fParticles = nparticles;
264 AliITSRiemannFit::AliPointtl *global = new AliPointtl[iMAX];
265 fPointRecs = new AliITSRiemannFit::AliPointtl*[iMAX];
267 for(Int_t j=0;j<iMAX;j++) {
268 fPointRecs[j] = &(global[j]);
271 Int_t ieta=0,ieta2=0;
272 Int_t i,id[4],idold[4];
273 Int_t track=0;// // track of hit
274 Float_t xpoint[3],errorPlus[3],errorMinus[3],globalError[3]; // position and error of the point
275 TLorentzVector oT,pE;
276 Float_t locals[3],localserror[3],localsplus[3],localsminus[3]; // local position and local errors
280 Int_t layer,ladder,detector;
281 Float_t xcluster,zcluster;
282 Int_t num=0,nspdi=0,nspdo=0,nsddi=0,nsddo=0,nssdi=0,nssdo=0;
284 for (mod=0; mod<nmodules; mod++) {
285 itsModule=(AliITSmodule*)iTSmodules->At(mod);
286 ITS->ResetRecPoints();
288 Int_t nrecp = iTSrec->GetEntries();
289 if (!nrecp) continue;
290 itsModule->GetID(layer,ladder,detector);
292 for (irec=0;irec<nrecp;irec++) {
293 recp = (AliITSRecPoint*)iTSrec->UncheckedAt(irec);
294 track=recp->fTracks[0];
295 if(track <0 ) continue;
296 xcluster=recp->GetX(); // x on cluster
297 zcluster=recp->GetZ(); // z on cluster
298 part = (TParticle*) gAlice->GetMCApp()->Particle(track);
299 part->ProductionVertex(oT); // set the vertex
300 part->Momentum(pE); // set the vertex momentum
301 name = part->GetName();
303 sprintf(nam2,"%s",name);
304 code = part->GetPdgCode();
310 locals[0]=xcluster; // x on cluster
311 locals[1]=0.0; // y on cluster
312 locals[2]=zcluster; // z on cluster
313 localserror[0]=sqrt(recp->GetSigmaX2());
315 localserror[2]=sqrt(recp->GetSigmaZ2());
316 localsplus[0]=xcluster+sqrt(recp->GetSigmaX2()); // x on cluster
317 if(layer==1||layer==2) localsplus[1]=0.0150/2; // y on cluster
318 else if(layer==3||layer==4) localsplus[1]=0.0280/2; // y on cluster
319 else if(layer==5||layer==6) localsplus[1]=0.0300/2; // y on cluster
320 localsplus[2]=zcluster+sqrt(recp->GetSigmaZ2()); // z on cluster
321 localsminus[0]=xcluster-sqrt(recp->GetSigmaX2()); // x on cluster
322 localsminus[1]=0.0; // y on cluster
323 localsminus[2]=zcluster-sqrt(recp->GetSigmaZ2()); // z on cluster
325 gm->LtoG(layer,ladder,detector,locals,xpoint);
326 gm->LtoG(layer,ladder,detector,localsplus,errorPlus);
327 gm->LtoG(layer,ladder,detector,localsminus,errorMinus);
328 globalError[0]=0.5*TMath::Abs(errorPlus[0]-errorMinus[0]);
329 globalError[1]=0.5*TMath::Abs(errorPlus[1]-errorMinus[1]);
330 globalError[2]=0.5*TMath::Abs(errorPlus[2]-errorMinus[2]);
332 if(TMath::Abs(part->Eta())<=1.0) ieta++;
333 if(TMath::Abs(part->Eta())<=0.5) ieta2++;
335 if(!(id[0]==idold[0]&&id[1]==idold[1]&&
336 id[2]==idold[2]&&id[3]==idold[3])) {
337 FillPoints(fPointRecs,num,xpoint,globalError,pE,oT,id,track,nam2,code,pPhi);
361 // FillPoints(fspdi,nspdi,xpoint,globalError,pE,oT,id,track,name,code,pPhi);
365 // FillPoints(fspdo,nspdo,xpoint,globalError,pE,oT,id,track,name,code,pPhi);
379 for(i=0;i<4;i++) idold[i] = id[i];
380 for(i=0;i<3;i++) xpoint[i] = 0.0;
381 } // end if id != idold
393 printf("%d primary tracks in eta=+-1\n",ieta);
394 printf("%d primary tracks#2 in eta=+-0.5\n",ieta2);
395 printf("\nInitPoints :\n\nPoints on Layer1 : %d on Layer2 : %d\n",nspdi,nspdo);
396 printf("Points on Layer3 : %d on Layer4 : %d\n",nsddi,nsddo);
397 printf("Points on Layer5 : %d on Layer6 : %d\n",nssdi,nssdo);
398 printf("Points on all Layers: %d\n",num);
399 printf("\n ************* Init Points Finished *************\n");
402 // ------------------------------------------------------------------------
403 ///////////////////////////////////////////////////////////
404 // Functions for sorting the fPointRecs array
405 ///////////////////////////////////////////////////////////
406 Bool_t SortZ(const AliITSRiemannFit::AliPointtl *s1,const AliITSRiemannFit::AliPointtl *s2){
407 // Z sorting function for qsort.
410 a = s1->GetZ() - s2->GetZ();
411 if(a<0.0) return kTRUE;
412 if(a>0.0) return kFALSE;
415 Bool_t SortTrack(const AliITSRiemannFit::AliPointtl *s1,const AliITSRiemannFit::AliPointtl *s2){
416 // track sorting function for qsort.
419 a = s1->GetTrack() - s2->GetTrack();
420 if(a<0.0) return kTRUE;
421 if(a>0.0) return kFALSE;
424 void hpsortTrack(AliITSRiemannFit::AliPointtl **ra,Int_t n){
426 AliITSRiemannFit::AliPointtl *rra;
430 l = ((n-1) >> 1) +1; // divide 2 + 1
434 rra = ra[--l]; // decrement first
438 if(--ir == 0){ // decrement first
446 if( j<ir && SortTrack(ra[j],ra[j+1]) ) j++;
447 if( SortTrack(rra,ra[j]) ){
458 void hpsortZ(AliITSRiemannFit::AliPointtl **ra,Int_t n){
460 AliITSRiemannFit::AliPointtl *rra;
464 l = ((n-1) >> 1) +1; // devide 2 + 1
468 rra = ra[--l]; // decrament first
472 if(--ir == 0){ // decrament first
480 if( j<ir && SortZ(ra[j],ra[j+1]) ) j++;
481 if( SortZ(rra,ra[j]) ){
492 //-----------------------------------------------------------------------
493 ////////////////////////////////////////////////////////////////////
495 ///////////////////////////////////////////////////////////////////
496 Int_t Partition(Int_t array[],Int_t left,Int_t right){
497 Int_t val = array[left];
499 Int_t rm = right + 1;
510 Int_t tempr = array[rm];
521 ///////////////////////////////////////////////////////////////////////
523 void AliITSRiemannFit::WritePoints(void) {
524 /////////////////////////////////////////////////////////////////////
525 // write the data in a file (temporary ascii)
526 /////////////////////////////////////////////////////////////////////
527 FILE *ascii= fopen("AsciiPoints.dat","w");
528 for(Int_t i=0;i<fPoints;i++) {
529 fprintf(ascii,"%d\t%d\t%f\t%f\t%f\n",fPointRecs[i]->GetLay(),
530 fPointRecs[i]->GetTrack(),fPointRecs[i]->GetX(),
531 fPointRecs[i]->GetY(),fPointRecs[i]->GetZ());
536 //-----------------------------------------------------------------------
538 void AliITSRiemannFit::ReadPoints(void) {
539 //////////////////////////////////////////////////////////////////////
540 // read the filled array
541 /////////////////////////////////////////////////////////////////////
542 hpsortTrack(fPointRecs,fPoints);
543 for(Int_t i=0;i<fPoints;i++)
544 printf("%d\t%d\t%d\t%f\t%f\t%f\t(%.0f,%.0f,%.0f)\t%.3f\t%s\n",
545 i,fPointRecs[i]->GetLay(),fPointRecs[i]->GetTrack(),
546 fPointRecs[i]->GetX(),fPointRecs[i]->GetY(),
547 fPointRecs[i]->GetZ(),fPointRecs[i]->GetOrigin()->X(),
548 fPointRecs[i]->GetOrigin()->Y(),fPointRecs[i]->GetOrigin()->Z(),
549 fPointRecs[i]->GetPt(),fPointRecs[i]->GetName());
552 //-----------------------------------------------------------------------
554 Int_t AliITSRiemannFit::SolveCubic(Double_t a,Double_t b,Double_t c,
555 Double_t &x1,Double_t &x2,Double_t &x3){
556 //////////////////////////////////////////////
557 /// Solve cubic equation:
558 /// x^3 + a*x^2 +b*x + c
560 /// returns x1 , x2 , x3
561 ////////////////////////////////////////
563 Double_t qQ = ((a*a - 3*b)/9);
564 Double_t rR = ((2*a*a*a - 9*a*b +27*c)/54);
566 Double_t aF = -2*sqrt(qQ);
568 Double_t pPI2 = TMath::Pi()*2;
570 if( rR*rR>qQ*qQ*qQ ) {
571 cout<<"\nTrack "<<"Determinant :\n\t\t No Real Solutions !!!\n"<<endl;
578 theta = TMath::ACos(rR/sqrt(qQ*qQ*qQ));
580 x1 = (aF*TMath::Cos(theta/3))-g;
581 x2 = (aF*TMath::Cos((theta+pPI2)/3))-g;
582 x3 = (aF*TMath::Cos((theta-pPI2)/3))-g;
586 //-----------------------------------------------------------------
588 void RiemannTransf(Int_t npoints,TVector3 **From,TVector3 **To) {
589 ///////////////////////////////////////////////////////////////////////
590 // This function apllies the transformation in the Riemann sphere
592 ///////////////////////////////////////////////////////////////////////
593 Float_t *rR = new Float_t[npoints];
594 Float_t *theta = new Float_t[npoints];
595 Float_t pPI2 = 2*TMath::Pi();
598 for(Int_t i=0;i<npoints;i++) {
599 rR[i] = sqrt(From[i]->X()*From[i]->X()+From[i]->Y()*From[i]->Y());
600 theta[i] = TMath::ATan2(From[i]->Y(),From[i]->X());
601 if(theta[i]<0) theta[i]+=pPI2;
602 x = rR[i]*cos(theta[i])/(1+rR[i]*rR[i]);
603 y = rR[i]*sin(theta[i])/(1+rR[i]*rR[i]);
604 z = rR[i]*rR[i]/(1+rR[i]*rR[i]);
605 To[i]->SetXYZ(x,y,z);
613 //---------------------------------------------------------------------
615 Int_t FitLinear(Int_t npoints, TVector3 **input, TVector3 **errors, Double_t omega,
616 Double_t &thu0, Double_t &thv0, Double_t &phi, TVector2 &zData, TVector3 &zError,
618 ///////////////////////////////////////////////////////////////////////
619 // Fit the points in the (z,s) plane - helix 3rd equation
621 ///////////////////////////////////////////////////////////////////////
623 //PH Double_t z[npoints],x[npoints],y[npoints],s[npoints];
624 //PH Double_t ez[npoints],ex[npoints],ey[npoints],es[npoints];
625 Double_t * z = new Double_t[npoints];
626 Double_t * x = new Double_t[npoints];
627 Double_t * y = new Double_t[npoints];
628 Double_t * s = new Double_t[npoints];
629 Double_t * ez = new Double_t[npoints];
630 Double_t * ex = new Double_t[npoints];
631 Double_t * ey = new Double_t[npoints];
632 Double_t * es = new Double_t[npoints];
633 Double_t z0=0.0,vpar=0.0,ez0=0.0,evpar=0.0, chisquare;
635 // Double_t chi=TMath::Pi()/2.0+phi;
636 Double_t chi=-TMath::Pi()-phi;
637 Double_t angold=0.0, tpang=0.0;
638 for(Int_t k = 0; k<npoints; k++) {
639 x[k] = 10.0*input[k]->X(); ex[k] = 10.0*errors[k]->X();
640 y[k] = 10.0*input[k]->Y(); ey[k] = 10.0*errors[k]->Y();
641 z[k] = 10.0*input[k]->Z(); ez[k] = 10.0*errors[k]->Z();
642 if(TMath::Abs(x[k]-thu0)<1.0e-5) { // should never happen, nor give troubles...
644 cerr<<"limit for x-x_0 "<<x[k]<<" "<<thu0<<endl;
655 Double_t ang1=TMath::ATan2((y[k]-thv0),(x[k]-thu0));
656 if( (x[k]-thu0)<0 ) {
658 tpang=ang1-TMath::Sign(TMath::Pi()*2.0,ang1);
663 if (k>0) direction+=(z[k]>z[k-1] ? 1 : -1);
664 s[k] = (ang1+chi)/omega;
665 es[k]=TMath::Sqrt(ey[k]*ey[k]+ex[k]*ex[k]/TMath::Power((x[k]-thu0),4))*TMath::Abs(s[k]);
667 if ( TMath::Abs(direction) != (npoints-1) ) {return 11;}
669 TGraphErrors *fitHist = new TGraphErrors(npoints,s,z,es,ez);
670 fitHist->Fit("pol1","qQ");
671 z0 = fitHist->GetFunction("pol1")->GetParameter(0);
672 vpar = fitHist->GetFunction("pol1")->GetParameter(1);
673 ez0 = fitHist->GetFunction("pol1")->GetParError(0);
674 evpar = fitHist->GetFunction("pol1")->GetParError(1);
675 chisquare = fitHist->GetFunction("pol1")->GetChisquare();
677 zError.SetXYZ(ez0,evpar,chisquare);
685 for(Int_t j = 0; j < npoints; j++) {
690 avs /= (Double_t)npoints;
691 avz /= (Double_t)npoints;
692 avsz /= (Double_t)npoints;
694 for(Int_t l = 0; l < npoints; l++) {
695 sigmas += (s[l]-avs)*(s[l]-avs);
696 sigmaz += (z[l]-avz)*(z[l]-avz);
698 sigmas /=(Double_t)npoints;
699 sigmaz /=(Double_t)npoints;
701 sigmas = sqrt(sigmas);
702 sigmaz = sqrt(sigmaz);
704 corrLin = (avsz-avs*avz)/(sigmas*sigmaz);
718 //-------------------------------------------------------------------
719 Int_t AliITSRiemannFit::FitHelix(Int_t tracknumber,Double_t Px,Double_t Py,Double_t Pz,Double_t& fd0,
720 Double_t& fphi,Double_t& u0, Double_t& v0, Double_t& rho,Double_t& omega, Double_t& z0,
721 Double_t& vpar,Double_t& chisql, Double_t& fCorrLin,Double_t& fFit,
722 Int_t first,Int_t second,Int_t third,Int_t fourth,Int_t fifth,Int_t sixth) {
723 ///////////////////////////////////////////////////////////////////////
724 // This function finds the helix paramenters
725 // d0 = impact parameter
726 // rho = radius of circle
729 // starting from the momentum and the outcome of
730 // the fit on the Riemann sphere (i.e. u0,v0,rho)
732 // MIND !!!! Here we assume both angular velocities be 1.0 (yes, one-dot-zero !)
735 ///////////////////////////////////////////////////////////////////////
737 // All this stuff relies on this hypothesis !!!
739 // FILE *pout=fopen("chisql.dat","a");
740 Int_t ierr = 0, ierrl=0;
743 Int_t bitlay[6]={1,1,1,1,1,1};
744 bitlay[0]*=first; bitlay[1]*=second; bitlay[2]*=third; bitlay[3]*=fourth; bitlay[4]*=fifth; bitlay[5]*=sixth;
745 fd0 = -9999; // No phisycs value
746 u0 = -9999.9999; // parameters of helix - strange value...
747 v0 = -9999.9999; // parameters of helix - strange value...
748 rho = -9999.9999; // parameters of helix -unphysical strange value...
750 const Char_t* name = 0;
754 Int_t npl[6]={0,0,0,0,0,0};
755 Double_t pP = sqrt(Px*Px+Py*Py+Pz*Pz);
756 Double_t pt = sqrt(Px*Px+Py*Py);
760 TVector3 *ori = new TVector3[iMAX];
761 TVector3 **original = new TVector3*[iMAX];
762 TVector3 *rie = new TVector3[iMAX];
763 TVector3 **riemann = new TVector3*[iMAX];
764 TVector3 *err = new TVector3[iMAX];
765 TVector3 **errors = new TVector3*[iMAX];
766 TVector3 *linerr = new TVector3[iMAX];
767 TVector3 **linerrors = new TVector3*[iMAX];
768 //PH Double_t weight[iMAX];
769 Double_t * weight = new Double_t[iMAX];
772 original[i] = &(ori[i]);
773 riemann[i] = &(rie[i]);
774 errors[i] = &(err[i]);
775 linerrors[i] = &(linerr[i]);
777 for(k =0;k<iMAX;k++) original[k]->SetXYZ(9999,9999,9999);
778 Double_t a11,a12,a13,a21,a22,a23,a31,a32,a33;
779 a11=0;a12=0;a13=0;a21=0;a22=0;a23=0;a31=0;a32=0;a33=0;
783 Double_t a,b,c,d; // cubic parameters
784 Double_t roots[3]= {0.0,0.0,0.0}; // cubic solutions
785 Double_t value = 0.0; // minimum eigenvalue
786 Double_t x1,x2,x3; // eigenvector component
787 Double_t n1,n2,n3,nr= 0;// unit eigenvector
788 Double_t radiusdm[7] = {0.3,0.4,0.7,1.49,2.38,3.91,4.36}; // beam pipe and layers radii [dm]
789 Double_t sigmaMS = 0;
790 TVector3 vVec,vVecNor;
792 // Select RecPoints belonging to the track
793 for(k =0;k<fPoints;k++){
794 if(fPointRecs[k]->GetTrack()==tracknumber) {
795 name = fPointRecs[k]->GetName();
796 pt = fPointRecs[k]->GetPt();
797 pLayer = fPointRecs[k]->GetLay();
798 Int_t ilay = pLayer-1;
799 if(npl[ilay]!=0) continue;
800 if(bitlay[ilay] == 1) {
801 original[nN]->SetXYZ(0.1*fPointRecs[k]->GetX(),0.1*fPointRecs[k]->GetY(),0.1*fPointRecs[k]->GetZ());
802 errors[nN]->SetXYZ(0.1*fPointRecs[k]->GetdX(),0.1*fPointRecs[k]->GetdY(),0.1*fPointRecs[k]->GetdZ());
803 sigmaMS = (radiusdm[pLayer]-radiusdm[0])*0.000724/pP;// beam pipe contribution
804 for(Int_t j=1;j<pLayer;j++) {
805 sigmaMS += (radiusdm[pLayer]-radiusdm[j])*0.00136/pP;
807 weight[nN] = ( 1 + original[nN]->Perp2() )*( 1+ original[nN]->Perp2() )/
808 ( errors[nN]->Perp2() + sigmaMS*sigmaMS );
809 linerrors[nN]->SetXYZ(errors[nN]->X(),errors[nN]->Y(),sqrt(errors[nN]->Z()*errors[nN]->Z()+sigmaMS*sigmaMS));
813 } //end if track==tracknumber
816 // 6 points, no more, no less
818 if(original[5]->X() == 9999 || original[6]->X() != 9999)
821 return 1; // not enough points
827 // FIT ON THE RIEMANN SPHERE FOR (x,y) PLANE
830 RiemannTransf(nN,original,riemann);
832 Double_t sumWeights = 0.0; // sum of weights factor
834 for(Int_t j=0;j<nN;j++){ // mean values for x[i],y[i],z[i]
835 xbar+=weight[j]*riemann[j]->X();
836 ybar+=weight[j]*riemann[j]->Y();
837 zbar+=weight[j]*riemann[j]->Z();
838 sumWeights+=weight[j];
845 for(Int_t j=0;j<nN;j++) { // Calculate the matrix elements
846 a11 += weight[j]*(riemann[j]->X() - xbar)*(riemann[j]->X() - xbar);
847 a12 += weight[j]*(riemann[j]->X() - xbar)*(riemann[j]->Y() - ybar);
848 a22 += weight[j]*(riemann[j]->Y() - ybar)*(riemann[j]->Y() - ybar);
849 a23 += weight[j]*(riemann[j]->Y() - ybar)*(riemann[j]->Z() - zbar);
850 a13 += weight[j]*(riemann[j]->X() - xbar)*(riemann[j]->Z() - zbar);
851 a33 += weight[j]*(riemann[j]->Z() - zbar)*(riemann[j]->Z() - zbar);
864 // ************** Determinant parameters ********************
865 // n.b. simplifications done keeping in mind symmetry of A
869 c = (a11*(a22+a33)+a33*a22-a12*a21-a13*a31-a23*a32);
870 d = (a31*a22*a13+(a12*a21-a11*a22)*a33-2.0*a23*a13*a12+a11*a23*a32);
872 // ************** Find the 3 eigenvalues *************************
873 Int_t checkCubic = SolveCubic(b,c,d,roots[0],roots[1],roots[2]);
876 printf("Track %d Has no real solution continuing ...\n",tracknumber);
881 // **************** Find the lowest eigenvalue *****************
882 if(roots[0]<=roots[1] && roots[0]<=roots[2]) value = roots[0];
883 if(roots[1]<=roots[0] && roots[1]<=roots[2]) value = roots[1];
884 if(roots[2]<=roots[0] && roots[2]<=roots[1]) value = roots[2];
886 // ************ Eigenvector relative to value **************
888 x2 = (a33*a21-a23*a31-value*a21)/(a22*a31-a32*a21-value*a31);
889 x1 = (value-a33-a32*x2)/a31;
890 vVec.SetXYZ(x1,x2,x3);
891 vVecNor = vVec.Unit();
895 nr = -n1*xbar-n2*ybar-n3*zbar;
897 u0 = -0.5*n1/(nr+n3);
898 v0 = -0.5*n2/(nr+n3);
899 rho = sqrt((n1*n1 + n2*n2 -4*nr*(nr+n3))/(4*(nr+n3)*(nr+n3)));
902 fFit += 10.*TMath::Abs(sqrt((original[0]->X()-u0)*(original[0]->X()-u0)+(original[0]->Y()-v0)*(original[0]->Y()-v0))-rho);
903 fFit += 10.*TMath::Abs(sqrt((original[1]->X()-u0)*(original[1]->X()-u0)+(original[1]->Y()-v0)*(original[1]->Y()-v0))-rho);
904 fFit += 10.*TMath::Abs(sqrt((original[2]->X()-u0)*(original[2]->X()-u0)+(original[2]->Y()-v0)*(original[2]->Y()-v0))-rho);
905 fFit += 10.*TMath::Abs(sqrt((original[3]->X()-u0)*(original[3]->X()-u0)+(original[3]->Y()-v0)*(original[3]->Y()-v0))-rho);
906 fFit += 10.*TMath::Abs(sqrt((original[4]->X()-u0)*(original[4]->X()-u0)+(original[4]->Y()-v0)*(original[4]->Y()-v0))-rho);
907 fFit += 10.*TMath::Abs(sqrt((original[5]->X()-u0)*(original[5]->X()-u0)+(original[5]->Y()-v0)*(original[5]->Y()-v0))-rho);
909 fd0 = 100000.*(TMath::Sqrt(u0*u0+v0*v0)-rho); // transverse impact parameter in microns
910 fphi = TMath::ATan2(v0,u0);
912 //**************************************************************************
913 // LINEAR FIT IN (z,s) PLANE: z = zData.X() + zData.Y()*s
914 // strictly linear (no approximation)
915 //**************************************************************************
917 ////////////////////////////////////////////////////////////////////////////////////////////////////////////
919 // REMEMBER, HERE STILL LENGHTS IN DM'S FOR ___INPUT___ BUT zDATA PARAMETERS ARE RETURNED IN CM'S //
920 // rho, u0, v0 parameters converted right now to cm's... it's a mess, I'll take care, sometimes... //
922 ////////////////////////////////////////////////////////////////////////////////////////////////////////////
927 ierrl=FitLinear(nN,original,linerrors,omega,u0,v0,fphi,zData,zError,corrLin);
929 // fprintf(pout,"%f \n",chisql);
933 ierr = (ierrl > ierr ? ierrl : ierr);
938 Int_t AliITSRiemannFit::FitHelix(Int_t NPoints, TVector3** fPointRecs,TVector3** fPointRecErrors,Float_t& f1, Float_t& f2, Float_t& f3) {
940 ///////////////////////////////////////////////////////////////////////
941 // This function finds the helix parameters
942 // d0 = impact parameter
943 // rho = radius of circle
946 // starting from the momentum and the outcome of
947 // the fit on the Riemann sphere (i.e. u0,v0,rho)
949 // MIND !!!! Here we assume both angular velocities be 1.0e-2 (yes, 0.01 !)
952 // Also linear fit in (z,s) is performed, so it's 3-D !
953 // z0 and vpar are calculated (intercept and z-component of velocity, but
954 // in units... you guess.
957 // Values calculated in addition:
959 // - transverse impact parameter fd0
960 // - sum of residuals in (x,y) plane fFit
961 // - chisquare of linear fit chisql
962 // - correlation coefficient fCorrLin
968 ///////////////////////////////////////////////////////////////////////
970 // All this stuff relies on this hypothesis !!!
972 Int_t ierr = 0, ierrl=0;
973 const Double_t kOmega = 1.0e-2;
978 Double_t fd0 = -9999; // fake values
979 Double_t u0 = -9999.9999; // for eventual
980 Double_t v0 = -9999.9999; // debugging
981 Double_t rho = -9999.9999; //
982 Double_t fphi, fFit, chisql, z0, vpar, fCorrLin;
985 // This info is no more there... to be re-considered... maybe
987 // Double_t pP = sqrt(Px*Px+Py*Py+Pz*Pz);
988 // Double_t pt = sqrt(Px*Px+Py*Py);
993 TVector3 *ori = new TVector3[NPoints];
994 TVector3 **original = new TVector3*[NPoints];
995 TVector3 *rie = new TVector3[NPoints];
996 TVector3 **riemann = new TVector3*[NPoints];
997 TVector3 *err = new TVector3[NPoints];
998 TVector3 **errors = new TVector3*[NPoints];
999 TVector3 *linerr = new TVector3[NPoints];
1000 TVector3 **linerrors = new TVector3*[NPoints];
1001 Double_t * weight = new Double_t[NPoints];
1003 for(Int_t i=0; i<NPoints; i++){
1005 original[i] = &(ori[i]);
1006 riemann[i] = &(rie[i]);
1007 errors[i] = &(err[i]);
1008 linerrors[i] = &(linerr[i]);
1010 original[i]->SetXYZ(9999,9999,9999);
1014 // Riemann fit parameters
1016 Double_t a11,a12,a13,a21,a22,a23,a31,a32,a33;
1017 a11=0;a12=0;a13=0;a21=0;a22=0;a23=0;a31=0;a32=0;a33=0;
1022 Double_t a,b,c,d; // cubic parameters
1023 Double_t roots[3]= {0.0,0.0,0.0}; // cubic solutions
1024 Double_t value = 0.0; // minimum eigenvalue
1025 Double_t x1,x2,x3; // eigenvector component
1026 Double_t n1,n2,n3,nr= 0; // unit eigenvector
1027 TVector3 vVec,vVecNor;
1029 for (Int_t ip=0; ip<NPoints; ip++) {
1030 original[ip]->SetXYZ(0.1*fPointRecs[ip]->X(),0.1*fPointRecs[ip]->Y(),0.1*fPointRecs[ip]->Z());
1032 errors[ip]->SetXYZ(0.1*fPointRecErrors[ip]->X(),0.1*fPointRecErrors[ip]->Y(),0.1*fPointRecErrors[ip]->Z());
1033 weight[ip] = (1+original[ip]->Perp2())*(1+original[ip]->Perp2())/(errors[ip]->Perp2());
1034 linerrors[ip]->SetXYZ(errors[ip]->X(),errors[ip]->Y(),errors[ip]->Z());
1040 // FIT ON THE RIEMANN SPHERE FOR (x,y) PLANE
1043 RiemannTransf(NPoints,original,riemann);
1045 Double_t sumWeights = 0.0; // sum of weights factor
1047 for(Int_t j=0;j<NPoints;j++){ // mean values for x[i],y[i],z[i]
1048 xbar+=weight[j]*riemann[j]->X();
1049 ybar+=weight[j]*riemann[j]->Y();
1050 zbar+=weight[j]*riemann[j]->Z();
1051 sumWeights+=weight[j];
1058 for(Int_t j=0;j<NPoints;j++) { // Calculate the matrix elements
1059 a11 += weight[j]*(riemann[j]->X() - xbar)*(riemann[j]->X() - xbar);
1060 a12 += weight[j]*(riemann[j]->X() - xbar)*(riemann[j]->Y() - ybar);
1061 a22 += weight[j]*(riemann[j]->Y() - ybar)*(riemann[j]->Y() - ybar);
1062 a23 += weight[j]*(riemann[j]->Y() - ybar)*(riemann[j]->Z() - zbar);
1063 a13 += weight[j]*(riemann[j]->X() - xbar)*(riemann[j]->Z() - zbar);
1064 a33 += weight[j]*(riemann[j]->Z() - zbar)*(riemann[j]->Z() - zbar);
1067 // this doesn't seem to work...
1069 // a11 /= sumWeights;
1070 // a12 /= sumWeights;
1071 // a22 /= sumWeights;
1072 // a23 /= sumWeights;
1073 // a13 /= sumWeights;
1074 // a33 /= sumWeights;
1086 // ************** Determinant parameters ********************
1087 // n.b. simplifications done keeping in mind symmetry of A
1091 c = (a11*(a22+a33)+a33*a22-a12*a21-a13*a31-a23*a32);
1092 d = (a31*a22*a13+(a12*a21-a11*a22)*a33-2.0*a23*a13*a12+a11*a23*a32);
1094 // ************** Find the 3 eigenvalues *************************
1095 Int_t checkCubic = SolveCubic(b,c,d,roots[0],roots[1],roots[2]);
1097 if(checkCubic !=1 ){
1098 printf("No real solution. Check data.\n");
1103 // **************** Find the lowest eigenvalue *****************
1104 if(roots[0]<=roots[1] && roots[0]<=roots[2]) value = roots[0];
1105 if(roots[1]<=roots[0] && roots[1]<=roots[2]) value = roots[1];
1106 if(roots[2]<=roots[0] && roots[2]<=roots[1]) value = roots[2];
1108 // ************ Eigenvector relative to value **************
1110 x2 = (a33*a21-a23*a31-value*a21)/(a22*a31-a32*a21-value*a31);
1111 x1 = (value-a33-a32*x2)/a31;
1112 vVec.SetXYZ(x1,x2,x3);
1113 vVecNor = vVec.Unit();
1117 nr = -n1*xbar-n2*ybar-n3*zbar;
1119 u0 = -0.5*n1/(nr+n3);
1120 v0 = -0.5*n2/(nr+n3);
1121 rho = sqrt((n1*n1 + n2*n2 -4*nr*(nr+n3))/(4*(nr+n3)*(nr+n3)));
1125 for (Int_t i=0; i<NPoints; i++) {
1126 fFit += 10.*TMath::Abs(sqrt((original[i]->X()-u0)*(original[i]->X()-u0)+(original[i]->Y()-v0)*(original[i]->Y()-v0))-rho);
1128 fd0 = 100000.*(TMath::Sqrt(u0*u0+v0*v0)-rho); // transverse impact parameter in microns
1129 fphi = TMath::ATan2(v0,u0);
1131 //**************************************************************************
1132 // LINEAR FIT IN (z,s) PLANE: z = zData.X() + zData.Y()*s
1133 // strictly linear (no approximation)
1134 //**************************************************************************
1136 ////////////////////////////////////////////////////////////////////////////////////////////////////////////
1138 // REMEMBER, HERE STILL LENGHTS IN DM'S FOR ___INPUT___ BUT zDATA PARAMETERS ARE RETURNED IN CM'S //
1139 // rho, u0, v0 parameters converted right now to cm's... it's a mess, I'll take care, sometimes... //
1141 ////////////////////////////////////////////////////////////////////////////////////////////////////////////
1147 ierrl=LinearFit(NPoints,original,linerrors,kOmega,u0,v0,fphi,zData,zError,corrLin);
1148 if(ierrl==33) return 0;
1150 // fprintf(pout,"%f \n",chisql);
1154 ierr = (ierrl > ierr ? ierrl : ierr);
1159 f2=vpar/(kOmega*TMath::Abs(rho));
1175 //____________________________________________________________
1177 Int_t AliITSRiemannFit::LinearFit(Int_t npoints, TVector3 **input,
1178 TVector3 **errors, Double_t omega,
1179 Double_t &thu0, Double_t &thv0, Double_t &phi,TVector2 &zData, TVector3 &zError,
1181 ///////////////////////////////////////////////////////////////////////
1182 // Fit the points in the (z,s) plane - helix 3rd equation
1184 ///////////////////////////////////////////////////////////////////////
1188 //PH Double_t z[npoints],x[npoints],y[npoints],s[npoints];
1189 //PH Double_t ez[npoints],ex[npoints],ey[npoints],es[npoints];
1190 Double_t * z = new Double_t[npoints];
1191 Double_t * x = new Double_t[npoints];
1192 Double_t * y = new Double_t[npoints];
1193 Double_t * s = new Double_t[npoints];
1194 Double_t * ez = new Double_t[npoints];
1195 Double_t * ex = new Double_t[npoints];
1196 Double_t * ey = new Double_t[npoints];
1197 Double_t * es = new Double_t[npoints];
1198 Double_t z0=0.0,vpar=0.0,ez0=0.0,evpar=0.0, chisquare;
1201 // Double_t chi=TMath::Pi()/2.0+phi;
1202 Double_t chi=-TMath::Pi()-phi;
1203 Double_t angold=0.0, tpang=0.0;
1204 for(Int_t k = 0; k<npoints; k++) {
1205 x[k] = 10.0*input[k]->X(); ex[k] = 10.0*errors[k]->X();
1206 y[k] = 10.0*input[k]->Y(); ey[k] = 10.0*errors[k]->Y();
1207 z[k] = 10.0*input[k]->Z(); ez[k] = 10.0*errors[k]->Z();
1208 if(TMath::Abs(x[k]-thu0)<1.0e-5) { // should never happen, nor give troubles...
1210 cerr<<"limit for x-x_0 "<<x[k]<<" "<<thu0<<endl;
1221 Double_t ang1=TMath::ATan2((y[k]-thv0),(x[k]-thu0));
1222 if( (x[k]-thu0)<0 ) {
1223 if (ang1*angold<0) {
1224 tpang=ang1-TMath::Sign(TMath::Pi()*2.0,ang1);
1229 if (k>0) direction+=(z[k]>z[k-1] ? 1 : -1);
1230 s[k] = (ang1+chi)/omega;
1231 es[k]=TMath::Sqrt(ey[k]*ey[k]+ex[k]*ex[k]/TMath::Power((x[k]-thu0),4))*TMath::Abs(s[k]);
1233 if ( TMath::Abs(direction) != (npoints-1) ) {return 11;}
1235 // if(s[0]>-636 && s[0]<-625) return 33;
1237 TGraph* fitHist = new TGraph(npoints,s,z);
1238 TF1* f1 = new TF1("f1",Fitfunction,-100,100,2);
1240 f1->SetParameter(0,1);
1241 f1->SetParameter(1,1);
1242 f1->SetLineColor(2);
1243 fitHist->Fit(f1,"qQ");
1245 z0 = f1->GetParameter(0);
1246 vpar = f1->GetParameter(1);
1247 ez0 = f1->GetParError(0);
1248 evpar= f1->GetParError(1);
1249 chisquare=f1->GetChisquare();
1251 zError.SetXYZ(ez0,evpar,chisquare);
1259 for(Int_t j = 0; j < npoints; j++) {
1264 avs /= (Double_t)npoints;
1265 avz /= (Double_t)npoints;
1266 avsz /= (Double_t)npoints;
1268 for(Int_t l = 0; l < npoints; l++) {
1269 sigmas += (s[l]-avs)*(s[l]-avs);
1270 sigmaz += (z[l]-avz)*(z[l]-avz);
1272 sigmas /=(Double_t)npoints;
1273 sigmaz /=(Double_t)npoints;
1275 sigmas = sqrt(sigmas);
1276 sigmaz = sqrt(sigmaz);
1278 corrLin = (avsz-avs*avz)/(sigmas*sigmaz);
1290 delete f1; delete fitHist;
1295 //_______________________________________________________
1297 Double_t AliITSRiemannFit::Fitfunction(Double_t *x, Double_t* par){
1298 // function used for fit
1299 return par[0]+(*x)*par[1];