]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliGenKinkInfo.cxx
Add parameterizations for "uniform" mag field a la former AliMagFC.
[u/mrichter/AliRoot.git] / PWG1 / AliGenKinkInfo.cxx
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 ///////////////////////////////////////////////////////////////////////////
18 /*
19
20 Origin: marian.ivanov@cern.ch
21 Container classes with MC infomation
22
23 */
24
25 #if !defined(__CINT__) || defined(__MAKECINT__)
26 #include <stdio.h>
27 #include <string.h>
28 //ROOT includes
29 #include "TROOT.h"
30 #include "Rtypes.h"
31 #include "TFile.h"
32 #include "TTree.h"
33 #include "TChain.h"
34 #include "TCut.h"
35 #include "TString.h"
36 #include "TStopwatch.h"
37 #include "TParticle.h"
38 #include "TSystem.h"
39 #include "TCanvas.h"
40 #include "TPolyLine3D.h"
41
42 //ALIROOT includes
43 #include "AliRun.h"
44 #include "AliStack.h"
45 #include "AliSimDigits.h"
46 #include "AliTPCParam.h"
47 #include "AliTPC.h"
48 #include "AliTPCLoader.h"
49 #include "AliDetector.h"
50 #include "AliTrackReference.h"
51 #include "AliTPCParamSR.h"
52 #include "AliTracker.h"
53 #include "AliMagF.h"
54 #include "AliHelix.h"
55 #include "AliTrackPointArray.h"
56
57 #endif
58 #include "AliGenKinkInfo.h"
59
60 //
61 // 
62
63 ClassImp(AliGenKinkInfo)
64
65
66
67 /////////////////////////////////////////////////////////////////////////////////
68
69 AliGenKinkInfo::AliGenKinkInfo():
70   fMCd(),       //info about daughter particle - second particle for V0
71   fMCm(),       //info about mother particle   - first particle for V0
72   fMCDist1(0),  //info about closest distance according closest MC - linear DCA
73   fMCDist2(0),  //info about closest distance parabolic DCA
74   fMCRr(0),     // rec position of the vertex 
75   fMCR(0)       //exact r position of the vertex
76 {
77   //
78   // default constructor
79   //
80   for (Int_t i=0;i<3;i++){
81     fMCPdr[i]=0;
82     fMCPd[i]=0;
83     fMCX[i]=0;
84     fMCPm[i]=0;
85     fMCAngle[i]=0;
86   }
87   for (Int_t i=0; i<2; i++) {
88     fPdg[i]= 0;
89     fLab[i]=0;
90   }
91 }
92
93 void AliGenKinkInfo::Update()
94 {
95   //
96   // Update information
97   //  some redundancy - faster acces to the values in analysis code
98   //
99   fMCPd[0] = fMCd.GetParticle().Px();
100   fMCPd[1] = fMCd.GetParticle().Py();
101   fMCPd[2] = fMCd.GetParticle().Pz();
102   fMCPd[3] = fMCd.GetParticle().P();
103   //
104   fMCX[0]  = fMCd.GetParticle().Vx();
105   fMCX[1]  = fMCd.GetParticle().Vy();
106   fMCX[2]  = fMCd.GetParticle().Vz();
107   fMCR       = TMath::Sqrt( fMCX[0]*fMCX[0]+fMCX[1]*fMCX[1]);
108   //
109   fPdg[0]    = fMCd.GetParticle().GetPdgCode();
110   fPdg[1]    = fMCm.GetParticle().GetPdgCode();
111   //
112   fLab[0]    = fMCd.GetParticle().GetUniqueID();
113   fLab[1]    = fMCm.GetParticle().GetUniqueID();
114   //
115   //
116   //
117   Double_t x1[3],p1[3];
118   TParticle& pdaughter = fMCd.GetParticle();
119   x1[0] = pdaughter.Vx();      
120   x1[1] = pdaughter.Vy();
121   x1[2] = pdaughter.Vz();
122   p1[0] = pdaughter.Px();
123   p1[1] = pdaughter.Py();
124   p1[2] = pdaughter.Pz();
125   Double_t sign = (pdaughter.GetPDG()->Charge()>0)? -1:1;
126   AliHelix dhelix1(x1,p1,sign);
127   //
128   //
129   Double_t x2[3],p2[3];            
130   //
131   TParticle& pmother = fMCm.GetParticle();
132   x2[0] = pmother.Vx();      
133   x2[1] = pmother.Vy();      
134   x2[2] = pmother.Vz();      
135   p2[0] = pmother.Px();
136   p2[1] = pmother.Py();
137   p2[2] = pmother.Pz();
138   //
139   const AliTrackReference & pdecay = fMCm.GetTRdecay();
140   x2[0] = pdecay.X();      
141   x2[1] = pdecay.Y();      
142   x2[2] = pdecay.Z();      
143   p2[0] = pdecay.Px();
144   p2[1] = pdecay.Py();
145   p2[2] = pdecay.Pz();
146   //
147   sign = (pmother.GetPDG()->Charge()>0) ? -1:1;
148   AliHelix mhelix(x2,p2,sign);
149   
150   //
151   //
152   //
153   //find intersection linear
154   //
155   Double_t distance1, distance2;
156   Double_t phase[2][2],radius[2];
157   Int_t  points = dhelix1.GetRPHIintersections(mhelix, phase, radius);
158   Double_t delta1=10000,delta2=10000;    
159   if (points>0){
160     dhelix1.LinearDCA(mhelix,phase[0][0],phase[0][1],radius[0],delta1);
161     dhelix1.LinearDCA(mhelix,phase[0][0],phase[0][1],radius[0],delta1);
162     dhelix1.LinearDCA(mhelix,phase[0][0],phase[0][1],radius[0],delta1);
163   }
164   if (points==2){    
165     dhelix1.LinearDCA(mhelix,phase[1][0],phase[1][1],radius[1],delta2);
166     dhelix1.LinearDCA(mhelix,phase[1][0],phase[1][1],radius[1],delta2);
167     dhelix1.LinearDCA(mhelix,phase[1][0],phase[1][1],radius[1],delta2);
168   }
169   distance1 = TMath::Min(delta1,delta2);
170   //
171   //find intersection parabolic
172   //
173   points = dhelix1.GetRPHIintersections(mhelix, phase, radius);
174   delta1=10000,delta2=10000;  
175   
176   if (points>0){
177     dhelix1.ParabolicDCA(mhelix,phase[0][0],phase[0][1],radius[0],delta1);
178   }
179   if (points==2){    
180     dhelix1.ParabolicDCA(mhelix,phase[1][0],phase[1][1],radius[1],delta2);
181   }
182   
183   distance2 = TMath::Min(delta1,delta2);
184   //
185   if (delta1<delta2){
186     //get V0 info
187     dhelix1.Evaluate(phase[0][0],fMCXr);
188     dhelix1.GetMomentum(phase[0][0],fMCPdr);
189     mhelix.GetMomentum(phase[0][1],fMCPm);
190     dhelix1.GetAngle(phase[0][0],mhelix,phase[0][1],fMCAngle);
191     fMCRr = TMath::Sqrt(radius[0]);
192   }
193   else{
194     dhelix1.Evaluate(phase[1][0],fMCXr);
195     dhelix1.GetMomentum(phase[1][0],fMCPdr);
196     mhelix.GetMomentum(phase[1][1],fMCPm);
197     dhelix1.GetAngle(phase[1][0],mhelix,phase[1][1],fMCAngle);
198     fMCRr = TMath::Sqrt(radius[1]);
199   }     
200   //            
201   //
202   fMCDist1 = TMath::Sqrt(distance1);
203   fMCDist2 = TMath::Sqrt(distance2);      
204     
205 }
206
207
208 Float_t AliGenKinkInfo::GetQt(){
209   //
210   //
211   Float_t momentumd = TMath::Sqrt(fMCPd[0]*fMCPd[0]+fMCPd[1]*fMCPd[1]+fMCPd[2]*fMCPd[2]);
212   return TMath::Sin(fMCAngle[2])*momentumd;
213 }
214
215