ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

//-------------------------------------------------------------------------
//               Implementation of the cascade vertexer class
//          Reads V0s and tracks, writes out cascade vertices
//                     Fills the ESD with the cascades 
//    Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
//-------------------------------------------------------------------------

//modified by R. Vernet 30/6/2006 : daughter label
//modified by R. Vernet  3/7/2006 : causality
//modified by I. Belikov 24/11/2006 : static setter for the default cuts

#include "AliESDEvent.h"
#include "AliESDcascade.h"
#include "AliCascadeVertexer.h"

ClassImp(AliCascadeVertexer)

//A set of loose cuts
Double_t 
  AliCascadeVertexer::fgChi2max=33.;   //maximal allowed chi2 
Double_t 
  AliCascadeVertexer::fgDV0min=0.01;   //min V0 impact parameter
Double_t 
  AliCascadeVertexer::fgMassWin=0.008; //"window" around the Lambda mass
Double_t 
  AliCascadeVertexer::fgDBachMin=0.01; //min bachelor impact parameter
Double_t 
  AliCascadeVertexer::fgDCAmax=2.0;    //max DCA between the V0 and the track 
Double_t 
  AliCascadeVertexer::fgCPAmin=0.98; //min cosine of the cascade pointing angle
Double_t 
  AliCascadeVertexer::fgRmin=0.2;      //min radius of the fiducial volume
Double_t 
  AliCascadeVertexer::fgRmax=100.;     //max radius of the fiducial volume
  

