FairRoot/PandaRoot
PndFsmTrack.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: PndFsmTrack.cc,v 1.16 2006/10/06 15:19:11 aida Exp $
4 //
5 // Description:
6 // Class PndFsmTrack
7 //
8 // Candidate "Tracks" or "Particles" for the Fast Simulation
9 //
10 // This software was developed for the PANDA collaboration. If you
11 // use all or part of it, please give an appropriate acknowledgement.
12 //
13 // Author List:
14 // Klaus Goetzen Original Author
15 //
16 // Copyright Information:
17 // Copyright (C) 2006 GSI
18 //
19 //------------------------------------------------------------------------
20 
21 //-----------------------
22 // This Class's Header --
23 //-----------------------
24 #include "PndFsmTrack.h"
25 #include "PndFsmResponse.h"
26 //#include "RhoBase/TRho.h"
27 //#include "FastSimApp/FsmHitMap.hh"
28 //#include "FsmDetTypes.hh"
29 #include "TDatabasePDG.h"
30 #include "TParticlePDG.h"
31 #include "TMatrixD.h"
32 #include "FairRunAna.h"
33 #include "FairRunSim.h"
34 #include "FairField.h"
35 //-------------
36 // C Headers --
37 //-------------
38 #include <iostream>
39 #include <cmath>
40 using std::cout;
41 using std::endl;
42 using std::ostream;
43 
44 
45 //---------------
46 // C++ Headers --
47 //---------------
48 
49 //-------------------------------
50 // Collaborating Class Headers --
51 //-------------------------------
52 //#include "PDT/PdtEntry.hh"
53 
54 //-----------------------------------------------------------------------
55 // Local Macros, Typedefs, Structures, Unions and Forward Declarations --
56 //-----------------------------------------------------------------------
57 
58 //----------------
59 // Constructors --
60 //----------------
61 
62 
64  setP4(TLorentzVector(0.,0.,0.,0.));
65  setStartVtx(TVector3(0.,0.,0.));
66  setStopVtx(TVector3(0.,0.,0.));
67  setPdt(0);
68  setMass2(0.);
69  setMvddEdX(0.);
70  setTpcdEdX(0.);
71  setSttdEdX(0.);
72  setCharge(0);
73  setGTrackId(0);
74  setDetResponse(0);
75  for (int i=0;i<15;i++)
76  fCov5[i][0]=0;
77  for (int i=0;i<5;i++)
78  fPar5[i]=0;
79  for (int i=0;i<28;i++)
80  fCov7[i][0]=0;
81 }
82 
83 PndFsmTrack::PndFsmTrack(TLorentzVector const inP4, TVector3 start, TVector3 stop, double inCharge, int inPdt, signed long trackId)
84 : fCov5(5,5), fCov7(7,7) {
85  setP4(inP4);
86  setStartVtx(start);
87  setStopVtx(stop);
88  setCharge(inCharge);
89  setPdt(inPdt);
90  setGTrackId(trackId);
91  setDetResponse(0);
92  setMass2( 0.0);
93  setMvddEdX(0.0);
94  setTpcdEdX(0.0);
95  setSttdEdX(0.0);
96 }
97 
98 void PndFsmTrack::HelixRep(TVector3 reference) {
99  fReference=reference;
100  if (fabs(charge())>1e-6) {
101  // calculate helix track representation (as in TFitParams.h)
102  reference=_startVtx-fReference;
103  double tandip=p4().Pz()/p4().Perp();
104  double pnt[3], Bf[3];
105  pnt[0]=fReference.X(); pnt[1]=fReference.Y(); pnt[2]=fReference.Z();
106  FairField* theField=0;
107  if(FairRun::Instance()->IsAna()){
108  theField=FairRunAna::Instance()->GetField();
109 // FairRunAna::Instance()->GetField()->GetFieldValue(pnt, Bf); //[kGs]
110  }else{
111  theField=FairRunSim::Instance()->GetField();
112 // FairRunSim::Instance()->GetField()->GetFieldValue(pnt, Bf); //[kGs]
113  }
114  if(theField==0) Fatal("PndFsmTrack::HelixRep()","Magnetic Field pointer missing. Set your field!");
115  theField->GetFieldValue(pnt, Bf); //[kGs]
116 
117 
118  double a=-2.99792458e-3*Bf[2]*charge();
119  double omega=a/p4().Perp();
120  // construct helix center (1/omega=R)
121  TVector3 center = reference + (1/omega)*(TVector3( -p4().Y(), p4().X(), 0).Unit());
122  // construct propagation length
123  double delta = (center-reference).Phi() - center.Phi();
124  double cosrs=cos(delta);
125  double sinrs=sin(delta);
126 
127  TVector3 p(p4().Vect());
128 
129  p.SetZ( p4().Pz() );
130  p.SetPhi( p4().Phi() - delta );
131 
132  reference.SetX( reference.X() - p.X()/a*sinrs + p.Y()/a*(1-cosrs) );
133  reference.SetY( reference.Y() - p.Y()/a*sinrs - p.X()/a*(1-cosrs) );
134  reference.SetZ( reference.Z() - tandip*delta/omega );
135 
136  fPar5[0]=reference.Cross(p).Z()<0 ? reference.Perp() : -reference.Perp();
137  fPar5[1]=p.Phi();
138  fPar5[2]=omega;
139  fPar5[3]=reference.Z();
140  fPar5[4]=tandip;
141  } else {
142  fPar5[0]=p4().X();
143  fPar5[1]=p4().Y();
144  fPar5[2]=p4().Z();
145  fPar5[3]=p4().T();
146  }
147 }
148 
149 void PndFsmTrack::Propagate(TVector3 origin, double deltaError) {
150  origin-=fReference;
151  // calculate p4 and start vertex at point
152  // on helix track closest to origin
153  double pnt[3], Bf[3];
154  pnt[0]=fReference.X(); pnt[1]=fReference.Y(); pnt[2]=fReference.Z();
155  FairField* theField=0;
156  if(FairRun::Instance()->IsAna()){
157  theField=FairRunAna::Instance()->GetField();
158 // FairRunAna::Instance()->GetField()->GetFieldValue(pnt, Bf); //[kGs]
159  }else{
160  theField=FairRunSim::Instance()->GetField();
161 // FairRunSim::Instance()->GetField()->GetFieldValue(pnt, Bf); //[kGs]
162  }
163  if(theField==0) Fatal("PndFsmTrack::HelixRep()","Magnetic Field pointer missing. Set your field!");
164  theField->GetFieldValue(pnt, Bf); //[kGs]
165  double a=2.99792458e-3*Bf[2];
166  double R=1/GetHelixOmega();
167  double pt=-a*R*charge();
168  double s0=sin(GetHelixPhi0());
169  double c0=cos(GetHelixPhi0());
170  // add some helix revolutions to get somewhere near origin
171  if (fabs(GetHelixTanDip()*R)>1e-9) {
172  double dz=2*M_PI*fabs(R*GetHelixTanDip());
173  fPar5[3]-=dz*floor((fPar5[3]-origin.Z())/dz+0.5);
174  }
175  // use gradient method to find POCA to origin
176  double s1;
177  double c1;
178  double ds=0;
179  double delta=0;
180  double alpha=1;
181  double distance2=1e150;
182  do {
183  delta+=ds;
184  s1=sin(GetHelixPhi0()+delta);
185  c1=cos(GetHelixPhi0()+delta);
186  // vertex setup
187  _startVtx.SetXYZ(-s0*(GetHelixD0()+R)+s1*R,
188  c0*(GetHelixD0()+R)-c1*R,
189  GetHelixZ0()+GetHelixTanDip()*R*delta );
190  TVector3 distance(_startVtx-origin);
191  if (distance2<distance.Mag2()) {
192  alpha*=0.5;
193  delta-=ds;
194  }
195  distance2=distance.Mag2();
196  // construct helix tangent at delta
197  TVector3 u(c1,s1,GetHelixTanDip());
198  ds=-u.Dot(distance)/(1+GetHelixTanDip()*GetHelixTanDip())*alpha*GetHelixOmega();
199  } while (fabs(ds)>1e-9);
200 
201  // momentum setup
202  _p4.SetXYZM( pt*c1, pt*s1, pt*GetHelixTanDip(),
203  TDatabasePDG::Instance()->GetParticle("pi-")->Mass());
204 
205  // calculate jacobian wrt d0..tandip
206  TMatrixD J_alpha(7,5);
207  J_alpha(0,0)=-s0;
208  J_alpha(0,1)=-_startVtx.Y();
209  J_alpha(0,2)=-R*R*(s1-s0);
210 
211  J_alpha(1,0)=+c0;
212  J_alpha(1,1)=+_startVtx.X();
213  J_alpha(1,2)=+R*R*(c1-c0);
214 
215  J_alpha(2,2)=-delta*GetHelixTanDip()*R*R;
216  J_alpha(2,3)=+1;
217  J_alpha(2,4)=+delta*R;
218 
219  J_alpha(3,1)=-_p4.Y();
220  J_alpha(3,2)=-_p4.X()*R;
221 
222  J_alpha(4,1)=+_p4.X();
223  J_alpha(4,2)=-_p4.Y()*R;
224 
225  J_alpha(5,2)=-_p4.Z()*R;
226  J_alpha(5,4)=+pt;
227 
228  J_alpha(6,2)=-_p4.Vect().Mag2()*R/_p4.T();
229  J_alpha(6,4)=+pt*pt*GetHelixTanDip()/_p4.T();
230 
231  if (deltaError>=0.001) {
232  // calculate jacobian wrt delta to allow fitter
233  // to move 1deg along the linearized trajectory
234  deltaError*=3.1416/180;
235  double covDelta=deltaError*deltaError;
236  TMatrixD J_delta(7,1);
237  J_delta(0,0)=+R*c1;
238  J_delta(1,0)=+R*s1;
239  J_delta(2,0)=+GetHelixTanDip()*R;
240  J_delta(3,0)=-pt*s1;
241  J_delta(4,0)=+pt*c1;
242  TMatrixD tmp2(J_delta, TMatrixD::kMultTranspose, J_delta);
243  tmp2*=covDelta;
244  fCov7+=tmp2;
245  }
246 
247  // calculate fCov7 = J_alpha * fCov5 * J_alpha.T +
248  // covDelta * J_delta * J_delta.T (covDelta is scalar)
249  TMatrixD tmp1(J_alpha, TMatrixD::kMult, fCov5);
250  fCov7.MultT(tmp1, J_alpha);
251 
253 }
254 
256 {
257  for (int i=0;i<7;i++)
258  for (int j=0;j<7;j++)
259  {
260  fCov7[i][j]=p7cov[i][j];
261  }
262 }
263 
265 {
266  for (int i=0;i<4;i++)
267  for (int j=0;j<4;j++)
268  {
269  fCov7[i+3][j+3]=p4cov[i][j];
270  }
271 }
272 
274 {
275  for (int i=0;i<3;i++)
276  for (int j=0;j<3;j++)
277  {
278  fCov7[i][j]=vcov[i][j];
279  }
280 }
281 
282 //--------------
283 // Destructor --
284 //--------------
285 
287 {
288  if (_detResponse)
289  delete _detResponse;
290 }
291 
292 //--------------
293 // Operations --
294 //--------------
295 
296 void
297 PndFsmTrack::setP4(TLorentzVector l)
298 {
299  _p4=l;
300 }
301 
302 void
304 {
305  _startVtx=v;
306 }
307 
308 void
310 {
311  _stopVtx=v;
312 }
313 
314 void
316 {
317  _charge=c;
318 }
319 
320 void
322 {
323  _Mass2=c;
324 }
325 
326 
327 void
329 {
330 _MvddEdX=c;
331 }
332 
333 void
335 {
336 _TpcdEdX=c;
337 }
338 
339 void
341 {
342 _SttdEdX=c;
343 }
344 
345 void
347 {
348  _pdt=n;
349 }
350 
351 void
353 {
354  _gTrackId=id;
355 }
356 
357 
358 void
360 {
361  _detResponse=resp;
362 }
363 
364 
365 
367 {
368  o<<"PndFsmTrack printout"<<endl;
369  o<<" P4 : < "<< _p4.Px() <<" / " <<_p4.Py() <<" / " <<_p4.Pz() <<" / "<<_p4.E()<<" >"<<endl;
370  o<<" Vtx1 : < "<< _startVtx.X()<<" / " <<_startVtx.Y()<<" / " <<_startVtx.Z() <<" > " << endl;
371  o<<" Vtx2 : < "<< _stopVtx.X() <<" / " <<_stopVtx.Y() <<" / " <<_stopVtx.Z( )<<" > " << endl;
372  o<<" charge = "<< _charge <<" / lundId = "<<_pdt <<" / gTrackId = "<<_gTrackId <<endl;
373  o<<" D0: "<<GetHelixD0();
374  o<<" Phi0: "<<GetHelixPhi0()*180/3.1416<<"deg";
375  if (GetHelixOmega()) o <<" 1/Omega: "<<1/GetHelixOmega();
376  else o <<" 1/Omega: " << "inf";
377  o<<" Z0: "<<GetHelixZ0();
378  o<<" TanDip: "<<GetHelixTanDip()<<endl;
379  // if (_detResponse) _detResponse->print(o);
380 }
381 
Double_t GetHelixZ0() const
Definition: PndFsmTrack.h:140
double _charge
Definition: PndFsmTrack.h:115
TMatrixD fCov7
Definition: PndFsmTrack.h:127
Double_t p
Definition: anasim.C:58
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
Double_t GetHelixPhi0() const
Definition: PndFsmTrack.h:138
void setMass2(double c)
Int_t i
Definition: run_full.C:25
void setCharge(double c)
void setStopVtx(TVector3 v)
void setPdt(int pdt)
TVector3 _startVtx
Definition: PndFsmTrack.h:113
Double_t GetHelixOmega() const
Definition: PndFsmTrack.h:139
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
void setDetResponse(PndFsmResponse *resp)
double _Mass2
Definition: PndFsmTrack.h:118
int n
TClonesArray * pnt
Double_t GetHelixD0() const
Definition: PndFsmTrack.h:137
void setMvddEdX(double c)
double charge()
Definition: PndFsmTrack.h:75
double _SttdEdX
Definition: PndFsmTrack.h:121
TVector3 fReference
Definition: PndFsmTrack.h:125
double Y
Definition: anaLmdDigi.C:68
__m128 v
Definition: P4_F32vec4.h:4
void SetP7Cov(TMatrixD &p7cov)
PndFsmResponse * _detResponse
Definition: PndFsmTrack.h:122
TString pt(TString pts, TString exts="px py pz")
Definition: invexp.C:133
Int_t a
Definition: anaLmdDigi.C:126
TLorentzVector p4()
Definition: PndFsmTrack.h:72
void SetVCov(TMatrixD &vcov)
signed long _gTrackId
Definition: PndFsmTrack.h:117
basic_ostream< char, char_traits< char > > ostream
double _TpcdEdX
Definition: PndFsmTrack.h:120
double _MvddEdX
Definition: PndFsmTrack.h:119
TMatrixD fCov5
Definition: PndFsmTrack.h:126
void setP4(TLorentzVector l)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
c1
Definition: plot_dirc.C:35
TLorentzVector _p4
Definition: PndFsmTrack.h:112
double fPar5[5]
Definition: PndFsmTrack.h:124
void setSttdEdX(double c)
double X
Definition: anaLmdDigi.C:68
void setTpcdEdX(double c)
void SetP4Cov(TMatrixD &p4cov)
void setStartVtx(TVector3 v)
virtual ~PndFsmTrack()
void HelixRep(TVector3 reference)
Definition: PndFsmTrack.cxx:98
double alpha
Definition: f_Init.h:9
Double_t GetHelixTanDip() const
Definition: PndFsmTrack.h:141
void setGTrackId(signed long id)
Double_t R
Definition: checkhelixhit.C:61
void Propagate(TVector3 origin, double deltaError=2.5)
void print(std::ostream &o)
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
TVector3 _stopVtx
Definition: PndFsmTrack.h:114