FairRoot/PandaRoot
Public Member Functions | Private Attributes | List of all members
PndFsmTrack Class Reference

#include <PndFsmTrack.h>

Public Member Functions

 PndFsmTrack ()
 
 PndFsmTrack (TLorentzVector p4, TVector3 start, TVector3 stop, double charge, int pdt, signed long trackId)
 
virtual ~PndFsmTrack ()
 
TLorentzVector p4 ()
 
TVector3 startVtx ()
 
TVector3 stopVtx ()
 
double charge ()
 
int pdt ()
 
signed long gTrackId ()
 
PndFsmResponsedetResponse ()
 
double Mass2 ()
 
double MvddEdX ()
 
double TpcdEdX ()
 
double SttdEdX ()
 
bool hitMapValid ()
 
bool hitMapResponse (unsigned int)
 
void setP4 (TLorentzVector l)
 
void setStartVtx (TVector3 v)
 
void setStopVtx (TVector3 v)
 
void setCharge (double c)
 
void setGTrackId (signed long id)
 
void setPdt (int pdt)
 
void setDetResponse (PndFsmResponse *resp)
 
void setMass2 (double c)
 
void setMvddEdX (double c)
 
void setTpcdEdX (double c)
 
void setSttdEdX (double c)
 
void print (std::ostream &o)
 
double * GetHelixParams ()
 
TMatrixDGetHelixCov ()
 
TMatrixDCov7 ()
 
void SetP7Cov (TMatrixD &p7cov)
 
void SetP4Cov (TMatrixD &p4cov)
 
void SetVCov (TMatrixD &vcov)
 
Double_t GetHelixD0 () const
 
Double_t GetHelixPhi0 () const
 
Double_t GetHelixOmega () const
 
Double_t GetHelixZ0 () const
 
Double_t GetHelixTanDip () const
 
void HelixRep (TVector3 reference)
 
void Propagate (TVector3 origin, double deltaError=2.5)
 

Private Attributes

TLorentzVector _p4
 
TVector3 _startVtx
 
TVector3 _stopVtx
 
double _charge
 
int _pdt
 
signed long _gTrackId
 
double _Mass2
 
double _MvddEdX
 
double _TpcdEdX
 
double _SttdEdX
 
PndFsmResponse_detResponse
 
double fPar5 [5]
 
TVector3 fReference
 
TMatrixD fCov5
 
TMatrixD fCov7
 

Detailed Description

Definition at line 49 of file PndFsmTrack.h.

Constructor & Destructor Documentation

PndFsmTrack::PndFsmTrack ( )

Definition at line 63 of file PndFsmTrack.cxx.

References fCov5, fCov7, fPar5, i, setCharge(), setDetResponse(), setGTrackId(), setMass2(), setMvddEdX(), setP4(), setPdt(), setStartVtx(), setStopVtx(), setSttdEdX(), and setTpcdEdX().

63  {
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 }
TMatrixD fCov7
Definition: PndFsmTrack.h:127
void setMass2(double c)
Int_t i
Definition: run_full.C:25
void setCharge(double c)
void setStopVtx(TVector3 v)
void setPdt(int pdt)
void setDetResponse(PndFsmResponse *resp)
void setMvddEdX(double c)
TMatrixD fCov5
Definition: PndFsmTrack.h:126
void setP4(TLorentzVector l)
double fPar5[5]
Definition: PndFsmTrack.h:124
void setSttdEdX(double c)
void setTpcdEdX(double c)
void setStartVtx(TVector3 v)
void setGTrackId(signed long id)
PndFsmTrack::PndFsmTrack ( TLorentzVector  p4,
TVector3  start,
TVector3  stop,
double  charge,
int  pdt,
signed long  trackId 
)

Definition at line 83 of file PndFsmTrack.cxx.

