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

Class for saving a FTS track cand. for Hough transform based FTS PR. More...

#include <PndFtsHoughTrackCand.h>

Inheritance diagram for PndFtsHoughTrackCand:

Public Member Functions

 PndFtsHoughTrackCand (PndFtsHoughTrackerTask *trackerTask=0)
 Set pointer to tracker task (super important as it provides an I/O interface to PandaRoot) More...
 
 ~PndFtsHoughTrackCand ()
 
void Print ()
 
Bool_t isComplete () const
 
PndTrackCand getPndTrackCand ()
 
PndTrack getPndTrack ()
 
FairTrackParP getTrackParPForHit (const UInt_t i)
 
Int_t getCharge () const
 
TVector3 getP (const Double_t zLabSys) const
 
TVector3 getPos (const Double_t zLabSys) const
 
Double_t getXLabSys (const Double_t zLabSys) const
 
Double_t getThetaZxRad (const Double_t zLabSys) const
 
Double_t getThetaZyRad (const Double_t) const
 
Double_t getZLineParabola ()
 
Double_t getZParabolaLine ()
 
void SetZxLineBeforeDipole (const PndFtsHoughTracklet zxLineParabola)
 
void SetZxParabola (const PndFtsHoughTracklet zxParabola)
 
void SetZxLineBehindDipole (const PndFtsHoughTracklet zxParabolaLine)
 
void SetZyLine (const PndFtsHoughTracklet zyLine)
 
 ClassDef (PndFtsHoughTrackCand, 1)
 

Private Member Functions

void throwError (const TString s) const
 For error reporting. More...
 
void throwIfZOutOfRange (const Double_t &zLabSys) const
 
void addUniqueTrackletHits (const PndFtsHoughTracklet inTracklet)
 
Double_t getQdivPzx () const
 
Double_t getPYLab () const
 
std::pair< Double_t, Double_tgetPZPXLabLine (const Double_t &zLabSys, const PndFtsHoughTracklet *const lineTracklet) const
 
std::pair< Double_t, Double_tgetPZPXLabParabola (const Double_t &zLabSys) const
 
Double_t getXOrYLabForLine (const Double_t &zLabSys, const PndFtsHoughTracklet *const lineTracklet) const
 
Double_t getXLabForParabola (const Double_t &zLabSys) const
 

Private Attributes

PndFtsHoughTrackerTaskfTrackerTask
 
Int_t fVerbose
 
PndFtsHoughTracklet fZxLineBeforeDipole
 
PndFtsHoughTracklet fZxParabola
 
PndFtsHoughTracklet fZxLineBehindDipole
 
PndFtsHoughTracklet fZyLine
 
PndTrackCand fIntTrackCand
 
Double_t fZCoordLineParabola
 
Double_t fZCoordParabolaLine
 

Detailed Description

Class for saving a FTS track cand. for Hough transform based FTS PR.

Author
Martin J. Galuska <martin [dot] j [dot] galuska (at) physik [dot] uni (minus) giessen [dot] de>

The track cand. consists of several tracklets: In zx-plane line + parabola + line In zy-plane line

This class takes care of geometric calculations for tracklets.

Loosely modeled after tools/riemannfit/PndRiemannTrack.h

TODO Maybe it would be better not to derive from PndTrackCand, but to contain an object of PndTrackCand

Created: 24.01.2014

Definition at line 47 of file PndFtsHoughTrackCand.h.

Constructor & Destructor Documentation

PndFtsHoughTrackCand::PndFtsHoughTrackCand ( PndFtsHoughTrackerTask trackerTask = 0)

Set pointer to tracker task (super important as it provides an I/O interface to PandaRoot)

Definition at line 18 of file PndFtsHoughTrackCand.cxx.

References fTrackerTask, fVerbose, and PndFtsHoughTrackerTask::GetVerbose().

18  :
19  fTrackerTask(trackerTask),
20 
21  fVerbose(0),
22 
23  fZxLineBeforeDipole(0., trackerTask), // TODO: It could be a problem here that I set the z reference value to 0.
24 
25  fZxParabola(0., trackerTask),
26 
27  fZxLineBehindDipole(0., trackerTask),
28 
29  fZyLine(0., trackerTask),
30 
31  fIntTrackCand(),
32 
35 {
36  if (0==fTrackerTask){
37  std::cout << "PndFtsHoughTrackCand FATAL ERROR Tracker task pointer not set.\n";
38  } else {
40  if(3<fVerbose) std::cout << "PndFtsHoughTrackCand called with tracker ptr " << fTrackerTask << '\n';
41  }
42 }
PndFtsHoughTracklet fZxLineBehindDipole
PndFtsHoughTrackerTask * fTrackerTask
PndFtsHoughTracklet fZxParabola
PndFtsHoughTracklet fZxLineBeforeDipole
PndFtsHoughTracklet fZyLine
PndFtsHoughTrackCand::~PndFtsHoughTrackCand ( )