Int_t AliCascadeVertexer::V0sTracks2CascadeVertices(AliESDEvent *event) {
  //--------------------------------------------------------------------
  // This function reconstructs cascade vertices
  //      Adapted to the ESD by I.Belikov (Jouri.Belikov@cern.ch)
  //--------------------------------------------------------------------
   const AliESDVertex *vtxT3D=event->GetPrimaryVertex();

   Double_t xPrimaryVertex=vtxT3D->GetX();
   Double_t yPrimaryVertex=vtxT3D->GetY();
   Double_t zPrimaryVertex=vtxT3D->GetZ();

   Double_t b=event->GetMagneticField();
   Int_t nV0=(Int_t)event->GetNumberOfV0s();

   //stores relevant V0s in an array
   TObjArray vtcs(nV0);
   Int_t i;
   for (i=0; i<nV0; i++) {
       AliESDv0 *v=event->GetV0(i);
       if (v->GetOnFlyStatus()) continue;
       if (v->GetD(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex)<fDV0min) continue;
       vtcs.AddLast(v);
   }
   nV0=vtcs.GetEntriesFast();

   // stores relevant tracks in another array
   Int_t nentr=(Int_t)event->GetNumberOfTracks();
   TArrayI trk(nentr); Int_t ntr=0;
   for (i=0; i<nentr; i++) {
       AliESDtrack *esdtr=event->GetTrack(i);
       ULong_t status=esdtr->GetStatus();

       if ((status&AliESDtrack::kITSrefit)==0)
	  if ((status&AliESDtrack::kTPCrefit)==0) continue;

       if (TMath::Abs(esdtr->GetD(xPrimaryVertex,yPrimaryVertex,b))<fDBachMin) continue;

       trk[ntr++]=i;
   }   

   Double_t massLambda=1.11568;
   Int_t ncasc=0;

   // Looking for the cascades...

   for (i=0; i<nV0; i++) { //loop on V0s

      AliESDv0 *v=(AliESDv0*)vtcs.UncheckedAt(i);
      AliESDv0 v0(*v);
      v0.ChangeMassHypothesis(kLambda0); // the v0 must be Lambda 
      if (TMath::Abs(v0.GetEffMass()-massLambda)>fMassWin) continue; 

      for (Int_t j=0; j<ntr; j++) {//loop on tracks
	 Int_t bidx=trk[j];
 	 //Bo:   if (bidx==v->GetNindex()) continue; //bachelor and v0's negative tracks must be different
         if (bidx==v0.GetIndex(0)) continue; //Bo:  consistency 0 for neg
	 AliESDtrack *btrk=event->GetTrack(bidx);
         if (btrk->GetSign()>0) continue;  // bachelor's charge 
          
    	 AliESDv0 *pv0=&v0;
         AliExternalTrackParam bt(*btrk), *pbt=&bt;

         Double_t dca=PropagateToDCA(pv0,pbt,b);
         if (dca > fDCAmax) continue;

         AliESDcascade cascade(*pv0,*pbt,bidx);//constucts a cascade candidate
	 //PH        if (cascade.GetChi2Xi() > fChi2max) continue;

	 Double_t x,y,z; cascade.GetXYZcascade(x,y,z); // Bo: bug correction
         Double_t r2=x*x + y*y; 
         if (r2 > fRmax*fRmax) continue;   // condition on fiducial zone
         if (r2 < fRmin*fRmin) continue;

	 Double_t pxV0,pyV0,pzV0;
	 pv0->GetPxPyPz(pxV0,pyV0,pzV0);
	 if (x*pxV0+y*pyV0+z*pzV0 < 0) continue; //causality

         Double_t x1,y1,z1; pv0->GetXYZ(x1,y1,z1);
         if (r2 > (x1*x1+y1*y1)) continue;

  	 if (cascade.GetCascadeCosineOfPointingAngle(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex) <fCPAmin) continue; //condition on the cascade pointing angle 
	 
         cascade.SetDcaXiDaughters(dca);
	 event->AddCascade(&cascade);
         ncasc++;
      } // end loop tracks
   } // end loop V0s

   // Looking for the anti-cascades...

   for (i=0; i<nV0; i++) { //loop on V0s
      AliESDv0 *v=(AliESDv0*)vtcs.UncheckedAt(i);
      AliESDv0 v0(*v);
      v0.ChangeMassHypothesis(kLambda0Bar); //the v0 must be anti-Lambda 
      if (TMath::Abs(v0.GetEffMass()-massLambda)>fMassWin) continue; 

      for (Int_t j=0; j<ntr; j++) {//loop on tracks
	 Int_t bidx=trk[j];
 	 //Bo:   if (bidx==v->GetPindex()) continue; //bachelor and v0's positive tracks must be different
         if (bidx==v0.GetIndex(1)) continue; //Bo:  consistency 1 for pos
	 AliESDtrack *btrk=event->GetTrack(bidx);
         if (btrk->GetSign()<0) continue;  // bachelor's charge 
          
	 AliESDv0 *pv0=&v0;
         AliExternalTrackParam bt(*btrk), *pbt=&bt;

         Double_t dca=PropagateToDCA(pv0,pbt,b);
         if (dca > fDCAmax) continue;

         AliESDcascade cascade(*pv0,*pbt,bidx); //constucts a cascade candidate
	 //PH         if (cascade.GetChi2Xi() > fChi2max) continue;

	 Double_t x,y,z; cascade.GetXYZcascade(x,y,z); // Bo: bug correction
         Double_t r2=x*x + y*y; 
         if (r2 > fRmax*fRmax) continue;   // condition on fiducial zone
         if (r2 < fRmin*fRmin) continue;

	 Double_t pxV0,pyV0,pzV0;
	 pv0->GetPxPyPz(pxV0,pyV0,pzV0);
	 if (x*pxV0+y*pyV0+z*pzV0 < 0) continue; //causality

         Double_t x1,y1,z1; pv0->GetXYZ(x1,y1,z1);
         if (r2 > (x1*x1+y1*y1)) continue;

	 if (cascade.GetCascadeCosineOfPointingAngle(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex) < fCPAmin) continue; //condition on the cascade pointing angle 

         cascade.SetDcaXiDaughters(dca);
	 event->AddCascade(&cascade);
         ncasc++;

      } // end loop tracks
   } // end loop V0s

Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc);

   return 0;
}


Double_t AliCascadeVertexer::Det(Double_t a00, Double_t a01, Double_t a10, Double_t a11) const {
  //--------------------------------------------------------------------
  // This function calculates locally a 2x2 determinant
  //--------------------------------------------------------------------
  return a00*a11 - a01*a10;
}

Double_t AliCascadeVertexer::Det(Double_t a00,Double_t a01,Double_t a02,
				 Double_t a10,Double_t a11,Double_t a12,
				 Double_t a20,Double_t a21,Double_t a22) const {
  //--------------------------------------------------------------------
  // This function calculates locally a 3x3 determinant
  //--------------------------------------------------------------------
  return  a00*Det(a11,a12,a21,a22)-a01*Det(a10,a12,a20,a22)+a02*Det(a10,a11,a20,a21);
}




