Public Member Functions | Public Attributes | Private Member Functions | Friends | List of all members
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 65 of file GFDetPlane.h.

Constructor & Destructor Documentation

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

Definition at line 44 of file GFDetPlane.cxx.

45  :fFinitePlane(finite)
46 {
47  static TRandom3 r(0);
48  fO.SetXYZ(0.,0.,0.);
49  fU.SetXYZ(r.Uniform(),r.Uniform(),0.);
50  fV.SetXYZ(r.Uniform(),r.Uniform(),0.);
51  sane();
52 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:162
genf::GFDetPlane::GFDetPlane ( const TVector3 &  o,
const TVector3 &  u,
const TVector3 &  v,
genf::GFAbsFinitePlane finite = NULL 
)

Definition at line 36 of file GFDetPlane.cxx.

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

Definition at line 54 of file GFDetPlane.cxx.

57  :fO(o),fFinitePlane(finite){
58  setNormal(n);
59 }
std::void_t< T > n
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:162
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:163
genf::GFDetPlane::~GFDetPlane ( )
virtual

Definition at line 61 of file GFDetPlane.cxx.

61  {
62  if(fFinitePlane!=NULL) delete fFinitePlane;
63 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:162
genf::GFDetPlane::GFDetPlane ( const genf::GFDetPlane rhs)

Definition at line 65 of file GFDetPlane.cxx.

65  : TObject(rhs) {
66  if(rhs.fFinitePlane != NULL) fFinitePlane = rhs.fFinitePlane->clone();
67  else fFinitePlane = NULL;
68  fO = rhs.fO;
69  fU = rhs.fU;
70  fV = rhs.fV;
71 }
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:162
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 214 of file GFDetPlane.cxx.

215 {
216  TVector2 p=LabToPlane(x);
217  TVector3 xplane=toLab(p);
218  return xplane-x;
219 }
TVector3 toLab(const TVector2 &x) const
transform from plane coordinates to lab system
Definition: GFDetPlane.cxx:203
TVector2 LabToPlane(const TVector3 &x) const
transform from Lab system into plane
Definition: GFDetPlane.cxx:196
p
Definition: test.py:223
list x
Definition: train.py:276
double genf::GFDetPlane::distance ( TVector3 &  v) const

Definition at line 354 of file GFDetPlane.cxx.

354  {
355  double s = (v - fO)*fU;
356  double t = (v - fO)*fV;
357  TVector3 distanceVector = v - fO - (s*fU) - (t*fV);
358  return distanceVector.Mag();
359 }
static QCString * s
Definition: config.cpp:1042
double genf::GFDetPlane::distance ( double  x,
double  y,
double  z 
) const

Definition at line 360 of file GFDetPlane.cxx.

360  {
361  TVector3 v(x,y,z);
362  double s = (v - fO)*fU;
363  double t = (v - fO)*fV;
364  TVector3 distanceVector = v - fO - (s*fU) - (t*fV);
365  return distanceVector.Mag();
366 }
list x
Definition: train.py:276
static QCString * s
Definition: config.cpp:1042
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 294 of file GFDetPlane.cxx.

294  {
295  *pl = new TPolyMarker3D(21*21,24);
296  (*pl)->SetMarkerSize(0.1);
297  (*pl)->SetMarkerColor(kBlue);
298  int nI=10;
299  int nJ=10;
300  *plLine = new TPolyLine3D(5);
301 
302  {
303  TVector3 linevec;
304  int i,j;
305  i=-1*nI;j=-1*nJ;
306  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
307  (*plLine)->SetPoint(0,linevec.X(),linevec.Y(),linevec.Z());
308  i=-1*nI;j=1*nJ;
309  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
310  (*plLine)->SetPoint(0,linevec.X(),linevec.Y(),linevec.Z());
311  i=1*nI;j=-1*nJ;
312  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
313  (*plLine)->SetPoint(2,linevec.X(),linevec.Y(),linevec.Z());
314  i=1*nI;j=1*nJ;
315  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
316  (*plLine)->SetPoint(1,linevec.X(),linevec.Y(),linevec.Z());
317  i=-1*nI;j=-1*nJ;
318  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
319  (*plLine)->SetPoint(4,linevec.X(),linevec.Y(),linevec.Z());
320 
321  }
322  for (int i=-1*nI;i<=nI;++i){
323  for (int j=-1*nJ;j<=nJ;++j){
324  TVector3 vec(fO+(mesh*i)*fU+(mesh*j)*fV);
325  int id=(i+10)*21+j+10;
326  (*pl)->SetPoint(id,vec.X(),vec.Y(),vec.Z());
327  }
328  }
329 
330 
331  *u = new TPolyLine3D(2);
332  (*u)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
333  (*u)->SetPoint(1,fO.X()+length*fU.X(),fO.Y()+length*fU.Y(),fO.Z()+length*fU.Z());
334  (*u)->SetLineWidth(2);
335  (*u)->SetLineColor(kGreen);
336 
337 
338  *v = new TPolyLine3D(2);
339  (*v)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
340  (*v)->SetPoint(1,fO.X()+length*fV.X(),fO.Y()+length*fV.Y(),fO.Z()+length*fV.Z());
341  (*v)->SetLineWidth(2);
342  (*v)->SetLineColor(kRed);
343 
344  if(n!=NULL){
345  *n = new TPolyLine3D(2);
346  TVector3 _n=getNormal();
347  (*n)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
348  (*n)->SetPoint(1,fO.X()+length*_n.X(),fO.Y()+length*_n.Y(),fO.Z()+length*_n.Z());
349  (*n)->SetLineWidth(2);
350  (*n)->SetLineColor(kBlue);
351  }
352 }
TVector3 getNormal() const
Definition: GFDetPlane.cxx:145
std::void_t< T > n
TVector3 genf::GFDetPlane::getNormal ( ) const

Definition at line 145 of file GFDetPlane.cxx.

146 {
147  TVector3 result=fU.Cross(fV);
148  result.SetMag(1.);
149  return result;
150 }
static QCString result
TVector3 genf::GFDetPlane::getO ( ) const
inline

Definition at line 81 of file GFDetPlane.h.

81 {return fO;}
TVector3 genf::GFDetPlane::getU ( ) const
inline

Definition at line 82 of file GFDetPlane.h.

82 {return fU;}
TVector3 genf::GFDetPlane::getV ( ) const
inline

Definition at line 83 of file GFDetPlane.h.

83 {return fV;}
bool genf::GFDetPlane::inActive ( const TVector3 &  point,
const TVector3 &  dir 
) const
inline

intersect in the active area? C.f. GFAbsFinitePlane

Definition at line 138 of file GFDetPlane.h.

138  {
139  return this->inActive( this->straightLineToPlane(point,dir));
140  }
string dir
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:368
bool inActive(const TVector3 &point, const TVector3 &dir) const
intersect in the active area? C.f. GFAbsFinitePlane
Definition: GFDetPlane.h:138
bool genf::GFDetPlane::inActive ( double  u,
double  v 
) const
inline

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

Definition at line 143 of file GFDetPlane.h.

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:162
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.

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

transform from Lab system into plane

Definition at line 196 of file GFDetPlane.cxx.

197 {
198  TVector3 d=x-fO;
199  return project(d);
200 }
TVector2 project(const TVector3 &x) const
projecting a direction onto the plane:
Definition: GFDetPlane.cxx:188
list x
Definition: train.py:276
genf::GFDetPlane & genf::GFDetPlane::operator= ( const genf::GFDetPlane rhs)

Definition at line 72 of file GFDetPlane.cxx.

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

Definition at line 247 of file GFDetPlane.cxx.

248 {
249  out<<"GFDetPlane: "
250  <<"O("<<fO.X()<<","<<fO.Y()<<","<<fO.Z()<<") "
251  <<"u("<<fU.X()<<","<<fU.Y()<<","<<fU.Z()<<") "
252  <<"v("<<fV.X()<<","<<fV.Y()<<","<<fV.Z()<<") "
253  <<"n("<<getNormal().X()<<","<<getNormal().Y()<<","<<getNormal().Z()<<") "
254  <<std::endl;
255  out << fFinitePlane << std::endl;
256  if(fFinitePlane!=NULL) fFinitePlane->Print(out);
257 
258 }
virtual void Print(std::ostream &out=std::cout) const =0
TVector3 getNormal() const
Definition: GFDetPlane.cxx:145
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:162
QTextStream & endl(QTextStream &s)
virtual void genf::GFDetPlane::Print ( Option_t *  ) const
inlineprivatevirtual

Definition at line 169 of file GFDetPlane.h.

170  { throw std::logic_error(std::string(__func__) + "::Print(Option_t*) not available"); }
std::string string
Definition: nybbler.cc:12
TVector2 genf::GFDetPlane::project ( const TVector3 &  x) const

projecting a direction onto the plane:

Definition at line 188 of file GFDetPlane.cxx.

189 {
190  Double_t xfU=fU*x;
191  Double_t xfV=fV*x;
192  return TVector2(xfU,xfV);
193 }
list x
Definition: train.py:276
void genf::GFDetPlane::sane ( )

Definition at line 224 of file GFDetPlane.cxx.

224  {
225  if (fU==fV)
226  throw GFException("genf::GFDetPlane::sane() sanity check failed", __LINE__, __FILE__).setFatal();
227  // ensure unit vectors
228  fU.SetMag(1.);
229  fV.SetMag(1.);
230  // ensure orthogonal system
231  // fV should be reached by
232  // rotating fU around _n in a counterclockwise direction
233 
234  TVector3 n=getNormal();
235  fV=n.Cross(fU);
236 
237  TVector3 v=fU;
238  v.Rotate(TMath::Pi()*0.5,n);
239  TVector3 null=v-fV;
240  if (null.Mag() >= 1E-6)
241  throw GFException("genf::GFDetPlane::sane(): non orthogonal!", __LINE__, __FILE__).setFatal();
242 
243 }
TVector3 getNormal() const
Definition: GFDetPlane.cxx:145
std::void_t< T > n
null
Definition: demo.py:50
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:48
GFException & setFatal(bool b=true)
set fatal flag. if this is true, the fit stops for this current track repr.
Definition: GFException.h:78
void genf::GFDetPlane::set ( const TVector3 &  o,
const TVector3 &  u,
const TVector3 &  v 
)

Definition at line 90 of file GFDetPlane.cxx.

93 {
94  fO=o;fU=u;fV=v;
95  sane();
96 }
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 102 of file GFDetPlane.h.

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

Definition at line 163 of file GFDetPlane.cxx.

163  {
164  n.SetMag(1.);
165  if( fabs(n.X()) > 0.1 ){
166  fU.SetXYZ(1./n.X()*(-1.*n.Y()-1.*n.Z()),1.,1.);
167  fU.SetMag(1.);
168  }
169  else {
170  if(fabs(n.Y()) > 0.1){
171  fU.SetXYZ(1.,1./n.Y()*(-1.*n.X()-1.*n.Z()),1.);
172  fU.SetMag(1.);
173  }
174  else{
175  fU.SetXYZ(1.,1.,1./n.Z()*(-1.*n.X()-1.*n.Y()));
176  fU.SetMag(1.);
177  }
178  }
179  fV=n.Cross(fU);
180 }
std::void_t< T > n
void genf::GFDetPlane::setNormal ( double  X,
double  Y,
double  Z 
)

Definition at line 158 of file GFDetPlane.cxx.

158  {
159  TVector3 N(X,Y,Z);
160  setNormal(N);
161 }
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:163
void genf::GFDetPlane::setNormal ( const double &  theta,
const double &  phi 
)

Definition at line 182 of file GFDetPlane.cxx.

182  {
183  TVector3 n(TMath::Sin(theta)*TMath::Cos(phi),TMath::Sin(theta)*TMath::Sin(phi),TMath::Cos(theta));
184  setNormal(n);
185 }
std::void_t< T > n
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:163
void genf::GFDetPlane::setO ( const TVector3 &  o)

Definition at line 99 of file GFDetPlane.cxx.

100 {
101  fO=o;
102  sane();
103 }
void genf::GFDetPlane::setO ( double  X,
double  Y,
double  Z 
)

Definition at line 105 of file GFDetPlane.cxx.

106 {
107  fO.SetXYZ(X,Y,Z);
108  sane();
109 }
void genf::GFDetPlane::setON ( const TVector3 &  o,
const TVector3 &  n 
)

Definition at line 152 of file GFDetPlane.cxx.

152  {
153  fO = o;
154  setNormal(n);
155 }
std::void_t< T > n
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:163
void genf::GFDetPlane::setU ( const TVector3 &  u)

Definition at line 112 of file GFDetPlane.cxx.

113 {
114  fU=u;
115  sane();
116 }
void genf::GFDetPlane::setU ( double  X,
double  Y,
double  Z 
)

Definition at line 118 of file GFDetPlane.cxx.

119 {
120  fU.SetXYZ(X,Y,Z);
121  sane();
122 }
void genf::GFDetPlane::setUV ( const TVector3 &  u,
const TVector3 &  v 
)

Definition at line 137 of file GFDetPlane.cxx.

138 {
139  fU=u;
140  fV=v;
141  sane();
142 }
void genf::GFDetPlane::setV ( const TVector3 &  v)

Definition at line 125 of file GFDetPlane.cxx.

126 {
127  fV=v;
128  sane();
129 }
void genf::GFDetPlane::setV ( double  X,
double  Y,
double  Z 
)

Definition at line 131 of file GFDetPlane.cxx.

132 {
133  fV.SetXYZ(X,Y,Z);
134  sane();
135 }
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 368 of file GFDetPlane.cxx.

368  {
369  TVector3 dirNorm(dir);
370  dirNorm.SetMag(1.);
371  TVector3 normal = getNormal();
372  double dirTimesN = dirNorm*normal;
373  if(fabs(dirTimesN)<1.E-6){//straight line is parallel to plane, so return infinity
374  //doesnt parallel mean that they cross in infinity ;-)
375  return TVector2(1.E100,1.E100);
376  }
377  double t = 1/dirTimesN * ((fO-point)*normal);
378  return project(point - fO + t * dirNorm);
379 }
string dir
TVector2 project(const TVector3 &x) const
projecting a direction onto the plane:
Definition: GFDetPlane.cxx:188
TVector3 getNormal() const
Definition: GFDetPlane.cxx:145
TVector3 genf::GFDetPlane::toLab ( const TVector2 &  x) const

transform from plane coordinates to lab system

Definition at line 203 of file GFDetPlane.cxx.

204 {
205  TVector3 d(fO);
206  d+=x.X()*fU;
207  d+=x.Y()*fV;
208  return d;
209 }
list x
Definition: train.py:276

Friends And Related Function Documentation

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

returns NOT ==

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.

Member Data Documentation

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

Definition at line 162 of file GFDetPlane.h.

TVector3 genf::GFDetPlane::fO

Definition at line 157 of file GFDetPlane.h.

TVector3 genf::GFDetPlane::fU

Definition at line 159 of file GFDetPlane.h.

TVector3 genf::GFDetPlane::fV

Definition at line 160 of file GFDetPlane.h.


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