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 <ROOT/REveBox.hxx>
15
16#include "TClass.h"
17
18#include <nlohmann/json.hpp>
19
20using namespace ROOT::Experimental;
21
22/** \class REveBox
23\ingroup REve
243D box with arbitrary vertices (cuboid).
25Vertices 0-3 specify the "bottom" rectangle in clockwise direction and
26vertices 4-7 the "top" rectangle so that 4 is above 0, 5 above 1 and so on.
27
28If vertices are provided some local coordinates the transformation matrix
29of the element should also be set (but then the memory usage is increased
30by the size of the REveTrans object).
31
32Currently only supports 3D -> 2D projections.
33*/
34
35
36////////////////////////////////////////////////////////////////////////////////
37/// Constructor.
38
39REveBox::REveBox(const char* n, const char* t) :
40 REveShape(n, t)
41{
42}
43
44////////////////////////////////////////////////////////////////////////////////
45/// Destructor.
46
50
51////////////////////////////////////////////////////////////////////////////////
52/// Set vertex 'i'.
53
55{
56 fVertices[i][0] = x;
57 fVertices[i][1] = y;
58 fVertices[i][2] = z;
59 ResetBBox();
60}
61
62////////////////////////////////////////////////////////////////////////////////
63/// Set vertex 'i'.
64
66{
67 fVertices[i][0] = v[0];
68 fVertices[i][1] = v[1];
69 fVertices[i][2] = v[2];
70 ResetBBox();
71}
72
73////////////////////////////////////////////////////////////////////////////////
74/// Set vertices.
75
77{
78 memcpy(fVertices, vs, sizeof(fVertices));
79 ResetBBox();
80}
81
82
83////////////////////////////////////////////////////////////////////////////////
84/// Compute bounding-box of the data.
85
87{
89
90 BBoxInit();
91 for (Int_t i=0; i<8; ++i)
92 {
94 }
95}
96
97
98////////////////////////////////////////////////////////////////////////////////
99/// Fill core part of JSON representation.
100
102{
104
105 j["fMainColor"] = GetFillColor();
106 j["fLineColor"] = GetLineColor();
107
108 return ret;
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// Crates 3D point array for rendering.
113
115{
116 int N = 8;
117 fRenderData = std::make_unique<REveRenderData>("makeBox", N*3);
118 for (Int_t i = 0; i < N; ++i)
119 {
120 fRenderData->PushV(fVertices[i][0], fVertices[i][1], fVertices[i][2]);
121 }
122}
123
124
125////////////////////////////////////////////////////////////////////////////////
126/// Virtual from REveProjectable, return REveBoxProjected class.
127
129{
130 return TClass::GetClass<REveBoxProjected>();
131}
132
133
134/** \class REveBoxProjected
135\ingroup REve
136Projection of REveBox.
137*/
138
139
140////////////////////////////////////////////////////////////////////////////////
141/// Constructor.
142
143REveBoxProjected::REveBoxProjected(const char* n, const char* t) :
144 REveShape(n, t),
145 fBreakIdx(0),
146 fDebugCornerPoints(false)
147{
148}
149
150////////////////////////////////////////////////////////////////////////////////
151/// Destructor.
152
156
157////////////////////////////////////////////////////////////////////////////////
158/// Compute bounding-box, virtual from TAttBBox.
159
161{
162 BBoxInit();
163 for (auto &pnt: fPoints)
164 BBoxCheckPoint(pnt.fX, pnt.fY, fDepth);
165}
166
167////////////////////////////////////////////////////////////////////////////////
168/// This is virtual method from base-class REveProjected.
169
174
175////////////////////////////////////////////////////////////////////////////////
176/// This is virtual method from base-class REveProjected.
177
183
184////////////////////////////////////////////////////////////////////////////////
185/// Re-project the box. Projects all points and finds 2D convex-hull.
186///
187/// The only issue is with making sure that initial conditions for
188/// hull-search are reasonable -- that is, there are no overlaps with the
189/// first point.
190
192{
193 REveBox *box = dynamic_cast<REveBox*>(fProjectable);
194
195 fDebugPoints.clear();
196
197 // Project points in global CS, remove overlaps.
198 vVector2_t pp[2];
199 {
201 REveTrans *trans = box->PtrMainTrans(kFALSE);
202
204 for (Int_t i = 0; i < 8; ++i)
205 {
206 projection->ProjectPointfv(trans, box->GetVertex(i), pbuf, fDepth);
207 vVector2_t& ppv = pp[projection->SubSpaceId(pbuf)];
208
211 for (auto &j: ppv)
212 {
213 if (p.SquareDistance(j) < REveProjection::fgEpsSqr)
214 {
215 overlap = kTRUE;
216 break;
217 }
218 }
219 if (! overlap)
220 {
221 ppv.push_back(p);
223 fDebugPoints.push_back(p);
224 }
225 }
226 }
227
228 fPoints.clear();
229 fBreakIdx = 0;
230
231 if ( ! pp[0].empty())
232 {
233 FindConvexHull(pp[0], fPoints, this);
234 }
235 if ( ! pp[1].empty())
236 {
237 fBreakIdx = fPoints.size();
238 FindConvexHull(pp[1], fPoints, this);
239 }
240}
241
242////////////////////////////////////////////////////////////////////////////////
243/// Get state of fgDebugCornerPoints static.
244
249
250////////////////////////////////////////////////////////////////////////////////
251/// Set state of fgDebugCornerPoints static.
252/// When this is true, points will be drawn at the corners of
253/// computed convex hull.
254
259
260////////////////////////////////////////////////////////////////////////////////
261/// Crates 3D point array for rendering.
262
264{
265 int N = fPoints.size();
266 fRenderData = std::make_unique<REveRenderData>("makeBoxProjected", N*3);
267 for (auto &v : fPoints)
268 {
269 fRenderData->PushV(v.fX);
270 fRenderData->PushV(v.fY);
271 fRenderData->PushV(fDepth);
272 }
273}
274
275////////////////////////////////////////////////////////////////////////////////
276/// Fill core part of JSON representation.
277
279{
281
282 j["fBreakIdx"] = fBreakIdx;
283
284 return ret;
285}
#define d(i)
Definition RSha256.hxx:102
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define N
winID h TVirtualViewer3D TVirtualGLPainter p
void BuildRenderData() override
Crates 3D point array for rendering.
Definition REveBox.cxx:263
Bool_t GetDebugCornerPoints()
Get state of fgDebugCornerPoints static.
Definition REveBox.cxx:245
~REveBoxProjected() override
Destructor.
Definition REveBox.cxx:153
void SetDepthLocal(Float_t d) override
This is virtual method from base-class REveProjected.
Definition REveBox.cxx:170
void UpdateProjection() override
Re-project the box.
Definition REveBox.cxx:191
void SetDebugCornerPoints(Bool_t d)
Set state of fgDebugCornerPoints static.
Definition REveBox.cxx:255
void ComputeBBox() override
Compute bounding-box, virtual from TAttBBox.
Definition REveBox.cxx:160
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveBox.cxx:278
void SetProjection(REveProjectionManager *mng, REveProjectable *model) override
This is virtual method from base-class REveProjected.
Definition REveBox.cxx:178
REveBoxProjected(const REveBoxProjected &)=delete
~REveBox() override
Destructor.
Definition REveBox.cxx:47
void SetVertex(Int_t i, Float_t x, Float_t y, Float_t z)
Set vertex 'i'.
Definition REveBox.cxx:54
TClass * ProjectedClass(const REveProjection *p) const override
Virtual from REveProjectable, return REveBoxProjected class.
Definition REveBox.cxx:128
void ComputeBBox() override
Compute bounding-box of the data.
Definition REveBox.cxx:86
void BuildRenderData() override
Crates 3D point array for rendering.
Definition REveBox.cxx:114
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveBox.cxx:101
void SetVertices(const Float_t *vs)
Set vertices.
Definition REveBox.cxx:76
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.
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:82
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.
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveShape.cxx:51
void BBoxCheckPoint(Float_t x, Float_t y, Float_t z)
Definition TAttBBox.h:69
void ResetBBox()
Definition TAttBBox.h:57
void BBoxInit(Float_t infinity=1e6)
Dynamic Float_t[6] X(min,max), Y(min,max), Z(min,max)
Definition TAttBBox.cxx:28
Float_t * fBBox
Definition TAttBBox.h:20
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:84
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