Definition at line 44 of file PndFtsHoughTrackCand.cxx.

45 {
46 
47 }

Member Function Documentation

void PndFtsHoughTrackCand::addUniqueTrackletHits ( const PndFtsHoughTracklet  inTracklet)
private

Definition at line 119 of file PndFtsHoughTrackCand.cxx.

References PndTrackCand::AddHit(), Double_t, fIntTrackCand, PndTrackCandHit::GetDetId(), PndTrackCandHit::GetHitId(), PndTrackCand::GetNHits(), PndTrackCandHit::GetRho(), PndTrackCand::GetSortedHit(), and PndTrackCand::HitInTrack().

Referenced by SetZxLineBeforeDipole(), SetZxLineBehindDipole(), SetZxParabola(), and SetZyLine().

120 {
121  // add all the hits from the tracklet which are not already in the hitId vector
122  for (UInt_t iHit = 0; iHit < inTracklet.GetNHits(); ++iHit)
123  {
124  PndTrackCandHit inHit = inTracklet.GetSortedHit(iHit);
125  const Int_t inHitId = inHit.GetHitId();
126  const Int_t inDetId = inHit.GetDetId();
127  // if hit is NOT in track -1 is returned by HitInTrack, otherwise the index (>=0) in the HitId vector is returned
128  if (-1==fIntTrackCand.HitInTrack(inDetId,inHitId)){
129  // add hit to track cand
130  const Double_t inRho = inHit.GetRho();
131  fIntTrackCand.AddHit(inDetId, inHitId, inRho);
132  }
133  }
134 }
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
void AddHit(UInt_t detId, UInt_t hitId, Double_t rho)
Double_t
Double_t GetRho() const
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
Int_t GetHitId() const
Int_t GetDetId() const
Int_t HitInTrack(UInt_t detId, UInt_t hitId) const
PndFtsHoughTrackCand::ClassDef ( PndFtsHoughTrackCand  ,
 
)
Int_t PndFtsHoughTrackCand::getCharge ( ) const
inline

Definition at line 163 of file PndFtsHoughTrackCand.h.

References getQdivPzx().

Referenced by getPZPXLabLine(), getPZPXLabParabola(), and getTrackParPForHit().

163  {
164  if (0 < getQdivPzx()){
165  return 1;
166  }
167  else{
168  return -1;
169  }
170 }
Double_t getQdivPzx() const
TVector3 PndFtsHoughTrackCand::getP ( const Double_t  zLabSys) const

Definition at line 265 of file PndFtsHoughTrackCand.cxx.

References Double_t, fVerbose, fZCoordLineParabola, fZCoordParabolaLine, fZxLineBeforeDipole, fZxLineBehindDipole, getPYLab(), getPZPXLabLine(), getPZPXLabParabola(), isComplete(), mom, and throwIfZOutOfRange().

Referenced by getTrackParPForHit().

265  {
266  throwIfZOutOfRange(zLabSys);
267 
268  TVector3 mom;
269 
270  if (kFALSE == isComplete())
271  {
272  std::cout << "getPforHit: Track cand. is not complete yet. Momentum will be calculated for incomplete track cand.\n";
273  }
274 
275  Double_t pZLabSys;
276  Double_t pXLabSys;
277  std::pair<Double_t, Double_t> pZPXLabSys;
278 
279  // track model is assumed to be line+parabola+line in zx and line in zy
280  if ( zLabSys <= fZCoordLineParabola ){
281  // use 1st line in zx plane
282  pZPXLabSys = getPZPXLabLine(zLabSys, &fZxLineBeforeDipole);
283  } else if ( zLabSys < fZCoordParabolaLine ){
284  // use "tangent" to parabola in zx plane
285  pZPXLabSys = getPZPXLabParabola(zLabSys);
286  } else {
287  // use 2nd line in zx plane
288  pZPXLabSys = getPZPXLabLine(zLabSys, &fZxLineBehindDipole);
289  }
290 
291  pZLabSys = pZPXLabSys.first;
292  pXLabSys = pZPXLabSys.second;
293  const Double_t pYLabSys = getPYLab();
294  mom.SetXYZ(pXLabSys, pYLabSys, pZLabSys);
295  if (fVerbose > 10) std::cout << "P-Vector for z=" << zLabSys << " : " << mom.X() << " " << mom.Y() << " " << mom.Z() << '\n';
296  return mom;
297 }
PndFtsHoughTracklet fZxLineBehindDipole
std::pair< Double_t, Double_t > getPZPXLabLine(const Double_t &zLabSys, const PndFtsHoughTracklet *const lineTracklet) const
Bool_t isComplete() const
std::pair< Double_t, Double_t > getPZPXLabParabola(const Double_t &zLabSys) const
Double_t mom
Definition: plot_dirc.C:14
Double_t
void throwIfZOutOfRange(const Double_t &zLabSys) const
PndFtsHoughTracklet fZxLineBeforeDipole
Double_t getPYLab() const
PndTrack PndFtsHoughTrackCand::getPndTrack ( )

