Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
REveBox.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 "TClass.h"
13
14#include "ROOT/REveBox.hxx"
17
18using namespace ROOT::Experimental;
19
20/** \class REveBox
21\ingroup REve
223D box with arbitrary vertices (cuboid).
23Vertices 0-3 specify the "bottom" rectangle in clockwise direction and
24vertices 4-7 the "top" rectangle so that 4 is above 0, 5 above 1 and so on.
25
26If vertices are provided some local coordinates the transformation matrix
27of the element should also be set (but then the memory usage is increased
28by the size of the REveTrans object).
29
30Currently only supports 3D -> 2D projections.
31*/
32
33//ClassImp(REveBox);
34
35////////////////////////////////////////////////////////////////////////////////
36/// Constructor.
37
38REveBox::REveBox(const char* n, const char* t) :
39 REveShape(n, t)
40{
41}
42
43////////////////////////////////////////////////////////////////////////////////
44/// Destructor.
45
47{
48}
49
50////////////////////////////////////////////////////////////////////////////////
51/// Set vertex 'i'.
52
54{
55 fVertices[i][0] = x;
56 fVertices[i][1] = y;
57 fVertices[i][2] = z;
58 ResetBBox();
59}
60
61////////////////////////////////////////////////////////////////////////////////
62/// Set vertex 'i'.
63
65{
66 fVertices[i][0] = v[0];
67 fVertices[i][1] = v[1];
68 fVertices[i][2] = v[2];
69 ResetBBox();
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Set vertices.
74
76{
77 memcpy(fVertices, vs, sizeof(fVertices));
78 ResetBBox();
79}
80
81
82////////////////////////////////////////////////////////////////////////////////
83/// Compute bounding-box of the data.
84
86{
88
89 BBoxInit();
90 for (Int_t i=0; i<8; ++i)
91 {
93 }
94}
95
96
97////////////////////////////////////////////////////////////////////////////////
98/// Fill core part of JSON representation.
99
100Int_t REveBox::WriteCoreJson(nlohmann::json &j, Int_t rnr_offset)
101{
102 Int_t ret = REveElement::WriteCoreJson(j, rnr_offset);
103
104 j["fMainColor"] = GetFillColor();
105 j["fLineColor"] = GetLineColor();
106
107 return ret;
108}
109
110////////////////////////////////////////////////////////////////////////////////
111/// Crates 3D point array for rendering.
112
114{
115 int N = 8;
116 fRenderData = std::make_unique<REveRenderData>("makeBox", N*3);
117 for (Int_t i = 0; i < N; ++i)
118 {
119 fRenderData->PushV(fVertices[i][0], fVertices[i][1], fVertices[i][2]);
120 }
121}
122
123
124////////////////////////////////////////////////////////////////////////////////
125/// Virtual from REveProjectable, return REveBoxProjected class.
126
128{
129 return TClass::GetClass<REveBoxProjected>();
130}
131
132
133/** \class REveBoxProjected
134\ingroup REve
135Projection of REveBox.
136*/
137
138
139////////////////////////////////////////////////////////////////////////////////
140/// Constructor.
141
142REveBoxProjected::REveBoxProjected(const char* n, const char* t) :
143 REveShape(n, t),
144 fBreakIdx(0),
145 fDebugCornerPoints(false)
146{
147}
148
149////////////////////////////////////////////////////////////////////////////////
150/// Destructor.
151
153{
154}
155
156////////////////////////////////////////////////////////////////////////////////
157/// Compute bounding-box, virtual from TAttBBox.
158
160{
161 BBoxInit();
162 for (auto &pnt: fPoints)
163 BBoxCheckPoint(pnt.fX, pnt.fY, fDepth);
164}
165
166////////////////////////////////////////////////////////////////////////////////
167/// This is virtual method from base-class REveProjected.
168
170{
171 SetDepthCommon(d, this, fBBox);
172}
173
174////////////////////////////////////////////////////////////////////////////////
175/// This is virtual method from base-class REveProjected.
176
178{
180 CopyVizParams(dynamic_cast<REveElement*>(model));
181}
182
183////////////////////////////////////////////////////////////////////////////////
184/// Re-project the box. Projects all points and finds 2D convex-hull.
185///
186/// The only issue is with making sure that initial conditions for
187/// hull-search are reasonable -- that is, there are no overlaps with the
188/// first point.
189
191{
192 REveBox *box = dynamic_cast<REveBox*>(fProjectable);
193
194 fDebugPoints.clear();
195
196 // Project points in global CS, remove overlaps.
197 vVector2_t pp[2];
198 {
199 REveProjection *projection = fManager->GetProjection();
200 REveTrans *trans = box->PtrMainTrans(kFALSE);
201
202 REveVector pbuf;
203 for (Int_t i = 0; i < 8; ++i)
204 {
205 projection->ProjectPointfv(trans, box->GetVertex(i), pbuf, fDepth);
206 vVector2_t& ppv = pp[projection->SubSpaceId(pbuf)];
207
208 REveVector2 p(pbuf);
209 Bool_t overlap = kFALSE;
210 for (auto &j: ppv)
211 {
213 {
214 overlap = kTRUE;
215 break;
216 }
217 }
218 if (! overlap)
219 {
220 ppv.push_back(p);
222 fDebugPoints.push_back(p);
223 }
224 }
225 }
226
227 fPoints.clear();
228 fBreakIdx = 0;
229
230 if ( ! pp[0].empty())
231 {
232 FindConvexHull(pp[0], fPoints, this);
233 }
234 if ( ! pp[1].empty())
235 {
236 fBreakIdx = fPoints.size();
237 FindConvexHull(pp[1], fPoints, this);
238 }
239}
240
241////////////////////////////////////////////////////////////////////////////////
242/// Get state of fgDebugCornerPoints static.
243
245{
246 return fDebugCornerPoints;
247}
248
249////////////////////////////////////////////////////////////////////////////////
250/// Set state of fgDebugCornerPoints static.
251/// When this is true, points will be drawn at the corners of
252/// computed convex hull.
253
255{
257}
258
259////////////////////////////////////////////////////////////////////////////////
260/// Crates 3D point array for rendering.
261
263{
264 int N = fPoints.size();
265 fRenderData = std::make_unique<REveRenderData>("makeBoxProjected", N*3);
266 for (auto &v : fPoints)
267 {
268 fRenderData->PushV(v.fX);
269 fRenderData->PushV(v.fY);
270 fRenderData->PushV(fDepth);
271 }
272}
273
274////////////////////////////////////////////////////////////////////////////////
275/// Fill core part of JSON representation.
276
277Int_t REveBoxProjected::WriteCoreJson(nlohmann::json &j, Int_t rnr_offset)
278{
279 Int_t ret = REveShape::WriteCoreJson(j, rnr_offset);
280
281 j["fBreakIdx"] = fBreakIdx;
282
283 return ret;
284}
#define d(i)
Definition RSha256.hxx:102
const Bool_t kFALSE
Definition RtypesCore.h:92
const Bool_t kTRUE
Definition RtypesCore.h:91
#define N
void BuildRenderData() override
Crates 3D point array for rendering.
Definition REveBox.cxx:262
Bool_t GetDebugCornerPoints()
Get state of fgDebugCornerPoints static.
Definition REveBox.cxx:244
virtual ~REveBoxProjected()
Destructor.
Definition REveBox.cxx:152
virtual void SetDepthLocal(Float_t d) override
This is virtual method from base-class REveProjected.
Definition REveBox.cxx:169
void UpdateProjection() override
Re-project the box.
Definition REveBox.cxx:190
void SetDebugCornerPoints(Bool_t d)
Set state of fgDebugCornerPoints static.
Definition REveBox.cxx:254
void ComputeBBox() override
Compute bounding-box, virtual from TAttBBox.
Definition REveBox.cxx:159
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveBox.cxx:277
void SetProjection(REveProjectionManager *mng, REveProjectable *model) override
This is virtual method from base-class REveProjected.
Definition REveBox.cxx:177
REveBoxProjected(const REveBoxProjected &)=delete
void SetVertex(Int_t i, Float_t x, Float_t y, Float_t z)
Set vertex 'i'.
Definition REveBox.cxx:53
virtual TClass * ProjectedClass(const REveProjection *p) const override
Virtual from REveProjectable, return REveBoxProjected class.
Definition REveBox.cxx:127
virtual ~REveBox()
Destructor.
Definition REveBox.cxx:46
virtual void ComputeBBox() override
Compute bounding-box of the data.
Definition REveBox.cxx:85
void BuildRenderData() override
Crates 3D point array for rendering.
Definition REveBox.cxx:113
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveBox.cxx:100
void SetVertices(const Float_t *vs)
Set vertices.
Definition REveBox.cxx:75
virtual Int_t WriteCoreJson(nlohmann::json &cj, Int_t rnr_offset)
Write core json.
std::unique_ptr< REveRenderData > fRenderData
Externally assigned and controlled user data.
virtual void SetProjection(REveProjectionManager *mng, REveProjectable *model)
Sets projection manager and reference in the projectable object.
void SetDepthCommon(Float_t d, REveElement *el, Float_t *bbox)
Utility function to update the z-values of the bounding-box.
REveProjectionManager Manager class for steering of projections and managing projected objects.
REveProjection Base for specific classes that implement non-linear projections.
void ProjectPointfv(Float_t *v, Float_t d)
Project float array.
virtual Int_t SubSpaceId(const REveVector &) const
static Int_t FindConvexHull(const vVector2_t &pin, vVector2_t &pout, REveElement *caller=nullptr)
Determines the convex-hull of points in pin.
void CopyVizParams(const REveElement *el) override
Copy visualization parameters from element el.
Definition REveShape.cxx:84
std::vector< REveVector2 > vVector2_t
Definition REveShape.hxx:37
static void CheckAndFixBoxOrientationFv(Float_t box[8][3])
Make sure box orientation is consistent with standard arrangement.
virtual Color_t GetFillColor() const
Definition REveShape.hxx:57
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveShape.cxx:55
virtual Color_t GetLineColor() const
Definition REveShape.hxx:58
REveVector2T A two-vector template without TObject inheritance and virtual functions.
TT SquareDistance(const REveVector2T &v) const
void BBoxCheckPoint(Float_t x, Float_t y, Float_t z)
Definition TAttBBox.h:58
void ResetBBox()
Definition TAttBBox.h:46
void BBoxInit(Float_t infinity=1e6)
Dynamic Float_t[6] X(min,max), Y(min,max), Z(min,max)
Definition TAttBBox.cxx:29
Float_t * fBBox
Definition TAttBBox.h:20
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:80
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16