Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TEveBox.cxx
Go to the documentation of this file.
1 // @(#)root/eve:$Id$
2 // Author: Matevz Tadel, 2010
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2007, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "TEveBox.h"
13 #include "TEveProjectionManager.h"
14 
15 /** \class TEveBox
16 \ingroup TEve
17 3D box with arbitrary vertices (cuboid).
18 Vertices 0-3 specify the "bottom" rectangle in clockwise direction and
19 vertices 4-7 the "top" rectangle so that 4 is above 0, 5 above 1 and so on.
20 
21 If vertices are provided some local coordinates the transformation matrix
22 of the element should also be set (but then the memory usage is increased
23 by the size of the TEveTrans object).
24 
25 Currently only supports 3D -> 2D projections.
26 */
27 
28 ClassImp(TEveBox);
29 
30 ////////////////////////////////////////////////////////////////////////////////
31 /// Constructor.
32 
33 TEveBox::TEveBox(const char* n, const char* t) :
34  TEveShape(n, t)
35 {
36 }
37 
38 ////////////////////////////////////////////////////////////////////////////////
39 /// Destructor.
40 
41 TEveBox::~TEveBox()
42 {
43 }
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 /// Set vertex 'i'.
47 
48 void TEveBox::SetVertex(Int_t i, Float_t x, Float_t y, Float_t z)
49 {
50  fVertices[i][0] = x;
51  fVertices[i][1] = y;
52  fVertices[i][2] = z;
53  ResetBBox();
54 }
55 
56 ////////////////////////////////////////////////////////////////////////////////
57 /// Set vertex 'i'.
58 
59 void TEveBox::SetVertex(Int_t i, const Float_t* v)
60 {
61  fVertices[i][0] = v[0];
62  fVertices[i][1] = v[1];
63  fVertices[i][2] = v[2];
64  ResetBBox();
65 }
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 /// Set vertices.
69 
70 void TEveBox::SetVertices(const Float_t* vs)
71 {
72  memcpy(fVertices, vs, sizeof(fVertices));
73  ResetBBox();
74 }
75 
76 
77 ////////////////////////////////////////////////////////////////////////////////
78 /// Compute bounding-box of the data.
79 
80 void TEveBox::ComputeBBox()
81 {
82  TEveShape::CheckAndFixBoxOrientationFv(fVertices);
83 
84  BBoxInit();
85  for (Int_t i=0; i<8; ++i)
86  {
87  BBoxCheckPoint(fVertices[i]);
88  }
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Virtual from TEveProjectable, return TEveBoxProjected class.
93 
94 TClass* TEveBox::ProjectedClass(const TEveProjection*) const
95 {
96  return TEveBoxProjected::Class();
97 }
98 
99 
100 /** \class TEveBoxProjected
101 \ingroup TEve
102 Projection of TEveBox.
103 */
104 
105 ClassImp(TEveBoxProjected);
106 
107 Bool_t TEveBoxProjected::fgDebugCornerPoints = kFALSE;
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 /// Constructor.
111 
112 TEveBoxProjected::TEveBoxProjected(const char* n, const char* t) :
113  TEveShape(n, t),
114  fBreakIdx(0)
115 {
116 }
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 /// Destructor.
120 
121 TEveBoxProjected::~TEveBoxProjected()
122 {
123 }
124 
125 ////////////////////////////////////////////////////////////////////////////////
126 /// Compute bounding-box, virtual from TAttBBox.
127 
128 void TEveBoxProjected::ComputeBBox()
129 {
130  BBoxInit();
131  for (vVector2_i i = fPoints.begin(); i != fPoints.end(); ++i)
132  {
133  BBoxCheckPoint(i->fX, i->fY, fDepth);
134  }
135 }
136 
137 ////////////////////////////////////////////////////////////////////////////////
138 /// This is virtual method from base-class TEveProjected.
139 
140 void TEveBoxProjected::SetDepthLocal(Float_t d)
141 {
142  SetDepthCommon(d, this, fBBox);
143 }
144 
145 ////////////////////////////////////////////////////////////////////////////////
146 /// This is virtual method from base-class TEveProjected.
147 
148 void TEveBoxProjected::SetProjection(TEveProjectionManager* mng, TEveProjectable* model)
149 {
150  TEveProjected::SetProjection(mng, model);
151  CopyVizParams(dynamic_cast<TEveElement*>(model));
152 }
153 
154 ////////////////////////////////////////////////////////////////////////////////
155 /// Re-project the box. Projects all points and finds 2D convex-hull.
156 ///
157 /// The only issue is with making sure that initial conditions for
158 /// hull-search are reasonable -- that is, there are no overlaps with the
159 /// first point.
160 
161 void TEveBoxProjected::UpdateProjection()
162 {
163  TEveBox *box = dynamic_cast<TEveBox*>(fProjectable);
164 
165  fDebugPoints.clear();
166 
167  // Project points in global CS, remove overlaps.
168  vVector2_t pp[2];
169  {
170  TEveProjection *projection = fManager->GetProjection();
171  TEveTrans *trans = box->PtrMainTrans(kFALSE);
172 
173  TEveVector pbuf;
174  for (Int_t i = 0; i < 8; ++i)
175  {
176  projection->ProjectPointfv(trans, box->GetVertex(i), pbuf, fDepth);
177  vVector2_t& ppv = pp[projection->SubSpaceId(pbuf)];
178 
179  TEveVector2 p(pbuf);
180  Bool_t overlap = kFALSE;
181  for (vVector2_i j = ppv.begin(); j != ppv.end(); ++j)
182  {
183  if (p.SquareDistance(*j) < TEveProjection::fgEpsSqr)
184  {
185  overlap = kTRUE;
186  break;
187  }
188  }
189  if (! overlap)
190  {
191  ppv.push_back(p);
192  if (fgDebugCornerPoints)
193  fDebugPoints.push_back(p);
194  }
195  }
196  }
197 
198  fPoints.clear();
199  fBreakIdx = 0;
200 
201  if ( ! pp[0].empty())
202  {
203  FindConvexHull(pp[0], fPoints, this);
204  }
205  if ( ! pp[1].empty())
206  {
207  fBreakIdx = fPoints.size();
208  FindConvexHull(pp[1], fPoints, this);
209  }
210 }
211 
212 ////////////////////////////////////////////////////////////////////////////////
213 /// Get state of fgDebugCornerPoints static.
214 
215 Bool_t TEveBoxProjected::GetDebugCornerPoints()
216 {
217  return fgDebugCornerPoints;
218 }
219 
220 ////////////////////////////////////////////////////////////////////////////////
221 /// Set state of fgDebugCornerPoints static.
222 /// When this is true, points will be drawn at the corners of
223 /// computed convex hull.
224 
225 void TEveBoxProjected::SetDebugCornerPoints(Bool_t d)
226 {
227  fgDebugCornerPoints = d;
228 }