source: trunk/packages/vizservers/vtkvis/PolyData.cpp @ 3680

Last change on this file since 3680 was 3680, checked in by ldelgass, 11 years ago

Improved cloud support in vtkvis: handle ugrids with no cells as well as
polydata clouds, add cloudstyle options to some graphics objects.

  • Property svn:eol-style set to native
File size: 6.6 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2004-2012  HUBzero Foundation, LLC
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cassert>
9
10#include <vtkDataSet.h>
11#include <vtkPolyData.h>
12#include <vtkUnstructuredGrid.h>
13#include <vtkPolyDataMapper.h>
14#include <vtkActor.h>
15#include <vtkProperty.h>
16#include <vtkTransform.h>
17#include <vtkDelaunay2D.h>
18#include <vtkDelaunay3D.h>
19#include <vtkDataSetSurfaceFilter.h>
20#include <vtkVertexGlyphFilter.h>
21
22#include "PolyData.h"
23#include "Trace.h"
24
25using namespace VtkVis;
26
27PolyData::PolyData() :
28    GraphicsObject(),
29    _cloudStyle(CLOUD_MESH)
30{
31}
32
33PolyData::~PolyData()
34{
35#ifdef WANT_TRACE
36    if (_dataSet != NULL)
37        TRACE("Deleting PolyData for %s", _dataSet->getName().c_str());
38    else
39        TRACE("Deleting PolyData with NULL DataSet");
40#endif
41}
42
43/**
44 * \brief Internal method to set up pipeline after a state change
45 */
46void PolyData::update()
47{
48    if (_dataSet == NULL) {
49        return;
50    }
51
52    if (_pdMapper == NULL) {
53        _pdMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
54        _pdMapper->SetResolveCoincidentTopologyToPolygonOffset();
55        _pdMapper->ScalarVisibilityOff();
56    }
57
58    vtkDataSet *ds = _dataSet->getVtkDataSet();
59    vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
60    if (pd) {
61        TRACE("Points: %d Verts: %d Lines: %d Polys: %d Strips: %d",
62              pd->GetNumberOfPoints(),
63              pd->GetNumberOfVerts(),
64              pd->GetNumberOfLines(),
65              pd->GetNumberOfPolys(),
66              pd->GetNumberOfStrips());
67    }
68    if (_dataSet->isCloud()) {
69        // DataSet is a point cloud
70        PrincipalPlane plane;
71        double offset;
72        if (_cloudStyle == CLOUD_POINTS ||
73            _dataSet->numDimensions() < 2 || ds->GetNumberOfPoints() < 3) { // 0D or 1D or not enough points to mesh
74            vtkSmartPointer<vtkVertexGlyphFilter> vgf = vtkSmartPointer<vtkVertexGlyphFilter>::New();
75#ifdef USE_VTK6
76            vgf->SetInputData(ds);
77#else
78            vgf->SetInput(ds);
79#endif
80            _pdMapper->SetInputConnection(vgf->GetOutputPort());
81        } else if (_dataSet->is2D(&plane, &offset)) {
82            vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
83            if (plane == PLANE_ZY) {
84                vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
85                trans->RotateWXYZ(90, 0, 1, 0);
86                if (offset != 0.0) {
87                    trans->Translate(-offset, 0, 0);
88                }
89                mesher->SetTransform(trans);
90            } else if (plane == PLANE_XZ) {
91                vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
92                trans->RotateWXYZ(-90, 1, 0, 0);
93                if (offset != 0.0) {
94                    trans->Translate(0, -offset, 0);
95                }
96                mesher->SetTransform(trans);
97            } else if (offset != 0.0) {
98                // XY with Z offset
99                vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
100                trans->Translate(0, 0, -offset);
101                mesher->SetTransform(trans);
102            }
103#ifdef USE_VTK6
104            mesher->SetInputData(ds);
105#else
106            mesher->SetInput(ds);
107#endif
108            mesher->ReleaseDataFlagOn();
109            mesher->Update();
110            vtkPolyData *outpd = mesher->GetOutput();
111            TRACE("Delaunay2D Verts: %d Lines: %d Polys: %d Strips: %d",
112                  outpd->GetNumberOfVerts(),
113                  outpd->GetNumberOfLines(),
114                  outpd->GetNumberOfPolys(),
115                  outpd->GetNumberOfStrips());
116            if (outpd->GetNumberOfPolys() == 0) {
117                WARN("Delaunay2D mesher failed");
118                vtkSmartPointer<vtkVertexGlyphFilter> vgf = vtkSmartPointer<vtkVertexGlyphFilter>::New();
119#ifdef USE_VTK6
120                vgf->SetInputData(ds);
121#else
122                vgf->SetInput(ds);
123#endif
124                _pdMapper->SetInputConnection(vgf->GetOutputPort());
125            } else {
126                _pdMapper->SetInputConnection(mesher->GetOutputPort());
127            }
128        } else {
129            vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
130#ifdef USE_VTK6
131            mesher->SetInputData(ds);
132#else
133            mesher->SetInput(ds);
134#endif
135            mesher->ReleaseDataFlagOn();
136            mesher->Update();
137            vtkUnstructuredGrid *grid = mesher->GetOutput();
138            TRACE("Delaunay3D Cells: %d",
139                  grid->GetNumberOfCells());
140            if (grid->GetNumberOfCells() == 0) {
141                WARN("Delaunay3D mesher failed");
142                vtkSmartPointer<vtkVertexGlyphFilter> vgf = vtkSmartPointer<vtkVertexGlyphFilter>::New();
143#ifdef USE_VTK6
144                vgf->SetInputData(ds);
145#else
146                vgf->SetInput(ds);
147#endif
148                _pdMapper->SetInputConnection(vgf->GetOutputPort());
149            } else {
150                // Delaunay3D returns an UnstructuredGrid, so feed it
151                // through a surface filter to get the grid boundary
152                // as a PolyData
153                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
154                gf->UseStripsOn();
155                gf->SetInputConnection(mesher->GetOutputPort());
156                gf->ReleaseDataFlagOn();
157                _pdMapper->SetInputConnection(gf->GetOutputPort());
158            }
159        }
160    } else if (pd) {
161        // DataSet is a vtkPolyData with cells
162#ifdef USE_VTK6
163        _pdMapper->SetInputData(pd);
164#else
165        _pdMapper->SetInput(pd);
166#endif
167    } else {
168        // DataSet is NOT a vtkPolyData
169        TRACE("DataSet is not a PolyData");
170        vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
171        gf->UseStripsOn();
172#ifdef USE_VTK6
173        gf->SetInputData(ds);
174#else
175        gf->SetInput(ds);
176#endif
177        gf->ReleaseDataFlagOn();
178        _pdMapper->SetInputConnection(gf->GetOutputPort());
179    }
180
181    initProp();
182    getActor()->SetMapper(_pdMapper);
183    _pdMapper->Update();
184#ifdef WANT_TRACE
185    double *b = getBounds();
186    TRACE("bounds: %g %g %g %g %g %g", b[0], b[1], b[2], b[3], b[4], b[5]);
187#endif
188}
189
190/**
191 * \brief Set a group of world coordinate planes to clip rendering
192 *
193 * Passing NULL for planes will remove all cliping planes
194 */
195void PolyData::setClippingPlanes(vtkPlaneCollection *planes)
196{
197    if (_pdMapper != NULL) {
198        _pdMapper->SetClippingPlanes(planes);
199    }
200}
201
202void PolyData::setCloudStyle(CloudStyle style)
203{
204    if (style != _cloudStyle) {
205        _cloudStyle = style;
206        if (_dataSet != NULL) {
207            update();
208        }
209    }
210}
Note: See TracBrowser for help on using the repository browser.