Definition at line 208 of file PndFtsHoughTrackCand.cxx.

References fIntTrackCand, PndTrackCand::GetNHits(), and getTrackParPForHit().

208  {
209  // calculates first and last parameter, but uses an empty PndTrackCand which has to be set lateron using SetTrackCandRef
210  FairTrackParP firstPar, lastPar;
211  PndTrackCand myCand;
212  if (fIntTrackCand.GetNHits() > 0){
213  firstPar = getTrackParPForHit(0);
215  }
216  return PndTrack(firstPar, lastPar, myCand);
217 }
FairTrackParP getTrackParPForHit(const UInt_t i)
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
PndTrackCand PndFtsHoughTrackCand::getPndTrackCand ( )

Definition at line 184 of file PndFtsHoughTrackCand.cxx.

References fIntTrackCand.

184  {
185  return fIntTrackCand;
186 
187  // // The following is what I used when I derived this class from PndTrackCand. Instead this class now contains a member of PndTrackCand
188  // // Maybe this is not necessary because PndFtsHoughTrackCand is derived from PndTrackCand
189  // PndTrackCand myCand;
190  // // copy all the hits from PndFtsHoughTrackCand *this to PndTrackCand myCand
191  // for (UInt_t iHit = 0; iHit < intTrackCand.GetNHits(); ++iHit)
192  // {
193  // PndTrackCandHit inHit = intTrackCand.GetSortedHit(iHit);
194  // const Int_t inHitId = inHit.GetHitId();
195  // const Int_t inDetId = inHit.GetDetId();
196  // const Double_t inRho = inHit.GetRho();
197  // myCand.AddHit(inDetId, inHitId, inRho);
198  // }
199  // return myCand;
200 }
TVector3 PndFtsHoughTrackCand::getPos ( const Double_t  zLabSys) const

Definition at line 331 of file PndFtsHoughTrackCand.cxx.

References Double_t, fVerbose, fZyLine, getXLabSys(), getXOrYLabForLine(), isComplete(), and throwIfZOutOfRange().

Referenced by getTrackParPForHit(), and plotHoughTrackCand().

331  {
332  // calculates the point on the track based on the results from the Hough transforms using my track model for the given z
333 
334  throwIfZOutOfRange(zLabSys);
335 
336  if (kFALSE == isComplete()) std::cout << "getPositionForHit: Track cand. is not complete yet. Position will be calculated for incomplete track cand.\n";
337 
338  TVector3 position;
339  // track model is assumed to be line+parabola+line in zx and line in zy
340  const Double_t yLabSys = getXOrYLabForLine(zLabSys, &fZyLine);
341 
342  Double_t xLabSys = getXLabSys(zLabSys);
343 
344  position.SetXYZ(xLabSys, yLabSys, zLabSys);
345  if (fVerbose > 10) std::cout << "Pos-Vector for z=" << zLabSys << " : " << position.X() << " " << position.Y() << " " << position.Z() << '\n';
346  return position;
347 }
Bool_t isComplete() const
Double_t getXOrYLabForLine(const Double_t &zLabSys, const PndFtsHoughTracklet *const lineTracklet) const
Double_t
void throwIfZOutOfRange(const Double_t &zLabSys) const
Double_t getXLabSys(const Double_t zLabSys) const
PndFtsHoughTracklet fZyLine
Double_t PndFtsHoughTrackCand::getPYLab ( ) const
inlineprivate

Definition at line 217 of file PndFtsHoughTrackCand.h.

References Double_t, fZxLineBeforeDipole, fZxParabola, fZyLine, getPZPXLabLine(), PndFtsHoughTracklet::getThetaRadVal(), and PndFtsHoughTracklet::getZRefLabSys().

Referenced by getP().

