source: trunk/packages/vizservers/vtkvis/RpPolyData.cpp @ 2402

Last change on this file since 2402 was 2402, checked in by ldelgass, 13 years ago
  • Let graphics objects handle DataSet? cumulative range changes, track vectors as well as scalars, also supply cumulative ranges in setDataSet()
  • Be more consistent about naming enums and commands for vectors
  • Add constructor arguments to some graphics objects to speed initialization (eliminates some pipeline changes)
  • Apply a scale factor to glyphs based on cell sizes
  • Add line glyph shape
  • Don't delete ColorMaps? in use
  • Update graphics objects when a ColorMap? is edited
  • Property svn:eol-style set to native
File size: 4.5 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2011, Purdue Research Foundation
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 <vtkDelaunay2D.h>
17#include <vtkDelaunay3D.h>
18#include <vtkDataSetSurfaceFilter.h>
19
20#include "RpPolyData.h"
21#include "Trace.h"
22
23using namespace Rappture::VtkVis;
24
25PolyData::PolyData() :
26    VtkGraphicsObject()
27{
28    _color[0] = 0.0f;
29    _color[1] = 0.0f;
30    _color[2] = 1.0f;
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 Create and initialize a VTK Prop to render a mesh
45 */
46void PolyData::initProp()
47{
48    if (_prop == NULL) {
49        _prop = vtkSmartPointer<vtkActor>::New();
50        vtkProperty *property = getActor()->GetProperty();
51        property->EdgeVisibilityOn();
52        property->SetColor(_color[0], _color[1], _color[2]);
53        property->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
54        property->SetLineWidth(_edgeWidth);
55        property->SetOpacity(_opacity);
56        property->SetAmbient(.2);
57        if (!_lighting)
58            property->LightingOff();
59    }
60}
61
62/**
63 * \brief Internal method to set up pipeline after a state change
64 */
65void PolyData::update()
66{
67    if (_dataSet == NULL) {
68        return;
69    }
70
71    if (_pdMapper == NULL) {
72        _pdMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
73        _pdMapper->SetResolveCoincidentTopologyToPolygonOffset();
74        _pdMapper->ScalarVisibilityOff();
75    }
76
77    vtkDataSet *ds = _dataSet->getVtkDataSet();
78    vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
79    if (pd) {
80        TRACE("Verts: %d Lines: %d Polys: %d Strips: %d",
81                  pd->GetNumberOfVerts(),
82                  pd->GetNumberOfLines(),
83                  pd->GetNumberOfPolys(),
84                  pd->GetNumberOfStrips());
85        // DataSet is a vtkPolyData
86        if (pd->GetNumberOfLines() == 0 &&
87            pd->GetNumberOfPolys() == 0 &&
88            pd->GetNumberOfStrips() == 0) {
89            // DataSet is a point cloud
90            if (_dataSet->is2D()) {
91                vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
92                mesher->SetInput(pd);
93                mesher->ReleaseDataFlagOn();
94                _pdMapper->SetInputConnection(mesher->GetOutputPort());
95#if defined(DEBUG) && defined(WANT_TRACE)
96                mesher->Update();
97                vtkPolyData *outpd = mesher->GetOutput();
98                TRACE("Delaunay2D Verts: %d Lines: %d Polys: %d Strips: %d",
99                      outpd->GetNumberOfVerts(),
100                      outpd->GetNumberOfLines(),
101                      outpd->GetNumberOfPolys(),
102                      outpd->GetNumberOfStrips());
103#endif
104            } else {
105                vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
106                mesher->SetInput(pd);
107                mesher->ReleaseDataFlagOn();
108                // Delaunay3D returns an UnstructuredGrid, so feed it through a surface filter
109                // to get the grid boundary as a PolyData
110                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
111                gf->UseStripsOn();
112                gf->SetInputConnection(mesher->GetOutputPort());
113                gf->ReleaseDataFlagOn();
114                _pdMapper->SetInputConnection(gf->GetOutputPort());
115            }
116        } else {
117            // DataSet is a vtkPolyData with lines and/or polygons
118            _pdMapper->SetInput(pd);
119        }
120    } else {
121        // DataSet is NOT a vtkPolyData
122        TRACE("DataSet is not a PolyData");
123        vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
124        gf->UseStripsOn();
125        gf->SetInput(ds);
126        gf->ReleaseDataFlagOn();
127        _pdMapper->SetInputConnection(gf->GetOutputPort());
128    }
129
130    initProp();
131    getActor()->SetMapper(_pdMapper);
132    _pdMapper->Update();
133}
134
135/**
136 * \brief Set a group of world coordinate planes to clip rendering
137 *
138 * Passing NULL for planes will remove all cliping planes
139 */
140void PolyData::setClippingPlanes(vtkPlaneCollection *planes)
141{
142    if (_pdMapper != NULL) {
143        _pdMapper->SetClippingPlanes(planes);
144    }
145}
146
Note: See TracBrowser for help on using the repository browser.