Double_t AliCascadeVertexer::PropagateToDCA(AliESDv0 *v, AliExternalTrackParam *t, Double_t b) {
  //--------------------------------------------------------------------
  // This function returns the DCA between the V0 and the track
  //--------------------------------------------------------------------
  Double_t alpha=t->GetAlpha(), cs1=TMath::Cos(alpha), sn1=TMath::Sin(alpha);
  Double_t r[3]; t->GetXYZ(r);
  Double_t x1=r[0], y1=r[1], z1=r[2];
  Double_t p[3]; t->GetPxPyPz(p);
  Double_t px1=p[0], py1=p[1], pz1=p[2];
  
  Double_t x2,y2,z2;     // position and momentum of V0
  Double_t px2,py2,pz2;
  
  v->GetXYZ(x2,y2,z2);
  v->GetPxPyPz(px2,py2,pz2);
 
// calculation dca
   
  Double_t dd= Det(x2-x1,y2-y1,z2-z1,px1,py1,pz1,px2,py2,pz2);
  Double_t ax= Det(py1,pz1,py2,pz2);
  Double_t ay=-Det(px1,pz1,px2,pz2);
  Double_t az= Det(px1,py1,px2,py2);

  Double_t dca=TMath::Abs(dd)/TMath::Sqrt(ax*ax + ay*ay + az*az);

//points of the DCA
  Double_t t1 = Det(x2-x1,y2-y1,z2-z1,px2,py2,pz2,ax,ay,az)/
                Det(px1,py1,pz1,px2,py2,pz2,ax,ay,az);
  
  x1 += px1*t1; y1 += py1*t1; //z1 += pz1*t1;
  

  //propagate track to the points of DCA

  x1=x1*cs1 + y1*sn1;
  if (!t->PropagateTo(x1,b)) {
    Error("PropagateToDCA","Propagation failed !");
    return 1.e+33;
  }  

  return dca;
}











 AliCascadeVertexer.cxx:1
 AliCascadeVertexer.cxx:2
 AliCascadeVertexer.cxx:3
 AliCascadeVertexer.cxx:4
 AliCascadeVertexer.cxx:5
 AliCascadeVertexer.cxx:6
 AliCascadeVertexer.cxx:7
 AliCascadeVertexer.cxx:8
 AliCascadeVertexer.cxx:9
 AliCascadeVertexer.cxx:10
 AliCascadeVertexer.cxx:11
 AliCascadeVertexer.cxx:12
 AliCascadeVertexer.cxx:13
 AliCascadeVertexer.cxx:14
 AliCascadeVertexer.cxx:15
 AliCascadeVertexer.cxx:16
 AliCascadeVertexer.cxx:17
 AliCascadeVertexer.cxx:18
 AliCascadeVertexer.cxx:19
 AliCascadeVertexer.cxx:20
 AliCascadeVertexer.cxx:21
 AliCascadeVertexer.cxx:22
 AliCascadeVertexer.cxx:23
 AliCascadeVertexer.cxx:24
 AliCascadeVertexer.cxx:25
 AliCascadeVertexer.cxx:26
 AliCascadeVertexer.cxx:27
 AliCascadeVertexer.cxx:28
 AliCascadeVertexer.cxx:29
 AliCascadeVertexer.cxx:30
 AliCascadeVertexer.cxx:31
 AliCascadeVertexer.cxx:32
 AliCascadeVertexer.cxx:33
 AliCascadeVertexer.cxx:34
 AliCascadeVertexer.cxx:35
 AliCascadeVertexer.cxx:36
 AliCascadeVertexer.cxx:37
 AliCascadeVertexer.cxx:38
 AliCascadeVertexer.cxx:39
 AliCascadeVertexer.cxx:40
 AliCascadeVertexer.cxx:41
 AliCascadeVertexer.cxx:42
 AliCascadeVertexer.cxx:43
 AliCascadeVertexer.cxx:44
 AliCascadeVertexer.cxx:45
 AliCascadeVertexer.cxx:46
 AliCascadeVertexer.cxx:47
 AliCascadeVertexer.cxx:48
 AliCascadeVertexer.cxx:49
 AliCascadeVertexer.cxx:50
 AliCascadeVertexer.cxx:51
 AliCascadeVertexer.cxx:52
 AliCascadeVertexer.cxx:53
 AliCascadeVertexer.cxx:54
 AliCascadeVertexer.cxx:55
 AliCascadeVertexer.cxx:56
 AliCascadeVertexer.cxx:57
 AliCascadeVertexer.cxx:58
 AliCascadeVertexer.cxx:59
 AliCascadeVertexer.cxx:60
 AliCascadeVertexer.cxx:61
 AliCascadeVertexer.cxx:62
 AliCascadeVertexer.cxx:63
 AliCascadeVertexer.cxx:64
 AliCascadeVertexer.cxx:65
 AliCascadeVertexer.cxx:66
 AliCascadeVertexer.cxx:67
 AliCascadeVertexer.cxx:68
 AliCascadeVertexer.cxx:69
 AliCascadeVertexer.cxx:70
 AliCascadeVertexer.cxx:71
 AliCascadeVertexer.cxx:72
 AliCascadeVertexer.cxx:73
 AliCascadeVertexer.cxx:74
 AliCascadeVertexer.cxx:75
 AliCascadeVertexer.cxx:76
 AliCascadeVertexer.cxx:77
 AliCascadeVertexer.cxx:78
 AliCascadeVertexer.cxx:79
 AliCascadeVertexer.cxx:80
 AliCascadeVertexer.cxx:81
 AliCascadeVertexer.cxx:82
 AliCascadeVertexer.cxx:83
 AliCascadeVertexer.cxx:84
 AliCascadeVertexer.cxx:85
 AliCascadeVertexer.cxx:86
 AliCascadeVertexer.cxx:87
 AliCascadeVertexer.cxx:88
 AliCascadeVertexer.cxx:89
 AliCascadeVertexer.cxx:90
 AliCascadeVertexer.cxx:91
 AliCascadeVertexer.cxx:92
 AliCascadeVertexer.cxx:93
 AliCascadeVertexer.cxx:94
 AliCascadeVertexer.cxx:95
 AliCascadeVertexer.cxx:96
 AliCascadeVertexer.cxx:97
 AliCascadeVertexer.cxx:98
 AliCascadeVertexer.cxx:99
 AliCascadeVertexer.cxx:100
 AliCascadeVertexer.cxx:101
 AliCascadeVertexer.cxx:102
 AliCascadeVertexer.cxx:103
 AliCascadeVertexer.cxx:104
 AliCascadeVertexer.cxx:105
 AliCascadeVertexer.cxx:106
 AliCascadeVertexer.cxx:107
 AliCascadeVertexer.cxx:108
 AliCascadeVertexer.cxx:109
 AliCascadeVertexer.cxx:110
 AliCascadeVertexer.cxx:111
 AliCascadeVertexer.cxx:112
 AliCascadeVertexer.cxx:113
 AliCascadeVertexer.cxx:114
 AliCascadeVertexer.cxx:115
 AliCascadeVertexer.cxx:116
 AliCascadeVertexer.cxx:117
 AliCascadeVertexer.cxx:118
 AliCascadeVertexer.cxx:119
 AliCascadeVertexer.cxx:120
 AliCascadeVertexer.cxx:121
 AliCascadeVertexer.cxx:122
 AliCascadeVertexer.cxx:123
 AliCascadeVertexer.cxx:124
 AliCascadeVertexer.cxx:125
 AliCascadeVertexer.cxx:126
 AliCascadeVertexer.cxx:127
 AliCascadeVertexer.cxx:128
 AliCascadeVertexer.cxx:129
 AliCascadeVertexer.cxx:130
 AliCascadeVertexer.cxx:131
 AliCascadeVertexer.cxx:132
 AliCascadeVertexer.cxx:133
 AliCascadeVertexer.cxx:134
 AliCascadeVertexer.cxx:135
 AliCascadeVertexer.cxx:136
 AliCascadeVertexer.cxx:137
 AliCascadeVertexer.cxx:138
 AliCascadeVertexer.cxx:139
 AliCascadeVertexer.cxx:140
 AliCascadeVertexer.cxx:141
 AliCascadeVertexer.cxx:142
 AliCascadeVertexer.cxx:143
 AliCascadeVertexer.cxx:144
 AliCascadeVertexer.cxx:145
 AliCascadeVertexer.cxx:146
 AliCascadeVertexer.cxx:147
 AliCascadeVertexer.cxx:148
 AliCascadeVertexer.cxx:149
 AliCascadeVertexer.cxx:150
 AliCascadeVertexer.cxx:151
 AliCascadeVertexer.cxx:152
 AliCascadeVertexer.cxx:153
 AliCascadeVertexer.cxx:154
 AliCascadeVertexer.cxx:155
 AliCascadeVertexer.cxx:156
 AliCascadeVertexer.cxx:157
 AliCascadeVertexer.cxx:158
 AliCascadeVertexer.cxx:159
 AliCascadeVertexer.cxx:160
 AliCascadeVertexer.cxx:161
 AliCascadeVertexer.cxx:162
 AliCascadeVertexer.cxx:163
 AliCascadeVertexer.cxx:164
 AliCascadeVertexer.cxx:165
 AliCascadeVertexer.cxx:166
 AliCascadeVertexer.cxx:167
 AliCascadeVertexer.cxx:168
 AliCascadeVertexer.cxx:169
 AliCascadeVertexer.cxx:170
 AliCascadeVertexer.cxx:171
 AliCascadeVertexer.cxx:172
 AliCascadeVertexer.cxx:173
 AliCascadeVertexer.cxx:174
 AliCascadeVertexer.cxx:175
 AliCascadeVertexer.cxx:176
 AliCascadeVertexer.cxx:177
 AliCascadeVertexer.cxx:178
 AliCascadeVertexer.cxx:179
 AliCascadeVertexer.cxx:180
 AliCascadeVertexer.cxx:181
 AliCascadeVertexer.cxx:182
 AliCascadeVertexer.cxx:183
 AliCascadeVertexer.cxx:184
 AliCascadeVertexer.cxx:185
 AliCascadeVertexer.cxx:186
 AliCascadeVertexer.cxx:187
 AliCascadeVertexer.cxx:188
 AliCascadeVertexer.cxx:189
 AliCascadeVertexer.cxx:190
 AliCascadeVertexer.cxx:191
 AliCascadeVertexer.cxx:192
 AliCascadeVertexer.cxx:193
 AliCascadeVertexer.cxx:194
 AliCascadeVertexer.cxx:195
 AliCascadeVertexer.cxx:196
 AliCascadeVertexer.cxx:197
 AliCascadeVertexer.cxx:198
 AliCascadeVertexer.cxx:199
 AliCascadeVertexer.cxx:200
 AliCascadeVertexer.cxx:201
 AliCascadeVertexer.cxx:202
 AliCascadeVertexer.cxx:203
 AliCascadeVertexer.cxx:204
 AliCascadeVertexer.cxx:205
 AliCascadeVertexer.cxx:206
 AliCascadeVertexer.cxx:207
 AliCascadeVertexer.cxx:208
 AliCascadeVertexer.cxx:209
 AliCascadeVertexer.cxx:210
 AliCascadeVertexer.cxx:211
 AliCascadeVertexer.cxx:212
 AliCascadeVertexer.cxx:213
 AliCascadeVertexer.cxx:214
 AliCascadeVertexer.cxx:215
 AliCascadeVertexer.cxx:216
 AliCascadeVertexer.cxx:217
 AliCascadeVertexer.cxx:218
 AliCascadeVertexer.cxx:219
 AliCascadeVertexer.cxx:220
 AliCascadeVertexer.cxx:221
 AliCascadeVertexer.cxx:222
 AliCascadeVertexer.cxx:223
 AliCascadeVertexer.cxx:224
 AliCascadeVertexer.cxx:225
 AliCascadeVertexer.cxx:226
 AliCascadeVertexer.cxx:227
 AliCascadeVertexer.cxx:228
 AliCascadeVertexer.cxx:229
 AliCascadeVertexer.cxx:230
 AliCascadeVertexer.cxx:231
 AliCascadeVertexer.cxx:232
 AliCascadeVertexer.cxx:233
 AliCascadeVertexer.cxx:234
 AliCascadeVertexer.cxx:235
 AliCascadeVertexer.cxx:236
 AliCascadeVertexer.cxx:237
 AliCascadeVertexer.cxx:238
 AliCascadeVertexer.cxx:239
 AliCascadeVertexer.cxx:240
 AliCascadeVertexer.cxx:241
 AliCascadeVertexer.cxx:242
 AliCascadeVertexer.cxx:243
 AliCascadeVertexer.cxx:244
 AliCascadeVertexer.cxx:245
 AliCascadeVertexer.cxx:246
 AliCascadeVertexer.cxx:247
 AliCascadeVertexer.cxx:248
 AliCascadeVertexer.cxx:249
 AliCascadeVertexer.cxx:250
 AliCascadeVertexer.cxx:251
 AliCascadeVertexer.cxx:252
 AliCascadeVertexer.cxx:253
 AliCascadeVertexer.cxx:254
 AliCascadeVertexer.cxx:255
 AliCascadeVertexer.cxx:256
 AliCascadeVertexer.cxx:257
 AliCascadeVertexer.cxx:258
 AliCascadeVertexer.cxx:259
 AliCascadeVertexer.cxx:260
 AliCascadeVertexer.cxx:261
 AliCascadeVertexer.cxx:262