217  {
218  // theta in radian in zy plane given at z where first zx line meets the parabola
219  const Double_t zRefLabSys = fZxParabola.getZRefLabSys();
220  const Double_t thetaZyRad = fZyLine.getThetaRadVal();
221  std::pair<Double_t, Double_t> pZPXLabSys = getPZPXLabLine(zRefLabSys, &fZxLineBeforeDipole);
222  const Double_t pZLabSys = pZPXLabSys.first;
223  Double_t pYLabSys = tan(thetaZyRad)*pZLabSys;
224  return pYLabSys;
225  }
Double_t getThetaRadVal() const
std::pair< Double_t, Double_t > getPZPXLabLine(const Double_t &zLabSys, const PndFtsHoughTracklet *const lineTracklet) const
Double_t getZRefLabSys() const
PndFtsHoughTracklet fZxParabola
Double_t
PndFtsHoughTracklet fZxLineBeforeDipole
PndFtsHoughTracklet fZyLine
std::pair< Double_t, Double_t > PndFtsHoughTrackCand::getPZPXLabLine ( const Double_t zLabSys,
const PndFtsHoughTracklet *const  lineTracklet 
) const
inlineprivate

Definition at line 203 of file PndFtsHoughTrackCand.h.

References cos(), Double_t, fZxParabola, getCharge(), PndFtsHoughTracklet::getSecondVal(), PndFtsHoughTracklet::getThetaRadVal(), and sin().

Referenced by getP(), and getPYLab().

203  { // zLabSys //[R.K.03/2017] unused variable(s)
204  // theta in radian in zx plane given at z = zRefLabSys
205  const Double_t thetaRad = lineTracklet->getThetaRadVal();
206  const Double_t qDivPzx = fZxParabola.getSecondVal(); // Q/pzx
207  const Double_t pZx = getCharge() / qDivPzx;
208 
209  const Double_t pZLabSys = pZx*cos(thetaRad);
210  const Double_t pXLabSys = pZx*sin(thetaRad);
211  std::pair<Double_t, Double_t> pZPXLabSys(pZLabSys, pXLabSys);
212 
213  return pZPXLabSys;
214  }
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
Double_t getThetaRadVal() const
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
PndFtsHoughTracklet fZxParabola
Double_t
Double_t getSecondVal() const
std::pair< Double_t, Double_t > PndFtsHoughTrackCand::getPZPXLabParabola ( const Double_t zLabSys) const
inlineprivate

Definition at line 188 of file PndFtsHoughTrackCand.h.

References cos(), Double_t, getCharge(), getQdivPzx(), getThetaZxRad(), and sin().

Referenced by getP().

188  {
189 
190  const Double_t currentThetaRad = getThetaZxRad(zLabSys);
191 
192  const Double_t qDivPzx = getQdivPzx(); // Q/pzx
193  const Double_t pZx = getCharge() / qDivPzx;
194 
195  const Double_t pZLabSys = pZx*cos(currentThetaRad);
196  const Double_t pXLabSys = pZx*sin(currentThetaRad);
197  std::pair<Double_t, Double_t> pZPXLabSys(pZLabSys, pXLabSys);
198 
199  return pZPXLabSys;
200 }
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
Double_t
Double_t getQdivPzx() const
Double_t getThetaZxRad(const Double_t zLabSys) const
Double_t PndFtsHoughTrackCand::getQdivPzx ( ) const
inlineprivate

Definition at line 103 of file PndFtsHoughTrackCand.h.

References fZxParabola, and PndFtsHoughTracklet::getSecondVal().

Referenced by getCharge(), getPZPXLabParabola(), and getXLabForParabola().

103 { return fZxParabola.getSecondVal(); };
PndFtsHoughTracklet fZxParabola
Double_t getSecondVal() const
Double_t PndFtsHoughTrackCand::getThetaZxRad ( const Double_t  zLabSys) const

Definition at line 300 of file PndFtsHoughTrackCand.cxx.

References atan2(), Double_t, getXLabSys(), and throwIfZOutOfRange().

Referenced by getPZPXLabParabola().

300  {
301  // estimate angle by calculating 2 points (zLabSys1, xLabSys1) and (zLabSys2, xLabSys2) which are close together
302 
303  throwIfZOutOfRange(zLabSys);
304 
305  const Double_t epsilon = 0.5; // small positive number
306 
307  const Double_t xLabSys1 = getXLabSys(zLabSys-epsilon);
308  const Double_t xLabSys2 = getXLabSys(zLabSys+epsilon);
309 
310  return atan2(xLabSys2-xLabSys1, 2*epsilon);
311 }
Double_t
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
void throwIfZOutOfRange(const Double_t &zLabSys) const
Double_t getXLabSys(const Double_t zLabSys) const
Double_t PndFtsHoughTrackCand::getThetaZyRad ( const Double_t  ) const
inline

Definition at line 75 of file PndFtsHoughTrackCand.h.

References fZyLine, and PndFtsHoughTracklet::getThetaRadVal().

Referenced by PndFtsHoughTrackFinder::LineBehindDipoleMatchesToLinePlusParabola().