References setCharge(), setDetResponse(), setGTrackId(), setMass2(), setMvddEdX(), setP4(), setPdt(), setStartVtx(), setStopVtx(), setSttdEdX(), and setTpcdEdX().

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 }
TMatrixD fCov7
Definition: PndFsmTrack.h:127
void setMass2(double c)
void setCharge(double c)
void setStopVtx(TVector3 v)
void setPdt(int pdt)
void setDetResponse(PndFsmResponse *resp)
void setMvddEdX(double c)
TMatrixD fCov5
Definition: PndFsmTrack.h:126
void setP4(TLorentzVector l)
void setSttdEdX(double c)
void setTpcdEdX(double c)
void setStartVtx(TVector3 v)
void setGTrackId(signed long id)
PndFsmTrack::~PndFsmTrack ( )
virtual

Definition at line 286 of file PndFsmTrack.cxx.

References _detResponse.

287 {
288  if (_detResponse)
289  delete _detResponse;
290 }
PndFsmResponse * _detResponse
Definition: PndFsmTrack.h:122

Member Function Documentation

double PndFsmTrack::charge ( )
inline
TMatrixD& PndFsmTrack::Cov7 ( )
inline

Definition at line 132 of file PndFsmTrack.h.

References fCov7.

Referenced by PndFastSim::Exec().

132 {return fCov7;}
TMatrixD fCov7
Definition: PndFsmTrack.h:127
PndFsmResponse* PndFsmTrack::detResponse ( )
inline

Definition at line 78 of file PndFsmTrack.h.

References _detResponse.

Referenced by PndFastSim::cutAndSmear(), and PndFastSim::Exec().

78 {return _detResponse;}
PndFsmResponse * _detResponse
Definition: PndFsmTrack.h:122
TMatrixD& PndFsmTrack::GetHelixCov ( )
inline

Definition at line 131 of file PndFsmTrack.h.

References fCov5.

Referenced by PndFastSim::cutAndSmear().

131 {return fCov5; }
TMatrixD fCov5
Definition: PndFsmTrack.h:126
Double_t PndFsmTrack::GetHelixD0 ( ) const
inline

Definition at line 137 of file PndFsmTrack.h.

Referenced by print(), and Propagate().

137 {return fPar5[0];}
double fPar5[5]
Definition: PndFsmTrack.h:124
Double_t PndFsmTrack::GetHelixOmega ( ) const
inline

Definition at line 139 of file PndFsmTrack.h.

Referenced by PndFastSim::cutAndSmear(), print(), and Propagate().

139 {return fPar5[2];}
double fPar5[5]
Definition: PndFsmTrack.h:124
double* PndFsmTrack::GetHelixParams ( )
inline

Definition at line 130 of file PndFsmTrack.h.

References fPar5.

Referenced by PndFastSim::cutAndSmear().

130 { return fPar5; }
double fPar5[5]
Definition: PndFsmTrack.h:124
Double_t PndFsmTrack::GetHelixPhi0 ( ) const
inline

Definition at line 138 of file PndFsmTrack.h.

Referenced by print(), and Propagate().

138 {return fPar5[1];}
double fPar5[5]
Definition: PndFsmTrack.h:124
Double_t PndFsmTrack::GetHelixTanDip ( ) const
inline

Definition at line 141 of file PndFsmTrack.h.

Referenced by PndFastSim::cutAndSmear(), print(), and Propagate().

141 {return fPar5[4];}
double fPar5[5]
Definition: PndFsmTrack.h:124
Double_t PndFsmTrack::GetHelixZ0 ( ) const
inline

Definition at line 140 of file PndFsmTrack.h.

Referenced by print(), and Propagate().

140 {return fPar5[3];}
double fPar5[5]
Definition: PndFsmTrack.h:124
signed long PndFsmTrack::gTrackId ( )
inline

Definition at line 77 of file PndFsmTrack.h.

References _gTrackId.

77 {return _gTrackId;}
signed long _gTrackId
Definition: PndFsmTrack.h:117
void PndFsmTrack::HelixRep ( TVector3  reference)

Definition at line 98 of file PndFsmTrack.cxx.

