CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TEveEllipsoidGL.cc
Go to the documentation of this file.
2 
3 #define protected public
4 #include "TEveProjections.h" // AMT missing getter for projection center / beam-spot
5 #undef protected
6 
9 #include "TEveProjectionManager.h"
10 
11 #include "TMath.h"
12 
13 #include "TGLRnrCtx.h"
14 
15 #include "TMatrixDEigen.h"
16 #include "TMatrixDSym.h"
17 
18 #include "TDecompSVD.h"
19 #include "TGLIncludes.h"
20 
21 //==============================================================================
22 // TEveEllipsoidGL
23 //==============================================================================
24 
25 //______________________________________________________________________________
26 // OpenGL renderer class for TEveEllipsoid.
27 //
28 
29 
30 
31 //______________________________________________________________________________
33  TGLObject(), fE(0)
34 {
35  // Constructor.
36 
37  // fDLCache = kFALSE; // Disable display list.
38 }
39 
40 //______________________________________________________________________________
41 Bool_t TEveEllipsoidGL::SetModel(TObject* obj, const Option_t* /*opt*/)
42 {
43  // Set model object.
44 
45  fE = SetModelDynCast<TEveEllipsoid>(obj);
46  return kTRUE;
47 }
48 
49 //______________________________________________________________________________
51 {
52  // Set bounding box.
53  ( (TEveEllipsoid*)fExternalObj)->ComputeBBox();
54  SetAxisAlignedBBox(((TEveEllipsoid*)fExternalObj)->AssertBBox());
55 }
56 
57 namespace {
58 GLUquadric* quad = 0; // !!!! AMT check why TGLQuadric crashes on mac
59 }
60 
61 //______________________________________________________________________________
62 void TEveEllipsoidGL::DirectDraw(TGLRnrCtx& /*rnrCtx*/) const
63 {
64  // Render with OpenGL.
65 
66  // printf("TEveEllipsoidGL::DirectDraw LOD %s\n", fE->GetName());
67 
68  glPushAttrib(GL_ENABLE_BIT | GL_POLYGON_BIT | GL_LIGHTING_BIT);
69  glEnable(GL_NORMALIZE );
70  if(!quad)
71  quad = gluNewQuadric();
72 
73  glPushMatrix();
74 
75  TMatrixDSym xxx(3);
76  for(int i=0;i<3;i++)
77  for(int j=0;j<3;j++)
78  {
79  xxx(i,j) = fE->RefEMtx()(i+1,j+1);
80  }
81  TMatrixDEigen eig(xxx);
82 
83 
84  // rewrite for multmatrix ....
85  TEveTrans x;
86  for(int i=0;i<3;i++)
87  for(int j=0;j<3;j++)
88  {
89  x(i+1, j+1) = eig.GetEigenVectors()(i,j);
90  }
91 
92  TVector3 a = x.GetBaseVec(1);
93  TVector3 c = a.Cross(x.GetBaseVec(2));
94  x.SetBaseVec(3, c);
95 
96  glTranslatef(fE->RefPos()[0], fE->RefPos()[1], fE->RefPos()[2]);
97  glMultMatrixd(x.Array());
98  glScalef(fE->RefExtent3D()[0] , fE->RefExtent3D()[1], fE->RefExtent3D()[2]);
99  gluSphere(quad,1. , 30, 30);
100 
101  glPopMatrix();
102  glPopAttrib();
103 
104 
105  // gluDeleteQuadric(quad);
106 }
107 
108 
109 //==============================================================================
110 // TEveEllipsoidProjectedGL
111 //==============================================================================
112 
113 //______________________________________________________________________________
114 // OpenGL renderer class for TEveEllipsoidProjected.
115 //
116 
117 
118 
119 //______________________________________________________________________________
121  fM(0)
122 {
123  // Constructor.
124 
125  // fDLCache = kFALSE; // Disable display list.
126  fMultiColor = kTRUE;
127 }
128 
129 //______________________________________________________________________________
130 Bool_t TEveEllipsoidProjectedGL::SetModel(TObject* obj, const Option_t* /*opt*/)
131 {
132  // Set model object.
133 
134  fM = SetModelDynCast<TEveEllipsoidProjected>(obj);
135  fE = dynamic_cast<TEveEllipsoid*>(fM->GetProjectable());
136  return fE != 0;
137 }
138 
139 //______________________________________________________________________________
141 {
142  // Set bounding box.
143 
144  SetAxisAlignedBBox(((TEveEllipsoidProjected*)fExternalObj)->AssertBBox());
145 }
146 
147 //______________________________________________________________________________
148 void TEveEllipsoidProjectedGL::DirectDraw(TGLRnrCtx& rnrCtx) const
149 {
150  // Render with OpenGL.
151 
152  TEveProjection *proj = fM->GetManager()->GetProjection();
153 
154 
155  glPushAttrib(GL_ENABLE_BIT| GL_POLYGON_BIT | GL_LINE_BIT | GL_POINT_BIT);
156  glDisable(GL_LIGHTING);
157  glDisable(GL_CULL_FACE);
158 
159  glPushMatrix();
160  if ( proj->GetType() == TEveProjection::kPT_RPhi)
161  DrawRhoPhi();
162  else
163  DrawRhoZ();
164 
165  glPopMatrix();
166  glPopAttrib();
167 }
168 
169 //-------------------------------------------------------------------------------
170 
171 void TEveEllipsoidProjectedGL::drawArch(float phiStart, float phiEnd, float phiStep, TEveVector& v0, TEveVector& v1, TEveVector& v2) const
172 {
173  TEveProjection *proj = fM->GetManager()->GetProjection();
174  float phi = phiStart;
175  while (phi < phiEnd ) {
176  TEveVector v = v0 + v1*((float)cos(phi)) + v2*((float)sin(phi));
177  proj->ProjectVector(v, fM->fDepth);
178  glVertex3fv(v.Arr());
179 
180  phi += phiStep;
181  }
182  TEveVector v = v0 + v1*((float)cos(phiEnd)) + v2*((float)sin(phiEnd));
183  proj->ProjectVector(v, fM->fDepth);
184  glVertex3fv(v.Arr());
185 }
186 
187  //-------------------------------------------------------------------------------
189 {
190  // printf("TEveEllipsoidProjectedGL::DirectDraw [%s ]\n", fE->GetName() );
191 
192  TMatrixDSym xxx(3);
193  for(int i=0;i<2;i++)
194  for(int j=0;j<2;j++)
195  {
196  xxx(i,j) = fE->RefEMtx()(i+1,j+1);
197  }
198 
199  TMatrixDEigen eig(xxx);
200  TVectorD xxxEig ( eig.GetEigenValuesRe());
201 
202 
203  // Projection supports only floats :(
204  TEveVector v0(fE->RefPos()[0], fE->RefPos()[1], 0);
205  TEveVector v1(eig.GetEigenVectors()(0, 0) , eig.GetEigenVectors()(0, 1), 0 );
206  v1 *= fE->fEScale * sqrt(TMath::Abs(xxxEig[0]));
207  TEveVector v2(eig.GetEigenVectors()(1, 0) , eig.GetEigenVectors()(1, 1), 0 );
208  v2 *= fE->fEScale * sqrt(TMath::Abs(xxxEig[1]));
209 
210  TEveProjection *proj = fM->GetManager()->GetProjection();
211 
212  // fill
213  glBegin(GL_POLYGON);
214  drawArch(0, TMath::TwoPi(), TMath::TwoPi()/20, v0, v1, v2);
215  glEnd();
216 
217  // frame
218  TGLUtil::LineWidth(fE->fLineWidth);
219  TGLUtil::Color(fE->fLineColor);
220 
221  glBegin(GL_LINE_LOOP);
222  drawArch(0, TMath::TwoPi(), TMath::TwoPi()/20, v0, v1, v2);
223  glEnd();
224 
225  glBegin(GL_LINES);
226  {
227  // glColor3f(1, 0, 0);
228  TEveVector p1 = v0 - v1;
229  TEveVector p2 = v0 + v1;
230  proj->ProjectVector(p1, fM->fDepth);
231  proj->ProjectVector(p2, fM->fDepth);
232  glVertex3fv(p1.Arr());
233  glVertex3fv(p2.Arr());
234  }
235  {
236  // glColor3f(0, 1, 0);
237  TEveVector p1 = v0 - v2;
238  TEveVector p2 = v0 + v2;
239  proj->ProjectVector(p1, fM->fDepth);
240  proj->ProjectVector(p2, fM->fDepth);
241  glVertex3fv(p1.Arr());
242  glVertex3fv(p2.Arr());
243  }
244  glEnd();
245 
246 
247 
248 }
249 
250 //--------------------------------------------------------------------
252 {
253  // printf("TEveEllipsoidProjectedGL::DirectDraw [%s ]\n", fE->GetTitle() );
254 
255  TEveVector v0(fE->RefPos()[0], fE->RefPos()[1], fE->RefPos()[2]);
256 
257  TMatrixDSym xxx(3);
258  xxx(0, 0) = 1;
259  for(int i=1;i<3;i++)
260  for(int j=1;j<3;j++)
261  {
262  xxx(i,j) = fE->RefEMtx()(i+1,j+1);
263  }
264  TMatrixDEigen eig(xxx);
265  TVectorD xxxEig ( eig.GetEigenValuesRe());
266 
267 
268  TEveVector v1(0, eig.GetEigenVectors()(1, 2), eig.GetEigenVectors()(2, 2));
269  v1 *= fE->fEScale * sqrt(TMath::Abs(xxxEig[2]));
270 
271  TEveVector v2(0, eig.GetEigenVectors()(1, 1), eig.GetEigenVectors()(2, 1));
272  v2 *= fE->fEScale * sqrt(TMath::Abs(xxxEig[1]));
273  if (v1[1]*v2[2] > v1[2]*v2[1])
274  v2 *= -1;
275 
276 
277  TEveProjection *proj = fM->GetManager()->GetProjection();
278 
279  // ellipse intersection with projection center
280  bool splitted = false;
281  int N = 20;
282  double phiStep = TMath::TwoPi()/N;
283 
284  // projection center can be moved in beam-spot
285  float bs = 0;
286  if (proj->GetDisplaceOrigin())
287  bs = proj->fCenter[1];
288 
289  float da = v2[1]*v2[1] + v1[1]*v1[1];
290  float db = 2 * v1[1] * (v0[1]-bs);
291  float dc = (v0[1]-bs)*(v0[1]-bs) - v2[1]*v2[1];
292  float disc = (db*db -4*da*dc);
293 
294  if ( disc > 0) {
295  disc = sqrt(disc);
296  float cosS1 = ( -db + disc)/(2 * da);
297  float cosS2 = ( -db - disc)/(2 * da);
298  if (TMath::Abs(cosS1) < 1) {
299  splitted = true;
300  // printf("splitted \n");
301 
302  double phi1 = acos(cosS1);
303  double phi2 = acos(cosS2);
304  TEveVector ps1 = v0 + v1*((float)cos(phi1)) + v2*((float)sin(phi1));
305  TEveVector ps2 = v0 + v1*((float)cos(phi2)) + v2*((float)sin(phi2));
306 
307 
308 
309  // acos has values [0, Pi] , check the symetry over x axis (mirroring)
310  if (TMath::Abs(ps1[1] - bs) > 1e-5)
311  phi1 = TMath::TwoPi() -phi1;
312 
313  if (TMath::Abs(ps2[1] - bs) > 1e-5)
314  phi2 = TMath::TwoPi() -phi2;
315 
316  int N = 20;
317  double phiStep = TMath::TwoPi()/N;
318  double phiOffset = phiStep*0.1;
319  double phiMin = TMath::Min(phi1, phi2);
320  double phiMax = TMath::Max(phi1, phi2);
321  // printf(" %f %f \n",phi1*TMath::RadToDeg(), phi2*TMath::RadToDeg() );
322 
323  // fill
324  // upper clothing
325  glBegin(GL_POLYGON);
326  drawArch(phiMin + phiOffset, phiMax - phiOffset, phiStep, v0, v1, v2);
327  glEnd();
328  // bottom clothing
329  glBegin(GL_POLYGON);
330  drawArch(phiMax + phiOffset, phiMax + TMath::TwoPi() - (phiMax -phiMin) - phiOffset, phiStep, v0, v1, v2);
331  glEnd();
332 
333  // frame
334  TGLUtil::LineWidth(fE->fLineWidth);
335  TGLUtil::Color(fE->fLineColor);
336  // upper clothing
337  glBegin(GL_LINE_LOOP);
338  drawArch(phiMin + phiOffset, phiMax - phiOffset, phiStep, v0, v1, v2);
339  glEnd();
340  // bottom clothing
341  glBegin(GL_LINE_LOOP);
342  drawArch(phiMax + phiOffset, phiMax + TMath::TwoPi() - (phiMax -phiMin) - phiOffset, phiStep, v0, v1, v2);
343  glEnd();
344 
345  }
346  }
347 
348 
349  if (!splitted) {
350  glBegin(GL_POLYGON);
351  drawArch(0, TMath::TwoPi(), phiStep, v0, v1, v2);
352  glEnd();
353  TGLUtil::LineWidth(fE->fLineWidth);
354  TGLUtil::Color(fE->fLineColor);
355  glBegin(GL_LINE_LOOP);
356  drawArch(0, TMath::TwoPi(), phiStep, v0, v1, v2);
357  glEnd();
358  }
359 
360  drawRhoZAxis(v0, v2);
361  drawRhoZAxis(v0, v1);
362 }
363 //______________________________________________________________________________
364 void TEveEllipsoidProjectedGL::drawRhoZAxis(TEveVector& v0, TEveVector& v2) const
365 {
366  glBegin(GL_LINES);
367  TEveProjection* proj = fM->GetManager()->GetProjection();
368 
369  float bs = 0;
370  if (proj->GetDisplaceOrigin())
371  bs = proj->fCenter[1];
372 
373  float off = (v2[1] > v0[1] ) ? 0.01 : -0.01;
374  TEveVector alu = v0 + v2;
375  proj->ProjectVector(alu, fM->fDepth);
376  glVertex3fv(alu.Arr());
377 
378  if (TMath::Abs(v0[1]/v2[1]) < 1 )
379  {
380  alu = v0 - ((float) ((1-off) *(v0[1]-bs)/v2[1])) * v2;
381  proj->ProjectVector(alu, fM->fDepth);
382  glVertex3fv(alu.Arr());
383 
384  //============================
385 
386  alu = v0 - ((float) ((1+off) * (v0[1]-bs)/v2[1])) * v2;
387  proj->ProjectVector(alu, fM->fDepth);
388  glVertex3fv(alu.Arr());
389  }
390 
391  alu = v0 - v2;
392  proj->ProjectVector(alu, fM->fDepth);
393  glVertex3fv(alu.Arr());
394 
395  glEnd();
396 }
const double TwoPi
int i
Definition: DBlmapReader.cc:9
TEveVector & RefPos()
Definition: TEveEllipsoid.h:38
void drawArch(float pStart, float pEnd, float phiStep, TEveVector &v0, TEveVector &v1, TEveVector &v2) const
virtual void DirectDraw(TGLRnrCtx &rnrCtx) const
virtual Bool_t SetModel(TObject *obj, const Option_t *opt=0)
TEveEllipsoid * fE
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
tuple db
Definition: EcalCondDB.py:151
virtual Bool_t SetModel(TObject *obj, const Option_t *opt=0)
void drawRhoZAxis(TEveVector &v, TEveVector &) const
T sqrt(T t)
Definition: SSEVec.h:46
TEveTrans & RefEMtx()
Definition: TEveEllipsoid.h:40
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int j
Definition: DBlmapReader.cc:9
virtual void SetBBox()
double p2[4]
Definition: TauolaWrapper.h:90
virtual void DirectDraw(TGLRnrCtx &rnrCtx) const
TEveVector & RefExtent3D()
Definition: TEveEllipsoid.h:39
#define N
Definition: blowfish.cc:9
TEveEllipsoidProjected * fM
double p1[4]
Definition: TauolaWrapper.h:89
double a
Definition: hdecay.h:121
x
Definition: VDTMath.h:216
mathSSE::Vec4< T > v
Definition: DDAxes.h:10