75 { return fZyLine.getThetaRadVal(); }; // zLabSys //[R.K.03/2017] unused variable(s)
Double_t getThetaRadVal() const
PndFtsHoughTracklet fZyLine
FairTrackParP PndFtsHoughTrackCand::getTrackParPForHit ( const UInt_t  i)

Definition at line 220 of file PndFtsHoughTrackCand.cxx.

References Double_t, fIntTrackCand, fTrackerTask, getCharge(), PndFtsHoughTrackerTask::GetFtsHit(), PndFtsHoughTrackerTask::GetFtsHitPosErrors(), PndTrackCandHit::GetHitId(), getP(), getPos(), PndTrackCand::GetSortedHit(), isComplete(), mom, and throwError().

Referenced by getPndTrack().

220  {
221  // index is the index of the hit in intTrackCand, not in the FTS hit array
222  // this method will sort the hitId vector of intTrackCand and therefore cannot be used for const track candidates
223 
224  // TODO: Check if all arguments are correct
225 
226  // Warn if we do not have a complete track candidate
227  if (!isComplete()) std::cout << "getHit: You try to access hits before we have a complete track candidate.\n";
228 
229  // not necessary, I check later if I get a hit or not
230  //if (index >= intTrackCand.GetNHits()) throwError("index in PndFtsHoughTrackCand::getTrackParPForHit is too large.");
231 
232  // Translate index in track candidate to hitId in FTS hit array
233  UInt_t hitId = fIntTrackCand.GetSortedHit(index).GetHitId();
234 
235  // get position of hit
236  const PndFtsHit *myHit = fTrackerTask->GetFtsHit(hitId);
237  if (0 == myHit) throwError("Cannot get hit, probably the tracking has not finished or the index is too large.");
238 
239  // Take z of hit and calculate the position and momentum for that z using the results from the Hough transforms and my track model
240  Double_t zLabSys = myHit->GetZ();
241 
242  // position should NOT come from hit, it should come from the pattern recognition track model
243  // position error can be large (like 1* or 2* tube size) as the Kalman filter will adjust it.
244  TVector3 hitPos = getPos(zLabSys);
245 
246 
247  TVector3 hitPosError = fTrackerTask->GetFtsHitPosErrors(myHit);
248 
249  // momentum comes from the pattern recognition track model
250  TVector3 mom = getP(zLabSys);
251  TVector3 momError = 0.1*mom; // TODO: add correct values here
252 
253  // set plane as detector plane in which the hit is, FTS planes are parallel to xy plane
254  TVector3 origin; // set origin of plane = position of hit
255  myHit->Position(origin);
256  TVector3 dj(1,0,0); // unit vector along x // TODO: Check if that is set correctly for FTS
257  TVector3 dk(0,1,0); // unit vector along y // TODO: Check if that is set correctly for FTS
258 
259  // ----- Constructor with track parameters in LAB -----------------------------------
260  // FairTrackParP::FairTrackParP(TVector3 pos, TVector3 Mom, TVector3 posErr, TVector3 MomErr, Int_t Q, TVector3 o, TVector3 dj, TVector3 dk)
261  FairTrackParP result(hitPos, mom, hitPosError, momError, getCharge(), origin, dj, dk);
262  return result;
263 }
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
Bool_t isComplete() const
PndFtsHoughTrackerTask * fTrackerTask
Double_t mom
Definition: plot_dirc.C:14
TVector3 getPos(const Double_t zLabSys) const
void throwError(const TString s) const
For error reporting.
Double_t
TVector3 getP(const Double_t zLabSys) const
const PndFtsHit * GetFtsHit(UInt_t hitId) const
Returns pointer to the hit with index hitId in the FTS hit array.
Int_t GetHitId() const
const TVector3 GetFtsHitPosErrors(const PndFtsHit *const ftsHit) const
Returns the position error (based on FTS straw geometry) for the hit with index hitId in the FTS hit ...
Double_t PndFtsHoughTrackCand::getXLabForParabola ( const Double_t zLabSys) const
private

Definition at line 138 of file PndFtsHoughTrackCand.cxx.

References a, Double_t, fabs(), fZxLineBeforeDipole, fZxParabola, getQdivPzx(), PndFtsHoughTracklet::getSecondVal(), PndFtsHoughTracklet::getThetaRadVal(), PndFtsHoughTracklet::getZRefLabSys(), sin(), and sqrt().

Referenced by getXLabSys().

