LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
GFDetPlane.cxx
Go to the documentation of this file.
1 /* Copyright 2008-2010, Technische Universitaet Muenchen,
2  Authors: Christian Hoeppner & Sebastian Neubert
3 
4  This file is part of GENFIT.
5 
6  GENFIT is free software: you can redistribute it and/or modify
7  it under the terms of the GNU Lesser General Public License as published
8  by the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  GENFIT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU Lesser General Public License for more details.
15 
16  You should have received a copy of the GNU Lesser General Public License
17  along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
18 */
21 
22 #include <iostream>
23 #include <cmath>
24 #include "TMath.h"
25 #include "TRandom3.h"
26 
27 //ClassImp(GFDetPlane)
28 
29 
30 
31 genf::GFDetPlane::GFDetPlane(const TVector3& o,
32  const TVector3& u,
33  const TVector3& v,
34  GFAbsFinitePlane* finite)
35 :fO(o),fU(u),fV(v),fFinitePlane(finite)
36 {
37  sane();
38 }
40  :fFinitePlane(finite)
41 {
42  static TRandom3 r(0);
43  fO.SetXYZ(0.,0.,0.);
44  fU.SetXYZ(r.Uniform(),r.Uniform(),0.);
45  fV.SetXYZ(r.Uniform(),r.Uniform(),0.);
46  sane();
47 }
48 
49 genf::GFDetPlane::GFDetPlane(const TVector3& o,
50  const TVector3& n,
51  GFAbsFinitePlane* finite)
52  :fO(o),fFinitePlane(finite){
53  setNormal(n);
54 }
55 
57  if(fFinitePlane!=NULL) delete fFinitePlane;
58 }
59 
61  if(rhs.fFinitePlane != NULL) fFinitePlane = rhs.fFinitePlane->clone();
62  else fFinitePlane = NULL;
63  fO = rhs.fO;
64  fU = rhs.fU;
65  fV = rhs.fV;
66 }
68  if (this == &rhs) return *this;
69  if(fFinitePlane!=NULL) {
70  delete fFinitePlane;
71  }
72  if(rhs.fFinitePlane != NULL){
74  }
75  else{
76  fFinitePlane = NULL;
77  }
78  fO = rhs.fO;
79  fU = rhs.fU;
80  fV = rhs.fV;
81  return *this;
82 }
83 
84 void
85 genf::GFDetPlane::set(const TVector3& o,
86  const TVector3& u,
87  const TVector3& v)
88 {
89  fO=o;fU=u;fV=v;
90  sane();
91 }
92 
93 void
94 genf::GFDetPlane::setO(const TVector3& o)
95 {
96  fO=o;
97  sane();
98 }
99 void
100 genf::GFDetPlane::setO(double X,double Y,double Z)
101 {
102  fO.SetXYZ(X,Y,Z);
103  sane();
104 }
105 
106 void
107 genf::GFDetPlane::setU(const TVector3& u)
108 {
109  fU=u;
110  sane();
111 }
112 void
113 genf::GFDetPlane::setU(double X,double Y,double Z)
114 {
115  fU.SetXYZ(X,Y,Z);
116  sane();
117 }
118 
119 void
120 genf::GFDetPlane::setV(const TVector3& v)
121 {
122  fV=v;
123  sane();
124 }
125 void
126 genf::GFDetPlane::setV(double X,double Y,double Z)
127 {
128  fV.SetXYZ(X,Y,Z);
129  sane();
130 }
131 void
132 genf::GFDetPlane::setUV(const TVector3& u,const TVector3& v)
133 {
134  fU=u;
135  fV=v;
136  sane();
137 }
138 
139 TVector3
141 {
142  TVector3 result=fU.Cross(fV);
143  result.SetMag(1.);
144  return result;
145 }
146 
147 void genf::GFDetPlane::setON(const TVector3& o,const TVector3& n){
148  fO = o;
149  setNormal(n);
150 }
151 
152 void
153 genf::GFDetPlane::setNormal(double X,double Y,double Z){
154  TVector3 N(X,Y,Z);
155  setNormal(N);
156 }
157 void
159  n.SetMag(1.);
160  if( fabs(n.X()) > 0.1 ){
161  fU.SetXYZ(1./n.X()*(-1.*n.Y()-1.*n.Z()),1.,1.);
162  fU.SetMag(1.);
163  }
164  else {
165  if(fabs(n.Y()) > 0.1){
166  fU.SetXYZ(1.,1./n.Y()*(-1.*n.X()-1.*n.Z()),1.);
167  fU.SetMag(1.);
168  }
169  else{
170  fU.SetXYZ(1.,1.,1./n.Z()*(-1.*n.X()-1.*n.Y()));
171  fU.SetMag(1.);
172  }
173  }
174  fV=n.Cross(fU);
175 }
176 
177 void genf::GFDetPlane::setNormal(const double& theta, const double& phi){
178  TVector3 n(TMath::Sin(theta)*TMath::Cos(phi),TMath::Sin(theta)*TMath::Sin(phi),TMath::Cos(theta));
179  setNormal(n);
180 }
181 
182 TVector2
183 genf::GFDetPlane::project(const TVector3& x)const
184 {
185  Double_t xfU=fU*x;
186  Double_t xfV=fV*x;
187  return TVector2(xfU,xfV);
188 }
189 
190 TVector2
191 genf::GFDetPlane::LabToPlane(const TVector3& x)const
192 {
193  TVector3 d=x-fO;
194  return project(d);
195 }
196 
197 TVector3
198 genf::GFDetPlane::toLab(const TVector2& x)const
199 {
200  TVector3 d(fO);
201  d+=x.X()*fU;
202  d+=x.Y()*fV;
203  return d;
204 }
205 
206 
207 
208 TVector3
209 genf::GFDetPlane::dist(const TVector3& x)const
210 {
211  TVector2 p=LabToPlane(x);
212  TVector3 xplane=toLab(p);
213  return xplane-x;
214 }
215 
216 
217 
218 void
220  if (fU==fV)
221  throw GFException("genf::GFDetPlane::sane() sanity check failed", __LINE__, __FILE__).setFatal();
222  // ensure unit vectors
223  fU.SetMag(1.);
224  fV.SetMag(1.);
225  // ensure orthogonal system
226  // fV should be reached by
227  // rotating fU around _n in a counterclockwise direction
228 
229  TVector3 n=getNormal();
230  fV=n.Cross(fU);
231 
232  TVector3 v=fU;
233  v.Rotate(TMath::Pi()*0.5,n);
234  TVector3 null=v-fV;
235  if (null.Mag() >= 1E-6)
236  throw GFException("genf::GFDetPlane::sane(): non orthogonal!", __LINE__, __FILE__).setFatal();
237 
238 }
239 
240 
241 void
242 genf::GFDetPlane::Print(std::ostream& out /* = std::cout */) const
243 {
244  out<<"GFDetPlane: "
245  <<"O("<<fO.X()<<","<<fO.Y()<<","<<fO.Z()<<") "
246  <<"u("<<fU.X()<<","<<fU.Y()<<","<<fU.Z()<<") "
247  <<"v("<<fV.X()<<","<<fV.Y()<<","<<fV.Z()<<") "
248  <<"n("<<getNormal().X()<<","<<getNormal().Y()<<","<<getNormal().Z()<<") "
249  <<std::endl;
250  out << fFinitePlane << std::endl;
251  if(fFinitePlane!=NULL) fFinitePlane->Print(out);
252 
253 }
254 
255 
256 /*
257  I could write pages of comments about correct equality checking for
258  floating point numbers, but: When two planes are as close as 10E-5 cm
259  in all nine numbers that define the plane, this will be enough for all
260  practical purposes
261  */
262 
263 // == and != are friends, not members, hence not "genf::GFDetPlane::operater==" below.
264 #define DETPLANE_EPSILON 1.E-5
265 bool genf::operator== (const GFDetPlane& lhs, const GFDetPlane& rhs){
266  if(
267  fabs( (lhs.fO.X()-rhs.fO.X()) ) > DETPLANE_EPSILON ||
268  fabs( (lhs.fO.Y()-rhs.fO.Y()) ) > DETPLANE_EPSILON ||
269  fabs( (lhs.fO.Z()-rhs.fO.Z()) ) > DETPLANE_EPSILON
270  ) return false;
271  else if(
272  fabs( (lhs.fU.X()-rhs.fU.X()) ) > DETPLANE_EPSILON ||
273  fabs( (lhs.fU.Y()-rhs.fU.Y()) ) > DETPLANE_EPSILON ||
274  fabs( (lhs.fU.Z()-rhs.fU.Z()) ) > DETPLANE_EPSILON
275  ) return false;
276  else if(
277  fabs( (lhs.fV.X()-rhs.fV.X()) ) > DETPLANE_EPSILON ||
278  fabs( (lhs.fV.Y()-rhs.fV.Y()) ) > DETPLANE_EPSILON ||
279  fabs( (lhs.fV.Z()-rhs.fV.Z()) ) > DETPLANE_EPSILON
280  ) return false;
281  return true;
282 }
283 
284 bool genf::operator!= (const GFDetPlane& lhs, const GFDetPlane& rhs){
285  return !(lhs == rhs);
286 }
287 
288 
289 void genf::GFDetPlane::getGraphics(double mesh, double length, TPolyMarker3D **pl,TPolyLine3D** plLine, TPolyLine3D **u, TPolyLine3D **v, TPolyLine3D**n){
290  *pl = new TPolyMarker3D(21*21,24);
291  (*pl)->SetMarkerSize(0.1);
292  (*pl)->SetMarkerColor(kBlue);
293  int nI=10;
294  int nJ=10;
295  *plLine = new TPolyLine3D(5);
296 
297  {
298  TVector3 linevec;
299  int i,j;
300  i=-1*nI;j=-1*nJ;
301  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
302  (*plLine)->SetPoint(0,linevec.X(),linevec.Y(),linevec.Z());
303  i=-1*nI;j=1*nJ;
304  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
305  (*plLine)->SetPoint(0,linevec.X(),linevec.Y(),linevec.Z());
306  i=1*nI;j=-1*nJ;
307  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
308  (*plLine)->SetPoint(2,linevec.X(),linevec.Y(),linevec.Z());
309  i=1*nI;j=1*nJ;
310  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
311  (*plLine)->SetPoint(1,linevec.X(),linevec.Y(),linevec.Z());
312  i=-1*nI;j=-1*nJ;
313  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
314  (*plLine)->SetPoint(4,linevec.X(),linevec.Y(),linevec.Z());
315 
316  }
317  for (int i=-1*nI;i<=nI;++i){
318  for (int j=-1*nJ;j<=nJ;++j){
319  TVector3 vec(fO+(mesh*i)*fU+(mesh*j)*fV);
320  int id=(i+10)*21+j+10;
321  (*pl)->SetPoint(id,vec.X(),vec.Y(),vec.Z());
322  }
323  }
324 
325 
326  *u = new TPolyLine3D(2);
327  (*u)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
328  (*u)->SetPoint(1,fO.X()+length*fU.X(),fO.Y()+length*fU.Y(),fO.Z()+length*fU.Z());
329  (*u)->SetLineWidth(2);
330  (*u)->SetLineColor(kGreen);
331 
332 
333  *v = new TPolyLine3D(2);
334  (*v)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
335  (*v)->SetPoint(1,fO.X()+length*fV.X(),fO.Y()+length*fV.Y(),fO.Z()+length*fV.Z());
336  (*v)->SetLineWidth(2);
337  (*v)->SetLineColor(kRed);
338 
339  if(n!=NULL){
340  *n = new TPolyLine3D(2);
341  TVector3 _n=getNormal();
342  (*n)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
343  (*n)->SetPoint(1,fO.X()+length*_n.X(),fO.Y()+length*_n.Y(),fO.Z()+length*_n.Z());
344  (*n)->SetLineWidth(2);
345  (*n)->SetLineColor(kBlue);
346  }
347 }
348 
349 double genf::GFDetPlane::distance(TVector3& v) const{
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 }
355 double genf::GFDetPlane::distance(double x,double y,double z) const{
356  TVector3 v(x,y,z);
357  double s = (v - fO)*fU;
358  double t = (v - fO)*fV;
359  TVector3 distanceVector = v - fO - (s*fU) - (t*fV);
360  return distanceVector.Mag();
361 }
362 
363 TVector2 genf::GFDetPlane::straightLineToPlane (const TVector3& point,const TVector3& dir) const{
364  TVector3 dirNorm(dir);
365  dirNorm.SetMag(1.);
366  TVector3 normal = getNormal();
367  double dirTimesN = dirNorm*normal;
368  if(fabs(dirTimesN)<1.E-6){//straight line is parallel to plane, so return infinity
369  //doesnt parallel mean that they cross in infinity ;-)
370  return TVector2(1.E100,1.E100);
371  }
372  double t = 1/dirTimesN * ((fO-point)*normal);
373  return project(point - fO + t * dirNorm);
374 }
375 
Float_t x
Definition: compare.C:6
Float_t s
Definition: plot.C:23
void set(const TVector3 &o, const TVector3 &u, const TVector3 &v)
Definition: GFDetPlane.cxx:85
void setON(const TVector3 &o, const TVector3 &n)
Definition: GFDetPlane.cxx:147
void setU(const TVector3 &u)
Definition: GFDetPlane.cxx:107
GFDetPlane & operator=(const genf::GFDetPlane &)
Definition: GFDetPlane.cxx:67
bool operator!=(const genf::GFDetPlane &, const genf::GFDetPlane &)
Definition: GFDetPlane.cxx:284
Float_t E
Definition: plot.C:23
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.
Definition: GFDetPlane.cxx:289
TVector3 dist(const TVector3 &point) const
Definition: GFDetPlane.cxx:209
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:279
Float_t Y
Definition: plot.C:39
void setO(const TVector3 &o)
Definition: GFDetPlane.cxx:94
void setUV(const TVector3 &u, const TVector3 &v)
Definition: GFDetPlane.cxx:132
void Print(std::ostream &out=std::cout) const
Definition: GFDetPlane.cxx:242
#define DETPLANE_EPSILON
Definition: GFDetPlane.cxx:264
TVector2 project(const TVector3 &x) const
projecting a direction onto the plane:
Definition: GFDetPlane.cxx:183
Float_t Z
Definition: plot.C:39
virtual void Print(std::ostream &out=std::cout) const =0
TVector3 getNormal() const
Definition: GFDetPlane.cxx:140
void setV(const TVector3 &v)
Definition: GFDetPlane.cxx:120
double distance(TVector3 &) const
Definition: GFDetPlane.cxx:349
bool operator==(const genf::GFDetPlane &, const genf::GFDetPlane &)
Definition: GFDetPlane.cxx:265
genf::GFAbsFinitePlane * fFinitePlane
Definition: GFDetPlane.h:157
TVector3 toLab(const TVector2 &x) const
transform from plane coordinates to lab system
Definition: GFDetPlane.cxx:198
TVector2 LabToPlane(const TVector3 &x) const
transform from Lab system into plane
Definition: GFDetPlane.cxx:191
Float_t d
Definition: plot.C:237
virtual ~GFDetPlane()
Definition: GFDetPlane.cxx:56
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:363
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:50
TDirectory * dir
Definition: macro.C:5
void setNormal(TVector3 n)
Definition: GFDetPlane.cxx:158
GFException & setFatal(bool b=true)
set fatal flag. if this is true, the fit stops for this current track repr.
Definition: GFException.h:80
Char_t n[5]
Float_t X
Definition: plot.C:39
virtual GFAbsFinitePlane * clone() const =0
Deep copy ctor for polymorphic class.
GFDetPlane(genf::GFAbsFinitePlane *finite=NULL)
Definition: GFDetPlane.cxx:39