EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GFDetPlane.cxx
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file GFDetPlane.cxx
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 */
19 #include "GFDetPlane.h"
20 
21 #include "assert.h"
22 #include <iostream>
23 #include <cmath>
24 #include "TMath.h"
25 #include "TRandom3.h"
26 
28 
29 GFDetPlane::GFDetPlane(const TVector3& o,
30  const TVector3& u,
31  const TVector3& v,
32  GFAbsFinitePlane* finite)
33 :fO(o),fU(u),fV(v),fFinitePlane(finite)
34 {
35  sane();
36 }
38  :fFinitePlane(finite)
39 {
40  static TRandom3 r(0);
41  fO.SetXYZ(0.,0.,0.);
42  fU.SetXYZ(r.Uniform(),r.Uniform(),0.);
43  fV.SetXYZ(r.Uniform(),r.Uniform(),0.);
44  sane();
45 }
46 
47 GFDetPlane::GFDetPlane(const TVector3& o,
48  const TVector3& n,
49  GFAbsFinitePlane* finite)
50  :fO(o),fFinitePlane(finite){
51  setNormal(n);
52 }
53 
55  if(fFinitePlane!=NULL) delete fFinitePlane;
56 }
57 
59  if(rhs.fFinitePlane != NULL) fFinitePlane = rhs.fFinitePlane->clone();
60  else fFinitePlane = NULL;
61  fO = rhs.fO;
62  fU = rhs.fU;
63  fV = rhs.fV;
64 }
66  assert(this!=&rhs);
67  if(fFinitePlane!=NULL) {
68  delete fFinitePlane;
69  }
70  if(rhs.fFinitePlane != NULL){
72  }
73  else{
74  fFinitePlane = NULL;
75  }
76  fO = rhs.fO;
77  fU = rhs.fU;
78  fV = rhs.fV;
79  return *this;
80 }
81 
82 void
83 GFDetPlane::set(const TVector3& o,
84  const TVector3& u,
85  const TVector3& v)
86 {
87  fO=o;fU=u;fV=v;
88  sane();
89 }
90 
91 void
92 GFDetPlane::setO(const TVector3& o)
93 {
94  fO=o;
95  sane();
96 }
97 void
98 GFDetPlane::setO(double X,double Y,double Z)
99 {
100  fO.SetXYZ(X,Y,Z);
101  sane();
102 }
103 
104 void
105 GFDetPlane::setU(const TVector3& u)
106 {
107  fU=u;
108  sane();
109 }
110 void
111 GFDetPlane::setU(double X,double Y,double Z)
112 {
113  fU.SetXYZ(X,Y,Z);
114  sane();
115 }
116 
117 void
118 GFDetPlane::setV(const TVector3& v)
119 {
120  fV=v;
121  sane();
122 }
123 void
124 GFDetPlane::setV(double X,double Y,double Z)
125 {
126  fV.SetXYZ(X,Y,Z);
127  sane();
128 }
129 void
130 GFDetPlane::setUV(const TVector3& u,const TVector3& v)
131 {
132  fU=u;
133  fV=v;
134  sane();
135 }
136 
137 TVector3
139 {
140  TVector3 result=fU.Cross(fV);
141  result.SetMag(1.);
142  return result;
143 }
144 
145 void GFDetPlane::setON(const TVector3& o,const TVector3& n){
146  fO = o;
147  setNormal(n);
148 }
149 
150 void
151 GFDetPlane::setNormal(double X,double Y,double Z){
152  TVector3 N(X,Y,Z);
153  setNormal(N);
154 }
155 void
157  n.SetMag(1.);
158  if( fabs(n.X()) > 0.1 ){
159  fU.SetXYZ(1./n.X()*(-1.*n.Y()-1.*n.Z()),1.,1.);
160  fU.SetMag(1.);
161  }
162  else {
163  if(fabs(n.Y()) > 0.1){
164  fU.SetXYZ(1.,1./n.Y()*(-1.*n.X()-1.*n.Z()),1.);
165  fU.SetMag(1.);
166  }
167  else{
168  fU.SetXYZ(1.,1.,1./n.Z()*(-1.*n.X()-1.*n.Y()));
169  fU.SetMag(1.);
170  }
171  }
172  fV=n.Cross(fU);
173 }
174 
175 void GFDetPlane::setNormal(const double& theta, const double& phi){
176  TVector3 n(TMath::Sin(theta)*TMath::Cos(phi),TMath::Sin(theta)*TMath::Sin(phi),TMath::Cos(theta));
177  setNormal(n);
178 }
179 
180 TVector2
181 GFDetPlane::project(const TVector3& x)const
182 {
183  Double_t xfU=fU*x;
184  Double_t xfV=fV*x;
185  return TVector2(xfU,xfV);
186 }
187 
188 TVector2
189 GFDetPlane::LabToPlane(const TVector3& x)const
190 {
191  TVector3 d=x-fO;
192  return project(d);
193 }
194 
195 TVector3
196 GFDetPlane::toLab(const TVector2& x)const
197 {
198  TVector3 d(fO);
199  d+=x.X()*fU;
200  d+=x.Y()*fV;
201  return d;
202 }
203 
204 
205 
206 TVector3
207 GFDetPlane::dist(const TVector3& x)const
208 {
209  TVector2 p=LabToPlane(x);
210  TVector3 xplane=toLab(p);
211  return xplane-x;
212 }
213 
214 
215 
216 void
218  assert(fU!=fV);
219  // ensure unit vectors
220  fU.SetMag(1.);
221  fV.SetMag(1.);
222  // ensure orthogonal system
223  // fV should be reached by
224  // rotating fU around _n in a counterclockwise direction
225 
226  TVector3 n=getNormal();
227  fV=n.Cross(fU);
228 
229  TVector3 v=fU;
230  v.Rotate(TMath::Pi()*0.5,n);
231  TVector3 null=v-fV;
232  assert(null.Mag()<1E-6);
233 
234 }
235 
236 
237 void
238 GFDetPlane::Print(const Option_t* option) const
239 {
240  std::cout<<"GFDetPlane: "
241  <<"O("<<fO.X()<<","<<fO.Y()<<","<<fO.Z()<<") "
242  <<"u("<<fU.X()<<","<<fU.Y()<<","<<fU.Z()<<") "
243  <<"v("<<fV.X()<<","<<fV.Y()<<","<<fV.Z()<<") "
244  <<"n("<<getNormal().X()<<","<<getNormal().Y()<<","<<getNormal().Z()<<") "
245  <<std::endl;
246  std::cout << fFinitePlane << std::endl;
247  if(fFinitePlane!=NULL) fFinitePlane->Print(option);
248 
249 }
250 
251 
252 /*
253  I could write pages of comments about correct equality checking for
254  floating point numbers, but: When two planes are as close as 10E-5 cm
255  in all nine numbers that define the plane, this will be enough for all
256  practical purposes
257  */
258 #define DETPLANE_EPSILON 1.E-5
259 bool operator== (const GFDetPlane& lhs, const GFDetPlane& rhs){
260  if(
261  fabs( (lhs.fO.X()-rhs.fO.X()) ) > DETPLANE_EPSILON ||
262  fabs( (lhs.fO.Y()-rhs.fO.Y()) ) > DETPLANE_EPSILON ||
263  fabs( (lhs.fO.Z()-rhs.fO.Z()) ) > DETPLANE_EPSILON
264  ) return false;
265  else if(
266  fabs( (lhs.fU.X()-rhs.fU.X()) ) > DETPLANE_EPSILON ||
267  fabs( (lhs.fU.Y()-rhs.fU.Y()) ) > DETPLANE_EPSILON ||
268  fabs( (lhs.fU.Z()-rhs.fU.Z()) ) > DETPLANE_EPSILON
269  ) return false;
270  else if(
271  fabs( (lhs.fV.X()-rhs.fV.X()) ) > DETPLANE_EPSILON ||
272  fabs( (lhs.fV.Y()-rhs.fV.Y()) ) > DETPLANE_EPSILON ||
273  fabs( (lhs.fV.Z()-rhs.fV.Z()) ) > DETPLANE_EPSILON
274  ) return false;
275  return true;
276 }
277 
278 bool operator!= (const GFDetPlane& lhs, const GFDetPlane& rhs){
279  return !(lhs==rhs);
280 }
281 
282 
283 void GFDetPlane::getGraphics(double mesh, double length, TPolyMarker3D **pl,TPolyLine3D** plLine, TPolyLine3D **u, TPolyLine3D **v, TPolyLine3D**n){
284  *pl = new TPolyMarker3D(21*21,24);
285  (*pl)->SetMarkerSize(0.1);
286  (*pl)->SetMarkerColor(kBlue);
287  int nI=10;
288  int nJ=10;
289  *plLine = new TPolyLine3D(5);
290 
291  {
292  TVector3 linevec;
293  int i,j;
294  i=-1*nI;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;j=1*nJ;
298  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
299  (*plLine)->SetPoint(0,linevec.X(),linevec.Y(),linevec.Z());
300  i=1*nI;j=-1*nJ;
301  linevec=(fO+(mesh*i)*fU+(mesh*j)*fV);
302  (*plLine)->SetPoint(2,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(1,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(4,linevec.X(),linevec.Y(),linevec.Z());
309 
310  }
311  for (int i=-1*nI;i<=nI;++i){
312  for (int j=-1*nJ;j<=nJ;++j){
313  TVector3 vec(fO+(mesh*i)*fU+(mesh*j)*fV);
314  int id=(i+10)*21+j+10;
315  (*pl)->SetPoint(id,vec.X(),vec.Y(),vec.Z());
316  }
317  }
318 
319 
320  *u = new TPolyLine3D(2);
321  (*u)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
322  (*u)->SetPoint(1,fO.X()+length*fU.X(),fO.Y()+length*fU.Y(),fO.Z()+length*fU.Z());
323  (*u)->SetLineWidth(2);
324  (*u)->SetLineColor(kGreen);
325 
326 
327  *v = new TPolyLine3D(2);
328  (*v)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
329  (*v)->SetPoint(1,fO.X()+length*fV.X(),fO.Y()+length*fV.Y(),fO.Z()+length*fV.Z());
330  (*v)->SetLineWidth(2);
331  (*v)->SetLineColor(kRed);
332 
333  if(n!=NULL){
334  *n = new TPolyLine3D(2);
335  TVector3 _n=getNormal();
336  (*n)->SetPoint(0,fO.X(),fO.Y(),fO.Z());
337  (*n)->SetPoint(1,fO.X()+length*_n.X(),fO.Y()+length*_n.Y(),fO.Z()+length*_n.Z());
338  (*n)->SetLineWidth(2);
339  (*n)->SetLineColor(kBlue);
340  }
341 }
342 
343 double GFDetPlane::distance(TVector3 v) const{
344  v -= fO;
345  double s = v*fU;
346  double t = v*fV;
347  TVector3 distanceVector = v - (s*fU) - (t*fV);
348  return distanceVector.Mag();
349 }
350 double GFDetPlane::distance(double x, double y, double z) const{
351  TVector3 v(x,y,z);
352  v -= fO;
353  double s = v*fU;
354  double t = v*fV;
355  TVector3 distanceVector = v - (s*fU) - (t*fV);
356  return distanceVector.Mag();
357 }
358 
359 TVector2 GFDetPlane::straightLineToPlane (const TVector3& point,const TVector3& dir) const{
360  TVector3 dirNorm(dir);
361  dirNorm.SetMag(1.);
362  TVector3 normal = getNormal();
363  double dirTimesN = dirNorm*normal;
364  if(fabs(dirTimesN)<1.E-6){//straight line is parallel to plane, so return infinity
365  //doesnt parallel mean that they cross in infinity ;-)
366  return TVector2(1.E100,1.E100);
367  }
368  double t = 1/dirTimesN * ((fO-point)*normal);
369  return project(point - fO + t * dirNorm);
370 }