138  {
139  // does not take option of varying B field into account
140  // calculate x in lab sys for a given z position in lab sys for which the parabola assumption holds
141  // theta in radian in zx plane given at z = zRefLabSys
142  const Double_t thetaRad = fZxParabola.getThetaRadVal();
143  const Double_t qDivPzx = getQdivPzx(); // Q/pzx
144  const Double_t zRefLabSys = fZxParabola.getZRefLabSys();
145  const Double_t zshifted = zLabSys-zRefLabSys;
146 
147  const Double_t interceptLine = fZxLineBeforeDipole.getSecondVal();
148 
149  const Double_t By = 1.;
150 
151  if (0==thetaRad){
152  Double_t xParabolaSys = qDivPzx /2. * By * zshifted;
153  // go from local parabola system to the lab system
154  Double_t xLabSys = xParabolaSys + interceptLine;
155  return xLabSys;
156  }
157 
158 
159  const Double_t tantheta = tan(thetaRad);
160  const Double_t pzstuff = 1. / qDivPzx / By / sin(thetaRad);
161  const Double_t ztantheta = zshifted / tantheta;
162  const Double_t a = -ztantheta + pzstuff / tantheta;
163 
164  const Double_t wurzel = sqrt(a * a - 2. * pzstuff * zshifted - ztantheta * ztantheta);
165  const Double_t x1ParabolaSys = a - wurzel;
166  const Double_t x2ParabolaSys = a + wurzel;
167 
168  // go from local parabola system to the lab system
169  Double_t x1LabSys = x1ParabolaSys + interceptLine;
170  Double_t x2LabSys = x2ParabolaSys + interceptLine;
171 
172  // there are two analytical solutions, take the one closer to 0
173  if (fabs(x1LabSys)<fabs(x2LabSys)) {
174  return x1LabSys;
175  } else {
176  return x2LabSys;
177  }
178 }
Double_t getThetaRadVal() const
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
Double_t getZRefLabSys() const
PndFtsHoughTracklet fZxParabola
Int_t a
Definition: anaLmdDigi.C:126
Double_t
Double_t getSecondVal() const
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
PndFtsHoughTracklet fZxLineBeforeDipole
Double_t getQdivPzx() const
Double_t PndFtsHoughTrackCand::getXLabSys ( const Double_t  zLabSys) const

Definition at line 313 of file PndFtsHoughTrackCand.cxx.

References fZCoordLineParabola, fZCoordParabolaLine, fZxLineBeforeDipole, fZxLineBehindDipole, getXLabForParabola(), getXOrYLabForLine(), throwError(), and throwIfZOutOfRange().

Referenced by PndFtsHoughSpace::CalculateHitPosFromIntersectionsWithZxTrackModel(), getPos(), and getThetaZxRad().

313  {
314  // calculates the corresponding x coordinate to a given z coordinate according to the track model
315 
316  throwIfZOutOfRange(zLabSys);
317 
318  if (zLabSys <= fZCoordLineParabola) {
319  // use 1st line in zx plane
320  return getXOrYLabForLine(zLabSys, &fZxLineBeforeDipole);
321  } else if (zLabSys < fZCoordParabolaLine) {
322  // use parabola in zx plane
323  return getXLabForParabola(zLabSys);
324  } else {
325  // use 2nd line in zx plane
326  return getXOrYLabForLine(zLabSys, &fZxLineBehindDipole);
327  }
328  throwError("Argument was not handled!");
329 }
PndFtsHoughTracklet fZxLineBehindDipole
Double_t getXLabForParabola(const Double_t &zLabSys) const
void throwError(const TString s) const
For error reporting.
Double_t getXOrYLabForLine(const Double_t &zLabSys, const PndFtsHoughTracklet *const lineTracklet) const
void throwIfZOutOfRange(const Double_t &zLabSys) const
PndFtsHoughTracklet fZxLineBeforeDipole
Double_t PndFtsHoughTrackCand::getXOrYLabForLine ( const Double_t zLabSys,
const PndFtsHoughTracklet *const  lineTracklet 
) const
inlineprivate

Definition at line 174 of file PndFtsHoughTrackCand.h.

References Double_t, PndFtsHoughTracklet::getSecondVal(), PndFtsHoughTracklet::getThetaRadVal(), and PndFtsHoughTracklet::getZRefLabSys().

Referenced by getPos(), and getXLabSys().

174  {
175  // calculate x or y in lab sys for a given z position in lab sys for which the line assumption holds
176  // theta in radian in zx plane given at z = zRefLabSys
177  const Double_t thetaRad = lineTracklet->getThetaRadVal();
178  const Double_t intercept = lineTracklet->getSecondVal();
179  const Double_t zRefLabSys = lineTracklet->getZRefLabSys();
180  const Double_t zshifted = zLabSys-zRefLabSys;
181 
182  Double_t xOrYLabSys = tan(thetaRad)*zshifted+intercept;
183 
184  return xOrYLabSys;
185 }
Double_t getThetaRadVal() const
Double_t getZRefLabSys() const
Double_t
Double_t getSecondVal() const
Double_t PndFtsHoughTrackCand::getZLineParabola ( )
inline

