ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 
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 
42 {
43 }
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 /// Set vertex 'i'.
47 
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 
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 
71 {
72  memcpy(fVertices, vs, sizeof(fVertices));
73  ResetBBox();
74 }
75 
76 
77 ////////////////////////////////////////////////////////////////////////////////
78 /// Compute bounding-box of the data.
79 
81 {
83 
84  BBoxInit();
85  for (Int_t i=0; i<8; ++i)
86  {
88  }
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Virtual from TEveProjectable, return TEveBoxProjected class.
93 
95 {
96  return TEveBoxProjected::Class();
97 }
98 
99 
100 /** \class TEveBoxProjected
101 \ingroup TEve
102 Projection of TEveBox.
103 */
104 
106 
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 
122 {
123 }
124 
125 ////////////////////////////////////////////////////////////////////////////////
126 /// Compute bounding-box, virtual from TAttBBox.
127 
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 
141 {
142  SetDepthCommon(d, this, fBBox);
143 }
144 
145 ////////////////////////////////////////////////////////////////////////////////
146 /// This is virtual method from base-class TEveProjected.
147 
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 
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  {
184  {
185  overlap = kTRUE;
186  break;
187  }
188  }
189  if (! overlap)
190  {
191  ppv.push_back(p);
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 
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 
226 {
228 }
TEveTrans is a 4x4 transformation matrix for homogeneous coordinates stored internally in a column-ma...
Definition: TEveTrans.h:26
static Bool_t fgDebugCornerPoints
Definition: TEveBox.h:72
virtual void ComputeBBox()
Compute bounding-box, virtual from TAttBBox.
Definition: TEveBox.cxx:128
Abstract base-class for 2D/3D shapes.
Definition: TEveShape.h:22
void SetVertex(Int_t i, Float_t x, Float_t y, Float_t z)
Set vertex 'i'.
Definition: TEveBox.cxx:48
virtual ~TEveBoxProjected()
Destructor.
Definition: TEveBox.cxx:121
vVector2_t fPoints
Definition: TEveBox.h:66
virtual Int_t SubSpaceId(const TEveVector &) const
std::vector< TEveVector2 > vVector2_t
Definition: TEveShape.h:33
float Float_t
Definition: RtypesCore.h:53
virtual ~TEveBox()
Destructor.
Definition: TEveBox.cxx:41
3D box with arbitrary vertices (cuboid).
Definition: TEveBox.h:21
Float_t fVertices[8][3]
Definition: TEveBox.h:30
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void SetProjection(TEveProjectionManager *mng, TEveProjectable *model)
Sets projection manager and reference in the projectable object.
virtual void SetDepthLocal(Float_t d)
This is virtual method from base-class TEveProjected.
Definition: TEveBox.cxx:140
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition: fillpatterns.C:1
Projection of TEveBox.
Definition: TEveBox.h:56
Float_t * fBBox
Definition: TAttBBox.h:22
static void CheckAndFixBoxOrientationFv(Float_t box[8][3])
Make sure box orientation is consistent with standard arrangement.
Definition: TEveShape.cxx:261
virtual void SetProjection(TEveProjectionManager *mng, TEveProjectable *model)
This is virtual method from base-class TEveProjected.
Definition: TEveBox.cxx:148
void BBoxCheckPoint(Float_t x, Float_t y, Float_t z)
Definition: TAttBBox.h:60
Double_t x[n]
Definition: legend1.C:17
void Class()
Definition: Class.C:29
int d
Definition: tornado.py:11
Base-class for non-linear projections.
Float_t z[5]
Definition: Ifit.C:16
TEveBox(const TEveBox &)
Manager class for steering of projections and managing projected objects.
vVector2_t fDebugPoints
Definition: TEveBox.h:68
Abstract base-class for non-linear projectable objects.
TThread * t[5]
Definition: threadsh1.C:13
SVector< double, 2 > v
Definition: Dict.h:5
std::vector< TEveVector2 >::iterator vVector2_i
Definition: TEveShape.h:34
TEveProjectable * fProjectable
TEveProjectionManager * fManager
void ResetBBox()
Definition: TAttBBox.h:48
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:81
const Float_t * GetVertex(Int_t i) const
Definition: TEveBox.h:40
Minimal, templated two-vector.
Definition: TEveVector.h:280
static Bool_t GetDebugCornerPoints()
Get state of fgDebugCornerPoints static.
Definition: TEveBox.cxx:215
void ProjectPointfv(Float_t *v, Float_t d)
Project float array.
TEveProjection * GetProjection()
void SetVertices(const Float_t *vs)
Set vertices.
Definition: TEveBox.cxx:70
static Float_t fgEpsSqr
virtual void CopyVizParams(const TEveElement *el)
Copy visualization parameters from element el.
Definition: TEveShape.cxx:70
Double_t y[n]
Definition: legend1.C:17
TT SquareDistance(const TEveVector2T &v) const
Definition: TEveVector.h:356
Int_t fBreakIdx
Definition: TEveBox.h:67
void SetDepthCommon(Float_t d, TEveElement *el, Float_t *bbox)
Utility function to update the z-values of the bounding-box.
static void SetDebugCornerPoints(Bool_t d)
Set state of fgDebugCornerPoints static.
Definition: TEveBox.cxx:225
virtual void ComputeBBox()
Compute bounding-box of the data.
Definition: TEveBox.cxx:80
ClassImp(TEveBox)
static Int_t FindConvexHull(const vVector2_t &pin, vVector2_t &pout, TEveElement *caller=0)
Determines the convex-hull of points in pin.
Definition: TEveShape.cxx:117
virtual TEveTrans * PtrMainTrans(Bool_t create=kTRUE)
Return pointer to main transformation.
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual void UpdateProjection()
Re-project the box.
Definition: TEveBox.cxx:161
const Int_t n
Definition: legend1.C:16
void BBoxInit(Float_t infinity=1e6)
Dynamic Float_t[6] X(min,max), Y(min,max), Z(min,max)
virtual TClass * ProjectedClass(const TEveProjection *p) const
Virtual from TEveProjectable, return TEveBoxProjected class.
Definition: TEveBox.cxx:94
TEveBoxProjected(const TEveBoxProjected &)