References _startVtx, a, charge(), cos(), fabs(), fPar5, fReference, p, p4(), pnt, sin(), X, and Y.

Referenced by PndFastSim::cutAndSmear().

98  {
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 }
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
TVector3 _startVtx
Definition: PndFsmTrack.h:113
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
TClonesArray * pnt
double charge()
Definition: PndFsmTrack.h:75
TVector3 fReference
Definition: PndFsmTrack.h:125
double Y
Definition: anaLmdDigi.C:68
Double_t p
Definition: anasim.C:58
Int_t a
Definition: anaLmdDigi.C:126
TLorentzVector p4()
Definition: PndFsmTrack.h:72
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
double fPar5[5]
Definition: PndFsmTrack.h:124
double X
Definition: anaLmdDigi.C:68
bool PndFsmTrack::hitMapResponse ( unsigned  int)
inline
bool PndFsmTrack::hitMapValid ( )
inline
double PndFsmTrack::Mass2 ( )
inline

Definition at line 79 of file PndFsmTrack.h.

References _Mass2.

79 {return _Mass2;}
double _Mass2
Definition: PndFsmTrack.h:118
double PndFsmTrack::MvddEdX ( )
inline

Definition at line 80 of file PndFsmTrack.h.

References _MvddEdX.

80 {return _MvddEdX;}
double _MvddEdX
Definition: PndFsmTrack.h:119
TLorentzVector PndFsmTrack::p4 ( )
inline

Definition at line 72 of file PndFsmTrack.h.

References _p4.

Referenced by PndFastSim::cutAndSmear(), PndFsmEmcFS::dE(), PndFsmEmcFwCap::dE(), PndFsmEmcBwCap::dE(), PndFsmEmcBarrel::dE(), PndFsmDetTemplate::detected(), PndFsmSimpleTracker::detected(), PndFsmSimpleVtx::detected(), PndFsmStt::detected(), PndFsmSttPid::detected(), PndFsmTof::detected(), PndFsmEmcFS::detected(), PndFsmEmcFwCap::detected(), PndFsmMdcFS::detected(), PndFsmMdcTS::detected(), PndFsmMvd::detected(), PndFsmMvd2::detected(), PndFsmEmcBwCap::detected(), PndFsmMvdPid::detected(), PndFsmRich::detected(), PndFsmDrcDisc::detected(), PndFsmMdtPid::detected(), PndFsmEmcBarrel::detected(), PndFsmDrcBarrel::detected(), PndFsmEffTracker::detected(), PndFsmEmcPid::detected(), PndFsmDetTemplate::dp(), PndFsmSimpleTracker::dp(), PndFsmStt::dp(), PndFsmMdcFS::dp(), PndFsmMdcTS::dp(), PndFsmMvd::dp(), PndFsmMvd2::dp(), PndFsmEffTracker::dp(), PndFsmSimpleTracker::dphi(), PndFsmEmcFS::dphi(), PndFsmEmcFwCap::dphi(), PndFsmEmcBwCap::dphi(), PndFsmEffTracker::dphi(), PndFsmSimpleTracker::dtheta(), PndFsmEmcBwCap::dtheta(), PndFsmEmcFS::dtheta(), PndFsmEmcFwCap::dtheta(), PndFsmEmcBarrel::dtheta(), PndFsmEffTracker::dtheta(), PndFastSim::Exec(), HelixRep(), PndFsmCmpDet::respond(), PndFsmCombiDet::respond(), PndFsmMvd2::respond(), PndFsmEmcFS::respond(), PndFsmMvd::respond(), PndFsmStt::respond(), PndFsmEmcFwCap::respond(), PndFsmEmcBarrel::respond(), PndFsmSttPid::respond(), PndFsmEmcBwCap::respond(), PndFsmRich::respond(), PndFsmMvdPid::respond(), PndFsmTof::respond(), PndFsmMdtPid::respond(), PndFsmDrcBarrel::respond(), PndFsmDrcDisc::respond(), PndFsmEmcPid::respond(), PndFsmEffTracker::respond(), PndFastSim::SetFlatCovMatrix(), PndFastSim::smearEnergy(), PndFastSim::smearM(), PndFastSim::smearMomentum(), PndFastSim::smearPhi(), PndFastSim::smearTheta(), and PndFastSim::UpdateGammaHit().