Definition at line 77 of file PndFtsHoughTrackCand.h.

References fZCoordLineParabola.

Referenced by plotHoughTrackCand().

77 { return fZCoordLineParabola; }; // gets z coordinate in laboratory system where I switch from line before dipole to parabola within dipole (in zx plane)
Double_t PndFtsHoughTrackCand::getZParabolaLine ( )
inline

Definition at line 78 of file PndFtsHoughTrackCand.h.

References fZCoordParabolaLine.

Referenced by plotHoughTrackCand().

78 { return fZCoordParabolaLine; }; // gets z coordinate in laboratory system where I switch from parabola within dipole to line behind dipole (in zx plane)
Bool_t PndFtsHoughTrackCand::isComplete ( ) const
inline

Definition at line 60 of file PndFtsHoughTrackCand.h.

References fZxLineBeforeDipole, fZxLineBehindDipole, fZxParabola, fZyLine, and PndFtsHoughTracklet::isSet().

Referenced by getP(), getPos(), and getTrackParPForHit().

PndFtsHoughTracklet fZxLineBehindDipole
PndFtsHoughTracklet fZxParabola
PndFtsHoughTracklet fZxLineBeforeDipole
PndFtsHoughTracklet fZyLine
void PndFtsHoughTrackCand::Print ( )

Definition at line 94 of file PndFtsHoughTrackCand.cxx.

References fIntTrackCand, fZxLineBeforeDipole, fZxLineBehindDipole, fZxParabola, fZyLine, PndFtsHoughTracklet::isSet(), PndFtsHoughTracklet::Print(), and PndTrackCand::Print().

Referenced by PndFtsHoughTrackerTask::Exec().

94  {
95  std::cout << "=========== PndFtsHoughTrackCand::Print() ==========\n";
97  if (kTRUE==fZxLineBeforeDipole.isSet()){
98  std::cout << "zx plane\n\n";
99  std::cout << "1st line: ";
101  }
102  if (kTRUE==fZxParabola.isSet()){
103  std::cout << "parabola: ";
104  fZxParabola.Print();
105  }
106  if (kTRUE==fZxLineBehindDipole.isSet()){
107  std::cout << "2nd line: ";
109  }
110  if (kTRUE==fZyLine.isSet()){
111  std::cout << "zy plane\n\n";
112  std::cout << "Line: ";
113  fZyLine.Print();
114  }
115 }
void Print() const
PndFtsHoughTracklet fZxLineBehindDipole
PndFtsHoughTracklet fZxParabola
PndFtsHoughTracklet fZxLineBeforeDipole
PndFtsHoughTracklet fZyLine
void PndFtsHoughTrackCand::SetZxLineBeforeDipole ( const PndFtsHoughTracklet  zxLineParabola)

Definition at line 50 of file PndFtsHoughTrackCand.cxx.

References addUniqueTrackletHits(), fZCoordLineParabola, fZxLineBeforeDipole, PndFtsHoughTracklet::getZRefLabSys(), and PndFtsHoughTracklet::isSet().

Referenced by PndFtsHoughTrackFinder::FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole().

50  {
51  if (fZxLineBeforeDipole.isSet()) {
52  std::cerr << "WARNING from PndFtsHoughTrackCand: Line before dipole in zx plane is set more than once! Ignore new values! Potentially FATAL ERROR!\n";
53  return;
54  }
55  fZxLineBeforeDipole = zxLineParabola;
58 }
void addUniqueTrackletHits(const PndFtsHoughTracklet inTracklet)
Double_t getZRefLabSys() const
PndFtsHoughTracklet fZxLineBeforeDipole
void PndFtsHoughTrackCand::SetZxLineBehindDipole ( const PndFtsHoughTracklet  zxParabolaLine)

Definition at line 71 of file PndFtsHoughTrackCand.cxx.

References addUniqueTrackletHits(), fZCoordParabolaLine, fZxLineBehindDipole, PndFtsHoughTracklet::getZRefLabSys(), and PndFtsHoughTracklet::isSet().

Referenced by PndFtsHoughTrackFinder::FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole().

71  {
72  if (fZxLineBehindDipole.isSet()) {
73  std::cerr << "WARNING from PndFtsHoughTrackCand: Line behind dipole in zx plane is set more than once! Ignore new values! Potentially FATAL ERROR!\n";
74  return;
75  }
76  fZxLineBehindDipole = zxParabolaLine;
79 }
PndFtsHoughTracklet fZxLineBehindDipole
void addUniqueTrackletHits(const PndFtsHoughTracklet inTracklet)
Double_t getZRefLabSys() const
void PndFtsHoughTrackCand::SetZxParabola ( const PndFtsHoughTracklet  zxParabola)

