LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
genf::GFDetPlane Class Reference

#include "GFDetPlane.h"

Inheritance diagram for genf::GFDetPlane:

Public Member Functions

 GFDetPlane (genf::GFAbsFinitePlane *finite=NULL)
 
 GFDetPlane (const TVector3 &o, const TVector3 &u, const TVector3 &v, genf::GFAbsFinitePlane *finite=NULL)
 
 GFDetPlane (const TVector3 &o, const TVector3 &n, genf::GFAbsFinitePlane *finite=NULL)
 
virtual ~GFDetPlane ()
 
 GFDetPlane (const genf::GFDetPlane &)
 
GFDetPlaneoperator= (const genf::GFDetPlane &)
 
TVector3 getO () const
 
TVector3 getU () const
 
TVector3 getV () const
 
void set (const TVector3 &o, const TVector3 &u, const TVector3 &v)
 
void setO (const TVector3 &o)
 
void setO (double, double, double)
 
void setU (const TVector3 &u)
 
void setU (double, double, double)
 
void setV (const TVector3 &v)
 
void setV (double, double, double)
 
void setUV (const TVector3 &u, const TVector3 &v)
 
void setON (const TVector3 &o, const TVector3 &n)
 
void setFinitePlane (genf::GFAbsFinitePlane *finite)
 
TVector3 getNormal () const
 
void setNormal (TVector3 n)
 
void setNormal (double, double, double)
 
void setNormal (const double &theta, const double &phi)
 
TVector2 project (const TVector3 &x) const
 projecting a direction onto the plane: More...
 
TVector2 LabToPlane (const TVector3 &x) const
 transform from Lab system into plane More...
 
TVector3 toLab (const TVector2 &x) const
 transform from plane coordinates to lab system More...
 
TVector3 dist (const TVector3 &point) const
 
TVector2 straightLineToPlane (const TVector3 &point, const TVector3 &dir) const
 gives u,v coordinates of the intersection point of a straight line with plane More...
 
void Print (std::ostream &out=std::cout) const
 
void getGraphics (double mesh, double length, TPolyMarker3D **pl, TPolyLine3D **plLine, TPolyLine3D **u, TPolyLine3D **v, TPolyLine3D **n=NULL)
 for poor attempts of making an event display. There is a lot of room for improvements. More...
 
double distance (TVector3 &) const
 
double distance (double, double, double) const
 
bool inActive (const TVector3 &point, const TVector3 &dir) const
 intersect in the active area? C.f. GFAbsFinitePlane More...
 
bool inActive (double u, double v) const
 inActive methods refer to finite plane. C.f. GFAbsFinitePlane More...
 
bool inActive (const TVector2 &v) const
 inActive methods refer to finite plane. C.f. GFAbsFinitePlane More...
 
void sane ()
 

Public Attributes

TVector3 fO
 
TVector3 fU
 
TVector3 fV
 
genf::GFAbsFinitePlanefFinitePlane
 

Private Member Functions

virtual void Print (Option_t *) const
 

Friends

bool operator== (const GFDetPlane &, const GFDetPlane &)
 
bool operator!= (const GFDetPlane &, const GFDetPlane &)
 returns NOT == More...
 

Detailed Description

Definition at line 64 of file GFDetPlane.h.

Constructor & Destructor Documentation

genf::GFDetPlane::GFDetPlane ( genf::GFAbsFinitePlane finite = NULL)

Definition at line 42 of file GFDetPlane.cxx.

References fO, fU, fV, r, and sane().

42  : fFinitePlane(finite)
43 {
44  static TRandom3 r(0);
45  fO.SetXYZ(0., 0., 0.);
46  fU.SetXYZ(r.Uniform(), r.Uniform(), 0.);
47  fV.SetXYZ(r.Uniform(), r.Uniform(), 0.);
48  sane();
49 }
TRandom r
Definition: spectrum.C:23
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
genf::GFDetPlane::GFDetPlane ( const TVector3 &  o,
const TVector3 &  u,
const TVector3 &  v,
genf::GFAbsFinitePlane finite = NULL 
)

Definition at line 34 of file GFDetPlane.cxx.

References sane().