72 {return _p4;} // 4-momentum
TLorentzVector _p4
Definition: PndFsmTrack.h:112
int PndFsmTrack::pdt ( )
inline
void PndFsmTrack::print ( std::ostream o)

Definition at line 366 of file PndFsmTrack.cxx.

References _charge, _gTrackId, _p4, _pdt, _startVtx, _stopVtx, GetHelixD0(), GetHelixOmega(), GetHelixPhi0(), GetHelixTanDip(), and GetHelixZ0().

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 }
Double_t GetHelixZ0() const
Definition: PndFsmTrack.h:140
double _charge
Definition: PndFsmTrack.h:115
Double_t GetHelixPhi0() const
Definition: PndFsmTrack.h:138
TVector3 _startVtx
Definition: PndFsmTrack.h:113
Double_t GetHelixOmega() const
Definition: PndFsmTrack.h:139
Double_t GetHelixD0() const
Definition: PndFsmTrack.h:137
signed long _gTrackId
Definition: PndFsmTrack.h:117
TLorentzVector _p4
Definition: PndFsmTrack.h:112
Double_t GetHelixTanDip() const
Definition: PndFsmTrack.h:141
TVector3 _stopVtx
Definition: PndFsmTrack.h:114
void PndFsmTrack::Propagate ( TVector3  origin,
double  deltaError = 2.5 
)

Definition at line 149 of file PndFsmTrack.cxx.

References _p4, _startVtx, a, alpha, c1, charge(), cos(), dz, fabs(), fCov5, fCov7, fPar5, fReference, GetHelixD0(), GetHelixOmega(), GetHelixPhi0(), GetHelixTanDip(), GetHelixZ0(), Mass, pnt, pt(), R, and sin().

Referenced by PndFastSim::cutAndSmear().

149  {
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 }
Double_t GetHelixZ0() const
Definition: PndFsmTrack.h:140
TMatrixD fCov7
Definition: PndFsmTrack.h:127
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
Double_t GetHelixPhi0() const
Definition: PndFsmTrack.h:138
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
TClonesArray * pnt
Double_t GetHelixD0() const
Definition: PndFsmTrack.h:137
double charge()
Definition: PndFsmTrack.h:75
TVector3 fReference
Definition: PndFsmTrack.h:125
TString pt(TString pts, TString exts="px py pz")
Definition: invexp.C:133
Int_t a
Definition: anaLmdDigi.C:126
TMatrixD fCov5
Definition: PndFsmTrack.h:126
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
double alpha
Definition: f_Init.h:9
Double_t GetHelixTanDip() const
Definition: PndFsmTrack.h:141
Double_t R
Definition: checkhelixhit.C:61
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
void PndFsmTrack::setCharge ( double  c)

Definition at line 315 of file PndFsmTrack.cxx.

References _charge, and c.

Referenced by PndFsmTrack().

316 {
317  _charge=c;
318 }
double _charge
Definition: PndFsmTrack.h:115
void PndFsmTrack::setDetResponse ( PndFsmResponse resp)

Definition at line 359 of file PndFsmTrack.cxx.

References _detResponse.

Referenced by PndFsmTrack(), and PndFastSim::smearTrack().

360 {
361  _detResponse=resp;
362 }
PndFsmResponse * _detResponse
Definition: PndFsmTrack.h:122
void PndFsmTrack::setGTrackId ( signed long  id)

Definition at line 352 of file PndFsmTrack.cxx.

References _gTrackId.

Referenced by PndFsmTrack().

353 {
354  _gTrackId=id;
355 }
signed long _gTrackId
Definition: PndFsmTrack.h:117
void PndFsmTrack::setMass2 ( double  c)

