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//ClassImp(REveBox);
36
37////////////////////////////////////////////////////////////////////////////////
38/// Constructor.
39
40REveBox::REveBox(const char* n, const char* t) :
41 REveShape(n, t)
42{
43}
44
45////////////////////////////////////////////////////////////////////////////////
46/// Destructor.
47
49{
50}
51
52////////////////////////////////////////////////////////////////////////////////
53/// Set vertex 'i'.
54
56{
57 fVertices[i][0] = x;
58 fVertices[i][1] = y;
59 fVertices[i][2] = z;
60 ResetBBox();
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Set vertex 'i'.
65
67{
68 fVertices[i][0] = v[0];
69 fVertices[i][1] = v[1];
70 fVertices[i][2] = v[2];
71 ResetBBox();
72}
73
74////////////////////////////////////////////////////////////////////////////////
75/// Set vertices.
76
78{
79 memcpy(fVertices, vs, sizeof(fVertices));
80 ResetBBox();
81}
82
83
84////////////////////////////////////////////////////////////////////////////////
85/// Compute bounding-box of the data.
86
88{
90
91 BBoxInit();
92 for (Int_t i=0; i<8; ++i)
93 {
95 }
96}
97
98
99////////////////////////////////////////////////////////////////////////////////
100/// Fill core part of JSON representation.
101
102Int_t REveBox::WriteCoreJson(nlohmann::json &j, Int_t rnr_offset)
103{
104 Int_t ret = REveElement::WriteCoreJson(j, rnr_offset);
105
106 j["fMainColor"] = GetFillColor();
107 j["fLineColor"] = GetLineColor();
108
109 return ret;
110}
111
112////////////////////////////////////////////////////////////////////////////////
113/// Crates 3D point array for rendering.
114
116{
117 int N = 8;
118 fRenderData = std::make_unique<REveRenderData>("makeBox", N*3);
119 for (Int_t i = 0; i < N; ++i)
120 {
121 fRenderData->PushV(fVertices[i][0], fVertices[i][1], fVertices[i][2]);
122 }
123}
124
125
126////////////////////////////////////////////////////////////////////////////////
127/// Virtual from REveProjectable, return REveBoxProjected class.
128
130{
131 return TClass::GetClass<REveBoxProjected>();
132}
133
134
135/** \class REveBoxProjected
136\ingroup REve
137Projection of REveBox.
138*/
139
140
141////////////////////////////////////////////////////////////////////////////////
142/// Constructor.
143
144REveBoxProjected::REveBoxProjected(const char* n, const char* t) :
145 REveShape(n, t),
146 fBreakIdx(0),
147 fDebugCornerPoints(false)
148{
149}
150
151////////////////////////////////////////////////////////////////////////////////
152/// Destructor.
153
155{
156}
157
158////////////////////////////////////////////////////////////////////////////////
159/// Compute bounding-box, virtual from TAttBBox.
160
162{
163 BBoxInit();
164 for (auto &pnt: fPoints)
165 BBoxCheckPoint(pnt.fX, pnt.fY, fDepth);
166}
167
168////////////////////////////////////////////////////////////////////////////////
169/// This is virtual method from base-class REveProjected.
170
172{
173 SetDepthCommon(d, this, fBBox);
174}
175
176////////////////////////////////////////////////////////////////////////////////
177/// This is virtual method from base-class REveProjected.
178
180{
182 CopyVizParams(dynamic_cast<REveElement*>(model));
183}
184
185////////////////////////////////////////////////////////////////////////////////
186/// Re-project the box. Projects all points and finds 2D convex-hull.
187///
188/// The only issue is with making sure that initial conditions for
189/// hull-search are reasonable -- that is, there are no overlaps with the
190/// first point.
191
193{
194 REveBox *box = dynamic_cast<REveBox*>(fProjectable);
195
196 fDebugPoints.clear();
197
198 // Project points in global CS, remove overlaps.
199 vVector2_t pp[2];
200 {
201 REveProjection *projection = fManager->GetProjection();
202 REveTrans *trans = box->PtrMainTrans(kFALSE);
203
204 REveVector pbuf;
205 for (Int_t i = 0; i < 8; ++i)
206 {
207 projection->ProjectPointfv(trans, box->GetVertex(i), pbuf, fDepth);
208 vVector2_t& ppv = pp[projection->SubSpaceId(pbuf)];
209
210 REveVector2 p(pbuf);
211 Bool_t overlap = kFALSE;
212 for (auto &j: ppv)
213 {
214 if (p.SquareDistance(j) < REveProjection::fgEpsSqr)
215 {
216 overlap = kTRUE;
217 break;
218 }
219 }
220 if (! overlap)
221 {
222 ppv.push_back(p);
224 fDebugPoints.push_back(p);
225 }
226 }
227 }
228
229 fPoints.clear();
230 fBreakIdx = 0;
231
232 if ( ! pp[0].empty())
233 {
234 FindConvexHull(pp[0], fPoints, this);
235 }
236 if ( ! pp[1].empty())
237 {
238 fBreakIdx = fPoints.size();
239 FindConvexHull(pp[1], fPoints, this);
240 }
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// Get state of fgDebugCornerPoints static.
245
247{
248 return fDebugCornerPoints;
249}
250
251////////////////////////////////////////////////////////////////////////////////
252/// Set state of fgDebugCornerPoints static.
253/// When this is true, points will be drawn at the corners of
254/// computed convex hull.
255
257{
259}
260
261////////////////////////////////////////////////////////////////////////////////
262/// Crates 3D point array for rendering.
263
265{
266 int N = fPoints.size();
267 fRenderData = std::make_unique<REveRenderData>("makeBoxProjected", N*3);
268 for (auto &v : fPoints)
269 {
270 fRenderData->PushV(v.fX);
271 fRenderData->PushV(v.fY);
272 fRenderData->PushV(fDepth);
273 }
274}
275
276////////////////////////////////////////////////////////////////////////////////
277/// Fill core part of JSON representation.
278
279Int_t REveBoxProjected::WriteCoreJson(nlohmann::json &j, Int_t rnr_offset)
280{
281 Int_t ret = REveShape::WriteCoreJson(j, rnr_offset);
282
283 j["fBreakIdx"] = fBreakIdx;
284
285 return ret;
286}
#define d(i)
Definition RSha256.hxx:102
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
#define N
winID h TVirtualViewer3D TVirtualGLPainter p
void BuildRenderData() override
Crates 3D point array for rendering.
Definition REveBox.cxx:264
Bool_t GetDebugCornerPoints()
Get state of fgDebugCornerPoints static.
Definition REveBox.cxx:246
~REveBoxProjected() override
Destructor.
Definition REveBox.cxx:154
void SetDepthLocal(Float_t d) override
This is virtual method from base-class REveProjected.
Definition REveBox.cxx:171
void UpdateProjection() override
Re-project the box.
Definition REveBox.cxx:192
void SetDebugCornerPoints(Bool_t d)
Set state of fgDebugCornerPoints static.
Definition REveBox.cxx:256
void ComputeBBox() override
Compute bounding-box, virtual from TAttBBox.
Definition REveBox.cxx:161
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveBox.cxx:279
void SetProjection(REveProjectionManager *mng, REveProjectable *model) override
This is virtual method from base-class REveProjected.
Definition REveBox.cxx:179
REveBoxProjected(const REveBoxProjected &)=delete
~REveBox() override
Destructor.
Definition REveBox.cxx:48
void SetVertex(Int_t i, Float_t x, Float_t y, Float_t z)
Set vertex 'i'.
Definition REveBox.cxx:55
TClass * ProjectedClass(const REveProjection *p) const override
Virtual from REveProjectable, return REveBoxProjected class.
Definition REveBox.cxx:129
void ComputeBBox() override
Compute bounding-box of the data.
Definition REveBox.cxx:87
void BuildRenderData() override
Crates 3D point array for rendering.
Definition REveBox.cxx:115
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
Definition REveBox.cxx:102
void SetVertices(const Float_t *vs)
Set vertices.
Definition REveBox.cxx:77
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:86
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:57
virtual Color_t GetLineColor() const
Definition REveShape.hxx:58
REveVector2T A two-vector template without TObject inheritance and virtual functions.
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:29
Float_t * fBBox
Definition TAttBBox.h:20
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
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