38  : fO(o), fU(u), fV(v), fFinitePlane(finite)
39 {
40  sane();
41 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
genf::GFDetPlane::GFDetPlane ( const TVector3 &  o,
const TVector3 &  n,
genf::GFAbsFinitePlane finite = NULL 
)

Definition at line 51 of file GFDetPlane.cxx.

References setNormal().

52  : fO(o), fFinitePlane(finite)
53 {
54  setNormal(n);
55 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:151
Char_t n[5]
genf::GFDetPlane::~GFDetPlane ( )
virtual

Definition at line 57 of file GFDetPlane.cxx.

References fFinitePlane.

58 {
59  if (fFinitePlane != NULL) delete fFinitePlane;
60 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
genf::GFDetPlane::GFDetPlane ( const genf::GFDetPlane rhs)

Definition at line 62 of file GFDetPlane.cxx.

References genf::GFAbsFinitePlane::clone(), fFinitePlane, fO, fU, and fV.

62  : TObject(rhs)
63 {
64  if (rhs.fFinitePlane != NULL)
66  else
67  fFinitePlane = NULL;
68  fO = rhs.fO;
69  fU = rhs.fU;
70  fV = rhs.fV;
71 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
virtual GFAbsFinitePlane * clone() const =0
Deep copy ctor for polymorphic class.

Member Function Documentation

TVector3 genf::GFDetPlane::dist ( const TVector3 &  point) const

Definition at line 199 of file GFDetPlane.cxx.

References LabToPlane(), toLab(), and x.

Referenced by genf::GFWireHitPolicy::detPlane(), genf::GFWirepointHitPolicy::detPlane(), genf::RKTrackRep::Extrap(), and setFinitePlane().

200 {
201  TVector2 p = LabToPlane(x);
202  TVector3 xplane = toLab(p);
203  return xplane - x;
204 }
Float_t x
Definition: compare.C:6
TVector3 toLab(const TVector2 &x) const
transform from plane coordinates to lab system
Definition: GFDetPlane.cxx:191
TVector2 LabToPlane(const TVector3 &x) const
transform from Lab system into plane
Definition: GFDetPlane.cxx:185
double genf::GFDetPlane::distance ( TVector3 &  v) const

Definition at line 340 of file GFDetPlane.cxx.

References fO, fU, and fV.

Referenced by genf::RKTrackRep::Extrap(), genf::RKTrackRep::RKutta(), and setFinitePlane().

341 {
342  double s = (v - fO) * fU;
343  double t = (v - fO) * fV;
344  TVector3 distanceVector = v - fO - (s * fU) - (t * fV);
345  return distanceVector.Mag();
346 }
double genf::GFDetPlane::distance ( double  x,
double  y,
double  z 
) const

Definition at line 347 of file GFDetPlane.cxx.

References fO, fU, and fV.

348 {
349  TVector3 v(x, y, z);
350  double s = (v - fO) * fU;
351  double t = (v - fO) * fV;
352  TVector3 distanceVector = v - fO - (s * fU) - (t * fV);
353  return distanceVector.Mag();
354 }
Float_t x
Definition: compare.C:6
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:276
void genf::GFDetPlane::getGraphics ( double  mesh,
double  length,
TPolyMarker3D **  pl,
TPolyLine3D **  plLine,
TPolyLine3D **  u,
TPolyLine3D **  v,
TPolyLine3D **  n = NULL 
)

for poor attempts of making an event display. There is a lot of room for improvements.

Definition at line 271 of file GFDetPlane.cxx.

References fO, fU, fV, and getNormal().

Referenced by setFinitePlane().

278 {
279  *pl = new TPolyMarker3D(21 * 21, 24);
280  (*pl)->SetMarkerSize(0.1);
281  (*pl)->SetMarkerColor(kBlue);
282  int nI = 10;
283  int nJ = 10;
284  *plLine = new TPolyLine3D(5);
285 
286  {
287  TVector3 linevec;
288  int i, j;
289  i = -1 * nI;
290  j = -1 * nJ;
291  linevec = (fO + (mesh * i) * fU + (mesh * j) * fV);
292  (*plLine)->SetPoint(0, linevec.X(), linevec.Y(), linevec.Z());
293  i = -1 * nI;
294  j = 1 * nJ;
295  linevec = (fO + (mesh * i) * fU + (mesh * j) * fV);
296  (*plLine)->SetPoint(0, linevec.X(), linevec.Y(), linevec.Z());
297  i = 1 * nI;
298  j = -1 * nJ;
299  linevec = (fO + (mesh * i) * fU + (mesh * j) * fV);
300  (*plLine)->SetPoint(2, linevec.X(), linevec.Y(), linevec.Z());
301  i = 1 * nI;
302  j = 1 * nJ;
303  linevec = (fO + (mesh * i) * fU + (mesh * j) * fV);
304  (*plLine)->SetPoint(1, linevec.X(), linevec.Y(), linevec.Z());
305  i = -1 * nI;
306  j = -1 * nJ;
307  linevec = (fO + (mesh * i) * fU + (mesh * j) * fV);
308  (*plLine)->SetPoint(4, linevec.X(), linevec.Y(), linevec.Z());
309  }
310  for (int i = -1 * nI; i <= nI; ++i) {
311  for (int j = -1 * nJ; j <= nJ; ++j) {
312  TVector3 vec(fO + (mesh * i) * fU + (mesh * j) * fV);
313  int id = (i + 10) * 21 + j + 10;
314  (*pl)->SetPoint(id, vec.X(), vec.Y(), vec.Z());
315  }
316  }
317 
318  *u = new TPolyLine3D(2);
319  (*u)->SetPoint(0, fO.X(), fO.Y(), fO.Z());
320  (*u)->SetPoint(1, fO.X() + length * fU.X(), fO.Y() + length * fU.Y(), fO.Z() + length * fU.Z());
321  (*u)->SetLineWidth(2);
322  (*u)->SetLineColor(kGreen);
323 
324  *v = new TPolyLine3D(2);
325  (*v)->SetPoint(0, fO.X(), fO.Y(), fO.Z());
326  (*v)->SetPoint(1, fO.X() + length * fV.X(), fO.Y() + length * fV.Y(), fO.Z() + length * fV.Z());
327  (*v)->SetLineWidth(2);
328  (*v)->SetLineColor(kRed);
329 
330  if (n != NULL) {
331  *n = new TPolyLine3D(2);
332  TVector3 _n = getNormal();
333  (*n)->SetPoint(0, fO.X(), fO.Y(), fO.Z());
334  (*n)->SetPoint(1, fO.X() + length * _n.X(), fO.Y() + length * _n.Y(), fO.Z() + length * _n.Z());
335  (*n)->SetLineWidth(2);
336  (*n)->SetLineColor(kBlue);
337  }
338 }
TVector3 getNormal() const
Definition: GFDetPlane.cxx:133
Char_t n[5]
TVector3 genf::GFDetPlane::getNormal ( ) const
bool genf::GFDetPlane::inActive ( const TVector3 &  point,
const TVector3 &  dir 
) const
inline

intersect in the active area? C.f. GFAbsFinitePlane

Definition at line 136 of file GFDetPlane.h.

References straightLineToPlane().

Referenced by genf::RKTrackRep::Extrap(), and genf::RKTrackRep::RKutta().

137  {
138  return this->inActive(this->straightLineToPlane(point, dir));
139  }
TVector2 straightLineToPlane(const TVector3 &point, const TVector3 &dir) const
gives u,v coordinates of the intersection point of a straight line with plane
Definition: GFDetPlane.cxx:356
TDirectory * dir
Definition: macro.C:5
bool inActive(const TVector3 &point, const TVector3 &dir) const
intersect in the active area? C.f. GFAbsFinitePlane
Definition: GFDetPlane.h:136
bool genf::GFDetPlane::inActive ( double  u,
double  v 
) const
inline

inActive methods refer to finite plane. C.f. GFAbsFinitePlane

Definition at line 142 of file GFDetPlane.h.

References fFinitePlane, and genf::GFAbsFinitePlane::inActive().

143  {
144  if (fFinitePlane == NULL) return true;
145  return fFinitePlane->inActive(u, v);
146  }
virtual bool inActive(const double &u, const double &v) const =0
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
bool genf::GFDetPlane::inActive ( const TVector2 &  v) const
inline

inActive methods refer to finite plane. C.f. GFAbsFinitePlane

Definition at line 149 of file GFDetPlane.h.

References inActive().

Referenced by inActive().

149 { return inActive(v.X(), v.Y()); }
bool inActive(const TVector3 &point, const TVector3 &dir) const
intersect in the active area? C.f. GFAbsFinitePlane
Definition: GFDetPlane.h:136
TVector2 genf::GFDetPlane::LabToPlane ( const TVector3 &  x) const

transform from Lab system into plane

Definition at line 185 of file GFDetPlane.cxx.

References d, fO, and project().

Referenced by dist(), and setFinitePlane().

186 {
187  TVector3 d = x - fO;
188  return project(d);
189 }
Float_t x
Definition: compare.C:6
TVector2 project(const TVector3 &x) const
projecting a direction onto the plane:
Definition: GFDetPlane.cxx:178
Float_t d
Definition: plot.C:235
genf::GFDetPlane & genf::GFDetPlane::operator= ( const genf::GFDetPlane rhs)

Definition at line 72 of file GFDetPlane.cxx.

References genf::GFAbsFinitePlane::clone(), fFinitePlane, fO, fU, and fV.

73 {
74  if (this == &rhs) return *this;
75  if (fFinitePlane != NULL) { delete fFinitePlane; }
76  if (rhs.fFinitePlane != NULL) { fFinitePlane = rhs.fFinitePlane->clone(); }
77  else {
78  fFinitePlane = NULL;
79  }
80  fO = rhs.fO;
81  fU = rhs.fU;
82  fV = rhs.fV;
83  return *this;
84 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
virtual GFAbsFinitePlane * clone() const =0
Deep copy ctor for polymorphic class.
void genf::GFDetPlane::Print ( std::ostream &  out = std::cout) const

Definition at line 228 of file GFDetPlane.cxx.

References fFinitePlane, fO, fU, fV, getNormal(), and genf::GFAbsFinitePlane::Print().

Referenced by genf::GFBookkeeping::Print(), genf::GFAbsTrackRep::Print(), trkf::Track3DKalman::produce(), trkf::Track3DKalmanSPS::produce(), and setFinitePlane().

229 {
230  out << "GFDetPlane: "
231  << "O(" << fO.X() << "," << fO.Y() << "," << fO.Z() << ") "
232  << "u(" << fU.X() << "," << fU.Y() << "," << fU.Z() << ") "
233  << "v(" << fV.X() << "," << fV.Y() << "," << fV.Z() << ") "
234  << "n(" << getNormal().X() << "," << getNormal().Y() << "," << getNormal().Z() << ") "
235  << std::endl;
236  out << fFinitePlane << std::endl;
237  if (fFinitePlane != NULL) fFinitePlane->Print(out);
238 }
virtual void Print(std::ostream &out=std::cout) const =0
TVector3 getNormal() const
Definition: GFDetPlane.cxx:133
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
virtual void genf::GFDetPlane::Print ( Option_t *  ) const
inlineprivatevirtual

Definition at line 167 of file GFDetPlane.h.

References operator!=, and operator==.

168  {
169  throw std::logic_error(std::string(__func__) + "::Print(Option_t*) not available");
170  }
TVector2 genf::GFDetPlane::project ( const TVector3 &  x) const

projecting a direction onto the plane:

Definition at line 178 of file GFDetPlane.cxx.

References fU, fV, and x.

Referenced by LabToPlane(), setFinitePlane(), and straightLineToPlane().

179 {
180  Double_t xfU = fU * x;
181  Double_t xfV = fV * x;
182  return TVector2(xfU, xfV);
183 }
Float_t x
Definition: compare.C:6
void genf::GFDetPlane::sane ( )

Definition at line 206 of file GFDetPlane.cxx.

References E, fU, fV, getNormal(), n, and GFException::setFatal().

Referenced by GFDetPlane(), set(), setO(), setU(), setUV(), and setV().

207 {
208  if (fU == fV)
209  throw GFException("genf::GFDetPlane::sane() sanity check failed", __LINE__, __FILE__)
210  .setFatal();
211  // ensure unit vectors
212  fU.SetMag(1.);
213  fV.SetMag(1.);
214  // ensure orthogonal system
215  // fV should be reached by
216  // rotating fU around _n in a counterclockwise direction
217 
218  TVector3 n = getNormal();
219  fV = n.Cross(fU);
220 
221  TVector3 v = fU;
222  v.Rotate(TMath::Pi() * 0.5, n);
223  TVector3 null = v - fV;
224  if (null.Mag() >= 1E-6)
225  throw GFException("genf::GFDetPlane::sane(): non orthogonal!", __LINE__, __FILE__).setFatal();
226 }
TVector3 getNormal() const
Definition: GFDetPlane.cxx:133
Float_t E
Definition: plot.C:20
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:47
GFException & setFatal(bool b=true)
set fatal flag. if this is true, the fit stops for this current track repr.
Definition: GFException.h:75
Char_t n[5]
void genf::GFDetPlane::set ( const TVector3 &  o,
const TVector3 &  u,
const TVector3 &  v 
)

Definition at line 86 of file GFDetPlane.cxx.

References fO, fU, fV, and sane().

Referenced by genf::GFAbsTrackRep::reset().

87 {
88  fO = o;
89  fU = u;
90  fV = v;
91  sane();
92 }
void genf::GFDetPlane::setFinitePlane ( genf::GFAbsFinitePlane finite)
inline

Optionally, set the finite plane definition. This is most important for avoiding fake intersection points in fitting of loopers. This should be implemented for silicon detectors most importantly.

Definition at line 96 of file GFDetPlane.h.

References dir, dist(), distance(), fFinitePlane, getGraphics(), getNormal(), LabToPlane(), operator!=, operator==, Print(), project(), setNormal(), straightLineToPlane(), toLab(), and x.

96 { fFinitePlane = finite; }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:160
void genf::GFDetPlane::setNormal ( TVector3  n)

Definition at line 151 of file GFDetPlane.cxx.

References fU, and fV.

Referenced by genf::GFSpacepointHitPolicy::detPlane(), GFDetPlane(), genf::RKTrackRep::RKTrackRep(), setFinitePlane(), setNormal(), and setON().

152 {
153  n.SetMag(1.);
154  if (fabs(n.X()) > 0.1) {
155  fU.SetXYZ(1. / n.X() * (-1. * n.Y() - 1. * n.Z()), 1., 1.);
156  fU.SetMag(1.);
157  }
158  else {
159  if (fabs(n.Y()) > 0.1) {
160  fU.SetXYZ(1., 1. / n.Y() * (-1. * n.X() - 1. * n.Z()), 1.);
161  fU.SetMag(1.);
162  }
163  else {
164  fU.SetXYZ(1., 1., 1. / n.Z() * (-1. * n.X() - 1. * n.Y()));
165  fU.SetMag(1.);
166  }
167  }
168  fV = n.Cross(fU);
169 }
Char_t n[5]
void genf::GFDetPlane::setNormal ( double  X,
double  Y,
double  Z 
)

Definition at line 146 of file GFDetPlane.cxx.

References setNormal().

147 {
148  TVector3 N(X, Y, Z);
149  setNormal(N);
150 }
Float_t Y
Definition: plot.C:37
Float_t Z
Definition: plot.C:37
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:151
Float_t X
Definition: plot.C:37
void genf::GFDetPlane::setNormal ( const double &  theta,
const double &  phi 
)

Definition at line 171 of file GFDetPlane.cxx.

References n, and setNormal().

172 {
173  TVector3 n(
174  TMath::Sin(theta) * TMath::Cos(phi), TMath::Sin(theta) * TMath::Sin(phi), TMath::Cos(theta));
175  setNormal(n);
176 }
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:151
Char_t n[5]
void genf::GFDetPlane::setO ( const TVector3 &  o)

Definition at line 94 of file GFDetPlane.cxx.

References fO, and sane().

Referenced by genf::GFSpacepointHitPolicy::detPlane(), genf::RKTrackRep::extrapolateToLine(), getV(), and genf::RKTrackRep::RKTrackRep().

95 {
96  fO = o;
97  sane();
98 }
void genf::GFDetPlane::setO ( double  X,
double  Y,
double  Z 
)

Definition at line 99 of file GFDetPlane.cxx.

References fO, and sane().

100 {
101  fO.SetXYZ(X, Y, Z);
102  sane();
103 }
Float_t Y
Definition: plot.C:37
Float_t Z
Definition: plot.C:37
Float_t X
Definition: plot.C:37
void genf::GFDetPlane::setON ( const TVector3 &  o,
const TVector3 &  n 
)

Definition at line 140 of file GFDetPlane.cxx.

References fO, and setNormal().

Referenced by genf::RKTrackRep::extrapolateToPoint(), and getV().

141 {
142  fO = o;
143  setNormal(n);
144 }
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:151
Char_t n[5]
void genf::GFDetPlane::setU ( const TVector3 &  u)

Definition at line 105 of file GFDetPlane.cxx.

References fU, and sane().

Referenced by genf::RKTrackRep::extrapolateToLine(), and getV().

106 {
107  fU = u;
108  sane();
109 }
void genf::GFDetPlane::setU ( double  X,
double  Y,
double  Z 
)

Definition at line 110 of file GFDetPlane.cxx.

References fU, and sane().

111 {
112  fU.SetXYZ(X, Y, Z);
113  sane();
114 }
Float_t Y
Definition: plot.C:37
Float_t Z
Definition: plot.C:37
Float_t X
Definition: plot.C:37
void genf::GFDetPlane::setUV ( const TVector3 &  u,
const TVector3 &  v 
)

Definition at line 126 of file GFDetPlane.cxx.

References fU, fV, and sane().

Referenced by getV().

127 {
128  fU = u;
129  fV = v;
130  sane();
131 }
void genf::GFDetPlane::setV ( const TVector3 &  v)

Definition at line 116 of file GFDetPlane.cxx.

References fV, and sane().

Referenced by genf::RKTrackRep::extrapolateToLine(), and getV().

117 {
118  fV = v;
119  sane();
120 }
void genf::GFDetPlane::setV ( double  X,
double  Y,
double  Z 
)

Definition at line 121 of file GFDetPlane.cxx.

References fV, and sane().

122 {
123  fV.SetXYZ(X, Y, Z);
124  sane();
125 }
Float_t Y
Definition: plot.C:37
Float_t Z
Definition: plot.C:37
Float_t X
Definition: plot.C:37
TVector2 genf::GFDetPlane::straightLineToPlane ( const TVector3 &  point,
const TVector3 &  dir 
) const

gives u,v coordinates of the intersection point of a straight line with plane

Definition at line 356 of file GFDetPlane.cxx.

References E, fO, getNormal(), and project().

Referenced by inActive(), and setFinitePlane().

357 {
358  TVector3 dirNorm(dir);
359  dirNorm.SetMag(1.);
360  TVector3 normal = getNormal();
361  double dirTimesN = dirNorm * normal;
362  if (fabs(dirTimesN) < 1.E-6) { //straight line is parallel to plane, so return infinity
363  //doesnt parallel mean that they cross in infinity ;-)
364  return TVector2(1.E100, 1.E100);
365  }
366  double t = 1 / dirTimesN * ((fO - point) * normal);
367  return project(point - fO + t * dirNorm);
368 }
TVector2 project(const TVector3 &x) const
projecting a direction onto the plane:
Definition: GFDetPlane.cxx:178
TVector3 getNormal() const
Definition: GFDetPlane.cxx:133
Float_t E
Definition: plot.C:20
TDirectory * dir
Definition: macro.C:5
TVector3 genf::GFDetPlane::toLab ( const TVector2 &  x) const

transform from plane coordinates to lab system

Definition at line 191 of file GFDetPlane.cxx.

References d, fO, fU, and fV.

Referenced by dist(), and setFinitePlane().

192 {
193  TVector3 d(fO);
194  d += x.X() * fU;
195  d += x.Y() * fV;
196  return d;
197 }
Float_t x
Definition: compare.C:6
Float_t d
Definition: plot.C:235

Friends And Related Function Documentation

bool operator!= ( const GFDetPlane ,
const GFDetPlane  
)
friend

returns NOT ==

Referenced by Print(), and setFinitePlane().

bool operator== ( const GFDetPlane ,
const GFDetPlane  
)
friend

this operator is called very often in Kalman filtering. It checks equality of planes by comparing the 9 double values that define them.

Referenced by Print(), and setFinitePlane().

Member Data Documentation

genf::GFAbsFinitePlane* genf::GFDetPlane::fFinitePlane

Definition at line 160 of file GFDetPlane.h.

Referenced by GFDetPlane(), inActive(), operator=(), Print(), setFinitePlane(), and ~GFDetPlane().

TVector3 genf::GFDetPlane::fO
TVector3 genf::GFDetPlane::fU
TVector3 genf::GFDetPlane::fV

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