Definition at line 321 of file PndFsmTrack.cxx.

References _Mass2, and c.

Referenced by PndFsmTrack(), and PndFastSim::smearM2().

322 {
323  _Mass2=c;
324 }
double _Mass2
Definition: PndFsmTrack.h:118
void PndFsmTrack::setMvddEdX ( double  c)

Definition at line 328 of file PndFsmTrack.cxx.

References _MvddEdX, and c.

Referenced by PndFsmTrack(), and PndFastSim::smearMvddEdx().

329 {
330 _MvddEdX=c;
331 }
double _MvddEdX
Definition: PndFsmTrack.h:119
void PndFsmTrack::setP4 ( TLorentzVector  l)

Definition at line 297 of file PndFsmTrack.cxx.

References _p4.

Referenced by PndFsmTrack(), PndFastSim::smearEnergy(), PndFastSim::smearM(), PndFastSim::smearMomentum(), PndFastSim::smearPhi(), and PndFastSim::smearTheta().

298 {
299  _p4=l;
300 }
TLorentzVector _p4
Definition: PndFsmTrack.h:112
void PndFsmTrack::SetP4Cov ( TMatrixD p4cov)

Definition at line 264 of file PndFsmTrack.cxx.

References fCov7, and i.

Referenced by PndFastSim::SetFlatCovMatrix().

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 }
TMatrixD fCov7
Definition: PndFsmTrack.h:127
Int_t i
Definition: run_full.C:25
void PndFsmTrack::SetP7Cov ( TMatrixD p7cov)

Definition at line 255 of file PndFsmTrack.cxx.

References fCov7, and i.

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 }
TMatrixD fCov7
Definition: PndFsmTrack.h:127
Int_t i
Definition: run_full.C:25
void PndFsmTrack::setPdt ( int  pdt)

Definition at line 346 of file PndFsmTrack.cxx.

References _pdt, and n.

Referenced by PndFsmTrack().

347 {
348  _pdt=n;
349 }
int n
void PndFsmTrack::setStartVtx ( TVector3  v)

Definition at line 303 of file PndFsmTrack.cxx.

References _startVtx, and v.

Referenced by PndFsmTrack(), and PndFastSim::smearVertex().

304 {
305  _startVtx=v;
306 }
TVector3 _startVtx
Definition: PndFsmTrack.h:113
__m128 v
Definition: P4_F32vec4.h:4
void PndFsmTrack::setStopVtx ( TVector3  v)

Definition at line 309 of file PndFsmTrack.cxx.

References _stopVtx, and v.

Referenced by PndFsmTrack(), PndFsmEmcBwCap::respond(), PndFsmEmcFwCap::respond(), PndFsmEmcFS::respond(), PndFsmEmcBarrel::respond(), and PndFastSim::UpdateGammaHit().

310 {
311  _stopVtx=v;
312 }
__m128 v
Definition: P4_F32vec4.h:4
TVector3 _stopVtx
Definition: PndFsmTrack.h:114
void PndFsmTrack::setSttdEdX ( double  c)

Definition at line 340 of file PndFsmTrack.cxx.

References _SttdEdX, and c.

Referenced by PndFsmTrack(), and PndFastSim::smearSttdEdx().

341 {
342 _SttdEdX=c;
343 }
double _SttdEdX
Definition: PndFsmTrack.h:121
void PndFsmTrack::setTpcdEdX ( double  c)

Definition at line 334 of file PndFsmTrack.cxx.

References _TpcdEdX, and c.

Referenced by PndFsmTrack().

335 {
336 _TpcdEdX=c;
337 }
double _TpcdEdX
Definition: PndFsmTrack.h:120
void PndFsmTrack::SetVCov ( TMatrixD vcov)

Definition at line 273 of file PndFsmTrack.cxx.

References fCov7, and i.