Definition at line 59 of file PndFtsHoughTrackCand.cxx.

References addUniqueTrackletHits(), fZCoordLineParabola, fZxParabola, PndFtsHoughTracklet::getZRefLabSys(), and PndFtsHoughTracklet::isSet().

Referenced by PndFtsHoughTrackFinder::FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole().

59  {
60  if (fZxParabola.isSet()) {
61  std::cerr << "WARNING from PndFtsHoughTrackCand: Parabola inside dipole in zx plane is set more than once! Ignore new values! Potentially FATAL ERROR!\n";
62  return;
63  }
64  // warn if line before dipole field and parabola within are not calculated wrt the same z reference value
65  if (zxParabola.getZRefLabSys() != fZCoordLineParabola){
66  std::cout << "WARNING from PndFtsHoughTrackCand: First line and parabola were not calculated wrt the same z position! Potentially FATAL ERROR!\n";
67  }
68  fZxParabola = zxParabola;
70 }
void addUniqueTrackletHits(const PndFtsHoughTracklet inTracklet)
Double_t getZRefLabSys() const
PndFtsHoughTracklet fZxParabola
void PndFtsHoughTrackCand::SetZyLine ( const PndFtsHoughTracklet  zyLine)

Definition at line 80 of file PndFtsHoughTrackCand.cxx.

References addUniqueTrackletHits(), fZyLine, and PndFtsHoughTracklet::isSet().

Referenced by PndFtsHoughTrackFinder::FindZyLineMatchingToLineParabolaLineInZx().

80  {
81  if (fZyLine.isSet()) {
82  std::cerr << "WARNING from PndFtsHoughTrackCand: Line in zy plane is set more than once! Ignore new values! Potentially FATAL ERROR!\n";
83  return;
84  }
85  fZyLine=zyLine;
87 }
void addUniqueTrackletHits(const PndFtsHoughTracklet inTracklet)
PndFtsHoughTracklet fZyLine
void PndFtsHoughTrackCand::throwError ( const TString  s) const
inlineprivate

For error reporting.

Definition at line 95 of file PndFtsHoughTrackCand.h.

Referenced by getTrackParPForHit(), getXLabSys(), and throwIfZOutOfRange().

95 { throw std::runtime_error(s.Data()); };
TLorentzVector s
Definition: Pnd2DStar.C:50
void PndFtsHoughTrackCand::throwIfZOutOfRange ( const Double_t zLabSys) const
inlineprivate

Definition at line 96 of file PndFtsHoughTrackCand.h.

References throwError().

Referenced by getP(), getPos(), getThetaZxRad(), and getXLabSys().

96  {
97  if ( zLabSys <= 100 ) throwError("zLabSys is too small, track model is not valid in that region.");
98  if ( zLabSys >= 1000 ) throwError("zLabSys is too big, track model is not valid in that region.");
99  };
void throwError(const TString s) const
For error reporting.

Member Data Documentation

PndTrackCand PndFtsHoughTrackCand::fIntTrackCand
private
PndFtsHoughTrackerTask* PndFtsHoughTrackCand::fTrackerTask
private

Definition at line 92 of file PndFtsHoughTrackCand.h.

Referenced by getTrackParPForHit(), and PndFtsHoughTrackCand().

Int_t PndFtsHoughTrackCand::fVerbose
private

Definition at line 115 of file PndFtsHoughTrackCand.h.

Referenced by getP(), getPos(), and PndFtsHoughTrackCand().

Double_t PndFtsHoughTrackCand::fZCoordLineParabola
private
Double_t PndFtsHoughTrackCand::fZCoordParabolaLine
private

Definition at line 150 of file PndFtsHoughTrackCand.h.

Referenced by getP(), getXLabSys(), getZParabolaLine(), and SetZxLineBehindDipole().

PndFtsHoughTracklet PndFtsHoughTrackCand::fZxLineBeforeDipole
private
PndFtsHoughTracklet PndFtsHoughTrackCand::fZxLineBehindDipole
private

Definition at line 131 of file PndFtsHoughTrackCand.h.

Referenced by getP(), getXLabSys(), isComplete(), Print(), and SetZxLineBehindDipole().

PndFtsHoughTracklet PndFtsHoughTrackCand::fZxParabola
private
PndFtsHoughTracklet PndFtsHoughTrackCand::fZyLine
private

Definition at line 138 of file PndFtsHoughTrackCand.h.

Referenced by getPos(), getPYLab(), getThetaZyRad(), isComplete(), Print(), and SetZyLine().


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