Referenced by PndFastSim::SetFlatCovMatrix().

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 }
TMatrixD fCov7
Definition: PndFsmTrack.h:127
Int_t i
Definition: run_full.C:25
TVector3 PndFsmTrack::startVtx ( )
inline

Definition at line 73 of file PndFsmTrack.h.

References _startVtx.

Referenced by PndFastSim::cutAndSmear(), PndFastSim::Exec(), and PndFastSim::smearVertex().

73 {return _startVtx;}
TVector3 _startVtx
Definition: PndFsmTrack.h:113
TVector3 PndFsmTrack::stopVtx ( )
inline

Definition at line 74 of file PndFsmTrack.h.

References _stopVtx.

Referenced by PndFastSim::Exec(), and PndFastSim::UpdateGammaHit().

74 {return _stopVtx;}
TVector3 _stopVtx
Definition: PndFsmTrack.h:114
double PndFsmTrack::SttdEdX ( )
inline

Definition at line 82 of file PndFsmTrack.h.

References _SttdEdX.

82 {return _SttdEdX;}
double _SttdEdX
Definition: PndFsmTrack.h:121
double PndFsmTrack::TpcdEdX ( )
inline

Definition at line 81 of file PndFsmTrack.h.

References _TpcdEdX.

81 {return _TpcdEdX;}
double _TpcdEdX
Definition: PndFsmTrack.h:120

Member Data Documentation

double PndFsmTrack::_charge
private

Definition at line 115 of file PndFsmTrack.h.

Referenced by charge(), print(), and setCharge().

PndFsmResponse* PndFsmTrack::_detResponse
private

Definition at line 122 of file PndFsmTrack.h.

Referenced by detResponse(), setDetResponse(), and ~PndFsmTrack().

signed long PndFsmTrack::_gTrackId
private

Definition at line 117 of file PndFsmTrack.h.

Referenced by gTrackId(), print(), and setGTrackId().

double PndFsmTrack::_Mass2
private

Definition at line 118 of file PndFsmTrack.h.

Referenced by Mass2(), and setMass2().

double PndFsmTrack::_MvddEdX
private

Definition at line 119 of file PndFsmTrack.h.

Referenced by MvddEdX(), and setMvddEdX().

TLorentzVector PndFsmTrack::_p4
private

Definition at line 112 of file PndFsmTrack.h.

Referenced by p4(), print(), Propagate(), and setP4().

int PndFsmTrack::_pdt
private

Definition at line 116 of file PndFsmTrack.h.

Referenced by pdt(), print(), and setPdt().

TVector3 PndFsmTrack::_startVtx
private

Definition at line 113 of file PndFsmTrack.h.

Referenced by HelixRep(), print(), Propagate(), setStartVtx(), and startVtx().

TVector3 PndFsmTrack::_stopVtx
private

Definition at line 114 of file PndFsmTrack.h.

Referenced by print(), setStopVtx(), and stopVtx().

double PndFsmTrack::_SttdEdX
private

Definition at line 121 of file PndFsmTrack.h.

Referenced by setSttdEdX(), and SttdEdX().

double PndFsmTrack::_TpcdEdX
private

Definition at line 120 of file PndFsmTrack.h.

Referenced by setTpcdEdX(), and TpcdEdX().

TMatrixD PndFsmTrack::fCov5
private

Definition at line 126 of file PndFsmTrack.h.

Referenced by GetHelixCov(), PndFsmTrack(), and Propagate().

TMatrixD PndFsmTrack::fCov7
private

Definition at line 127 of file PndFsmTrack.h.

Referenced by Cov7(), PndFsmTrack(), Propagate(), SetP4Cov(), SetP7Cov(), and SetVCov().

double PndFsmTrack::fPar5[5]
private

Definition at line 124 of file PndFsmTrack.h.

Referenced by GetHelixParams(), HelixRep(), PndFsmTrack(), and Propagate().

TVector3 PndFsmTrack::fReference
private

Definition at line 125 of file PndFsmTrack.h.

Referenced by HelixRep(), and Propagate().


The documentation for this class was generated from the following files: