source: branches/blt4/packages/vizservers/vtkvis/RpHeightMap.cpp @ 2681

Last change on this file since 2681 was 2681, checked in by gah, 12 years ago
File size: 27.0 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 <vtkPointData.h>
12#include <vtkCellData.h>
13#include <vtkCellDataToPointData.h>
14#include <vtkDataSetMapper.h>
15#include <vtkPolyDataMapper.h>
16#include <vtkUnstructuredGrid.h>
17#include <vtkProperty.h>
18#include <vtkImageData.h>
19#include <vtkLookupTable.h>
20#include <vtkTransform.h>
21#include <vtkDelaunay2D.h>
22#include <vtkDelaunay3D.h>
23#include <vtkGaussianSplatter.h>
24#include <vtkExtractVOI.h>
25#include <vtkDataSetSurfaceFilter.h>
26#include <vtkContourFilter.h>
27#include <vtkStripper.h>
28#include <vtkWarpScalar.h>
29#include <vtkPropAssembly.h>
30#include <vtkCutter.h>
31#include <vtkPlane.h>
32
33#include "RpHeightMap.h"
34#include "RpVtkRenderer.h"
35#include "Trace.h"
36
37using namespace Rappture::VtkVis;
38
39#define MESH_POINT_CLOUDS
40
41HeightMap::HeightMap(int numContours, double heightScale) :
42    VtkGraphicsObject(),
43    _numContours(numContours),
44    _colorMap(NULL),
45    _contourEdgeWidth(1.0),
46    _warpScale(heightScale),
47    _sliceAxis(Z_AXIS),
48    _pipelineInitialized(false)
49{
50    _contourEdgeColor[0] = 1.0f;
51    _contourEdgeColor[1] = 0.0f;
52    _contourEdgeColor[2] = 0.0f;
53}
54
55HeightMap::HeightMap(const std::vector<double>& contours, double heightScale) :
56    VtkGraphicsObject(),
57    _numContours(contours.size()),
58    _contours(contours),
59    _colorMap(NULL),
60    _contourEdgeWidth(1.0),
61    _warpScale(heightScale),
62    _sliceAxis(Z_AXIS),
63    _pipelineInitialized(false)
64{
65    _contourEdgeColor[0] = 1.0f;
66    _contourEdgeColor[1] = 0.0f;
67    _contourEdgeColor[2] = 0.0f;
68}
69
70HeightMap::~HeightMap()
71{
72#ifdef WANT_TRACE
73    if (_dataSet != NULL)
74        TRACE("Deleting HeightMap for %s", _dataSet->getName().c_str());
75    else
76        TRACE("Deleting HeightMap with NULL DataSet");
77#endif
78}
79
80void HeightMap::setDataSet(DataSet *dataSet,
81                           Renderer *renderer)
82{
83    if (_dataSet != dataSet) {
84        _dataSet = dataSet;
85
86        if (_dataSet != NULL) {
87            if (renderer->getUseCumulativeRange()) {
88                renderer->getCumulativeDataRange(_dataRange,
89                                                 _dataSet->getActiveScalarsName(),
90                                                 1);
91            } else {
92                _dataSet->getScalarRange(_dataRange);
93            }
94
95            // Compute a data scaling factor to make maximum
96            // height (at default _warpScale) equivalent to
97            // largest dimension of data set
98            double bounds[6];
99            double boundsRange = 0.0;
100            _dataSet->getBounds(bounds);
101            for (int i = 0; i < 6; i+=2) {
102                double r = bounds[i+1] - bounds[i];
103                if (r > boundsRange)
104                    boundsRange = r;
105            }
106            double datRange = _dataRange[1] - _dataRange[0];
107            if (datRange < 1.0e-10) {
108                _dataScale = 1.0;
109            } else {
110                _dataScale = boundsRange / datRange;
111            }
112            TRACE("Bounds range: %g data range: %g scaling: %g",
113                  boundsRange, datRange, _dataScale);
114        }
115
116        update();
117    }
118}
119
120/**
121 * \brief Create and initialize VTK Props to render the colormapped dataset
122 */
123void HeightMap::initProp()
124{
125    if (_dsActor == NULL) {
126        _dsActor = vtkSmartPointer<vtkActor>::New();
127        _dsActor->GetProperty()->SetOpacity(_opacity);
128        _dsActor->GetProperty()->SetEdgeColor(_edgeColor[0],
129                                              _edgeColor[1],
130                                              _edgeColor[2]);
131        _dsActor->GetProperty()->SetLineWidth(_edgeWidth);
132        _dsActor->GetProperty()->EdgeVisibilityOff();
133        _dsActor->GetProperty()->SetAmbient(.2);
134        _dsActor->GetProperty()->LightingOn();
135    }
136    if (_contourActor == NULL) {
137        _contourActor = vtkSmartPointer<vtkActor>::New();
138        _contourActor->GetProperty()->SetOpacity(_opacity);
139        _contourActor->GetProperty()->SetColor(_contourEdgeColor[0],
140                                               _contourEdgeColor[1],
141                                               _contourEdgeColor[2]);
142        _contourActor->GetProperty()->SetEdgeColor(_contourEdgeColor[0],
143                                                   _contourEdgeColor[1],
144                                                   _contourEdgeColor[2]);
145        _contourActor->GetProperty()->SetLineWidth(_contourEdgeWidth);
146        _contourActor->GetProperty()->EdgeVisibilityOff();
147        _contourActor->GetProperty()->SetAmbient(.2);
148        _contourActor->GetProperty()->LightingOff();
149    }
150    if (_prop == NULL) {
151        _prop = vtkSmartPointer<vtkAssembly>::New();
152        getAssembly()->AddPart(_dsActor);
153        getAssembly()->AddPart(_contourActor);
154    }
155}
156
157/**
158 * \brief Internal method to set up pipeline after a state change
159 */
160void HeightMap::update()
161{
162    if (_dataSet == NULL)
163        return;
164
165    vtkDataSet *ds = _dataSet->getVtkDataSet();
166
167    TRACE("DataSet type: %s", _dataSet->getVtkType());
168
169    // Contour filter to generate isolines
170    if (_contourFilter == NULL) {
171        _contourFilter = vtkSmartPointer<vtkContourFilter>::New();
172    }
173
174    // Mapper, actor to render color-mapped data set
175    if (_dsMapper == NULL) {
176        _dsMapper = vtkSmartPointer<vtkDataSetMapper>::New();
177        // Map scalars through lookup table regardless of type
178        _dsMapper->SetColorModeToMapScalars();
179    }
180
181    if (!_pipelineInitialized) {
182        vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
183
184        if (ds->GetPointData() == NULL ||
185            ds->GetPointData()->GetScalars() == NULL) {
186            WARN("No scalar point data in dataset %s", _dataSet->getName().c_str());
187            if (ds->GetCellData() != NULL &&
188                ds->GetCellData()->GetScalars() != NULL) {
189                cellToPtData =
190                    vtkSmartPointer<vtkCellDataToPointData>::New();
191                cellToPtData->SetInput(ds);
192                //cellToPtData->PassCellDataOn();
193                cellToPtData->Update();
194                ds = cellToPtData->GetOutput();
195            } else {
196                ERROR("No scalar cell data in dataset %s", _dataSet->getName().c_str());
197                return;
198            }
199        }
200
201        vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
202        if (pd != NULL) {
203            // DataSet is a vtkPolyData
204            if (pd->GetNumberOfLines() == 0 &&
205                pd->GetNumberOfPolys() == 0 &&
206                pd->GetNumberOfStrips() == 0) {
207                // DataSet is a point cloud
208                PrincipalPlane plane;
209                double offset;
210                if (_dataSet->is2D(&plane, &offset)) {
211#ifdef MESH_POINT_CLOUDS
212                    // Result of Delaunay2D is a PolyData
213                    vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
214                    if (plane == PLANE_ZY) {
215                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
216                        trans->RotateWXYZ(90, 0, 1, 0);
217                        if (offset != 0.0) {
218                            trans->Translate(-offset, 0, 0);
219                        }
220                        mesher->SetTransform(trans);
221                        _sliceAxis = X_AXIS;
222                    } else if (plane == PLANE_XZ) {
223                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
224                        trans->RotateWXYZ(-90, 1, 0, 0);
225                        if (offset != 0.0) {
226                            trans->Translate(0, -offset, 0);
227                        }
228                        mesher->SetTransform(trans);
229                        _sliceAxis = Y_AXIS;
230                    } else if (offset != 0.0) {
231                        // XY with Z offset
232                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
233                        trans->Translate(0, 0, -offset);
234                        mesher->SetTransform(trans);
235                    }
236                    mesher->SetInput(pd);
237                    vtkAlgorithmOutput *warpOutput = initWarp(mesher->GetOutputPort());
238                    _dsMapper->SetInputConnection(warpOutput);
239                    _contourFilter->SetInputConnection(warpOutput);
240#else
241                    if (_pointSplatter == NULL)
242                        _pointSplatter = vtkSmartPointer<vtkGaussianSplatter>::New();
243                    if (_volumeSlicer == NULL)
244                        _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
245                    _pointSplatter->SetInput(pd);
246                    int dims[3];
247                    _pointSplatter->GetSampleDimensions(dims);
248                    TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
249                    if (plane == PLANE_ZY) {
250                        dims[0] = 3;
251                        _volumeSlicer->SetVOI(1, 1, 0, dims[1]-1, 0, dims[1]-1);
252                        _sliceAxis = X_AXIS;
253                    } else if (plane == PLANE_XZ) {
254                        dims[1] = 3;
255                        _volumeSlicer->SetVOI(0, dims[0]-1, 1, 1, 0, dims[2]-1);
256                        _sliceAxis = Y_AXIS;
257                    } else {
258                        dims[2] = 3;
259                        _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
260                    }
261                    _pointSplatter->SetSampleDimensions(dims);
262                    double bounds[6];
263                    _pointSplatter->Update();
264                    _pointSplatter->GetModelBounds(bounds);
265                    TRACE("Model bounds: %g %g %g %g %g %g",
266                          bounds[0], bounds[1],
267                          bounds[2], bounds[3],
268                          bounds[4], bounds[5]);
269                    _volumeSlicer->SetInputConnection(_pointSplatter->GetOutputPort());
270                    _volumeSlicer->SetSampleRate(1, 1, 1);
271                    vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
272                    gf->UseStripsOn();
273                    gf->SetInputConnection(_volumeSlicer->GetOutputPort());
274                    vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
275                    _dsMapper->SetInputConnection(warpOutput);
276                    _contourFilter->SetInputConnection(warpOutput);
277#endif
278                } else {
279#ifdef MESH_POINT_CLOUDS
280                    // Data Set is a 3D point cloud
281                    // Result of Delaunay3D mesher is unstructured grid
282                    vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
283                    mesher->SetInput(pd);
284                    // Run the mesher
285                    mesher->Update();
286                    // Get bounds of resulting grid
287                    double bounds[6];
288                    mesher->GetOutput()->GetBounds(bounds);
289                    // Sample a plane within the grid bounding box
290                    vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
291                    cutter->SetInputConnection(mesher->GetOutputPort());
292                    if (_cutPlane == NULL) {
293                        _cutPlane = vtkSmartPointer<vtkPlane>::New();
294                    }
295                    _cutPlane->SetNormal(0, 0, 1);
296                    _cutPlane->SetOrigin(0,
297                                         0,
298                                         bounds[4] + (bounds[5]-bounds[4])/2.);
299                    cutter->SetCutFunction(_cutPlane);
300                    vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
301                    gf->UseStripsOn();
302                    gf->SetInputConnection(cutter->GetOutputPort());
303#else
304                    if (_pointSplatter == NULL)
305                        _pointSplatter = vtkSmartPointer<vtkGaussianSplatter>::New();
306                    _pointSplatter->SetInput(pd);
307                    int dims[3];
308                    _pointSplatter->GetSampleDimensions(dims);
309                    TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
310                    dims[2] = 3;
311                    _pointSplatter->SetSampleDimensions(dims);
312                    double bounds[6];
313                    _pointSplatter->Update();
314                    _pointSplatter->GetModelBounds(bounds);
315                    TRACE("Model bounds: %g %g %g %g %g %g",
316                          bounds[0], bounds[1],
317                          bounds[2], bounds[3],
318                          bounds[4], bounds[5]);
319                    if (_volumeSlicer == NULL)
320                        _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
321                    _volumeSlicer->SetInputConnection(_pointSplatter->GetOutputPort());
322                    _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
323                    _volumeSlicer->SetSampleRate(1, 1, 1);
324                    vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
325                    gf->UseStripsOn();
326                    gf->SetInputConnection(_volumeSlicer->GetOutputPort());
327#endif
328                    vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
329                    _dsMapper->SetInputConnection(warpOutput);
330                    _contourFilter->SetInputConnection(warpOutput);
331                }
332            } else {
333                // DataSet is a vtkPolyData with lines and/or polygons
334                vtkAlgorithmOutput *warpOutput = initWarp(pd);
335                if (warpOutput != NULL) {
336                    _dsMapper->SetInputConnection(warpOutput);
337                    _contourFilter->SetInputConnection(warpOutput);
338                } else {
339                    _dsMapper->SetInput(pd);
340                    _contourFilter->SetInput(pd);
341                }
342            }
343        } else {
344            // DataSet is NOT a vtkPolyData
345            // Can be: image/volume/uniform grid, structured grid, unstructured grid, rectilinear grid
346            vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
347            gf->UseStripsOn();
348            vtkImageData *imageData = vtkImageData::SafeDownCast(ds);
349            if (!_dataSet->is2D() && imageData != NULL) {
350                // 3D image/volume/uniform grid
351                if (_volumeSlicer == NULL)
352                    _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
353                int dims[3];
354                imageData->GetDimensions(dims);
355                TRACE("Image data dimensions: %d %d %d", dims[0], dims[1], dims[2]);
356                _volumeSlicer->SetInput(ds);
357                _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, (dims[2]-1)/2, (dims[2]-1)/2);
358                _volumeSlicer->SetSampleRate(1, 1, 1);
359                gf->SetInputConnection(_volumeSlicer->GetOutputPort());
360            } else if (!_dataSet->is2D() && imageData == NULL) {
361                // 3D structured grid, unstructured grid, or rectilinear grid
362                double bounds[6];
363                ds->GetBounds(bounds);
364                // Sample a plane within the grid bounding box
365                vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
366                cutter->SetInput(ds);
367                if (_cutPlane == NULL) {
368                    _cutPlane = vtkSmartPointer<vtkPlane>::New();
369                }
370                _cutPlane->SetNormal(0, 0, 1);
371                _cutPlane->SetOrigin(0,
372                                     0,
373                                     bounds[4] + (bounds[5]-bounds[4])/2.);
374                cutter->SetCutFunction(_cutPlane);
375                gf->SetInputConnection(cutter->GetOutputPort());
376            } else {
377                // 2D data
378                gf->SetInput(ds);
379            }
380            vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
381            _dsMapper->SetInputConnection(warpOutput);
382            _contourFilter->SetInputConnection(warpOutput);
383        }
384    }
385
386    _pipelineInitialized = true;
387
388    if (_lut == NULL) {
389        setColorMap(ColorMap::getDefault());
390    }
391    //_dsMapper->InterpolateScalarsBeforeMappingOn();
392
393    initProp();
394
395    _contourFilter->ComputeNormalsOff();
396    _contourFilter->ComputeGradientsOff();
397
398    // Speed up multiple contour computation at cost of extra memory use
399    if (_numContours > 1) {
400        _contourFilter->UseScalarTreeOn();
401    } else {
402        _contourFilter->UseScalarTreeOff();
403    }
404
405    _contourFilter->SetNumberOfContours(_numContours);
406
407    if (_contours.empty()) {
408        // Evenly spaced isovalues
409        _contourFilter->GenerateValues(_numContours, _dataRange[0], _dataRange[1]);
410    } else {
411        // User-supplied isovalues
412        for (int i = 0; i < _numContours; i++) {
413            _contourFilter->SetValue(i, _contours[i]);
414        }
415    }
416    if (_contourMapper == NULL) {
417        _contourMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
418        _contourMapper->ScalarVisibilityOff();
419        _contourMapper->SetResolveCoincidentTopologyToPolygonOffset();
420        vtkSmartPointer<vtkStripper> stripper = vtkSmartPointer<vtkStripper>::New();
421        stripper->SetInputConnection(_contourFilter->GetOutputPort());
422        _contourMapper->SetInputConnection(stripper->GetOutputPort());
423        _contourActor->SetMapper(_contourMapper);
424    }
425
426    _dsActor->SetMapper(_dsMapper);
427
428    _dsMapper->Update();
429    _contourMapper->Update();
430}
431
432vtkAlgorithmOutput *HeightMap::initWarp(vtkAlgorithmOutput *input)
433{
434    TRACE("Warp scale: %g", _warpScale);
435    if (_warpScale == 0.0) {
436        _warp = NULL;
437        return input;
438    } else if (input == NULL) {
439        ERROR("NULL input");
440        return input;
441    } else {
442        if (_warp == NULL)
443            _warp = vtkSmartPointer<vtkWarpScalar>::New();
444        switch (_sliceAxis) {
445        case X_AXIS:
446            _warp->SetNormal(1, 0, 0);
447            break;
448        case Y_AXIS:
449            _warp->SetNormal(0, 1, 0);
450            break;
451        default:
452            _warp->SetNormal(0, 0, 1);
453        }
454        _warp->UseNormalOn();
455        _warp->SetScaleFactor(_warpScale * _dataScale);
456        _warp->SetInputConnection(input);
457        return _warp->GetOutputPort();
458    }
459}
460
461vtkAlgorithmOutput *HeightMap::initWarp(vtkPolyData *pdInput)
462{
463    TRACE("Warp scale: %g", _warpScale);
464    if (_warpScale == 0.0) {
465        _warp = NULL;
466        return NULL;
467    } else if (pdInput == NULL) {
468        ERROR("NULL input");
469        return NULL;
470    } else {
471        if (_warp == NULL)
472            _warp = vtkSmartPointer<vtkWarpScalar>::New();
473        switch (_sliceAxis) {
474        case X_AXIS:
475            _warp->SetNormal(1, 0, 0);
476            break;
477        case Y_AXIS:
478            _warp->SetNormal(0, 1, 0);
479            break;
480        default:
481            _warp->SetNormal(0, 0, 1);
482        }
483        _warp->UseNormalOn();
484        _warp->SetScaleFactor(_warpScale * _dataScale);
485        _warp->SetInput(pdInput);
486        return _warp->GetOutputPort();
487    }
488}
489
490/**
491 * \brief Controls relative scaling of height of mountain plot
492 */
493void HeightMap::setHeightScale(double scale)
494{
495    if (_warpScale == scale)
496        return;
497
498    _warpScale = scale;
499    if (_warp == NULL) {
500        vtkAlgorithmOutput *warpOutput = initWarp(_dsMapper->GetInputConnection(0, 0));
501        _dsMapper->SetInputConnection(warpOutput);
502        _contourFilter->SetInputConnection(warpOutput);
503    } else if (scale == 0.0) {
504        vtkAlgorithmOutput *warpInput = _warp->GetInputConnection(0, 0);
505        _dsMapper->SetInputConnection(warpInput);
506        _contourFilter->SetInputConnection(warpInput);
507        _warp = NULL;
508    } else {
509        _warp->SetScaleFactor(_warpScale * _dataScale);
510    }
511
512    if (_dsMapper != NULL)
513        _dsMapper->Update();
514    if (_contourMapper != NULL)
515        _contourMapper->Update();
516}
517
518/**
519 * \brief Select a 2D slice plane from a 3D DataSet
520 *
521 * \param[in] axis Axis of slice plane
522 * \param[in] ratio Position [0,1] of slice plane along axis
523 */
524void HeightMap::selectVolumeSlice(Axis axis, double ratio)
525{
526    if (_dataSet->is2D()) {
527        WARN("DataSet not 3D, returning");
528        return;
529    }
530
531    if (_volumeSlicer == NULL &&
532        _cutPlane == NULL) {
533        WARN("Called before update() or DataSet is not a volume");
534        return;
535    }
536
537    _sliceAxis = axis;
538    if (_warp != NULL) {
539        switch (axis) {
540        case X_AXIS:
541            _warp->SetNormal(1, 0, 0);
542            break;
543        case Y_AXIS:
544            _warp->SetNormal(0, 1, 0);
545            break;
546        case Z_AXIS:
547            _warp->SetNormal(0, 0, 1);
548            break;
549        default:
550            ERROR("Invalid Axis");
551            return;
552        }
553    }
554
555    if (_cutPlane != NULL) {
556        double bounds[6];
557        _dataSet->getBounds(bounds);
558        switch (axis) {
559        case X_AXIS:
560            _cutPlane->SetNormal(1, 0, 0);
561            _cutPlane->SetOrigin(bounds[0] + (bounds[1]-bounds[0]) * ratio,
562                                 0,
563                                 0);
564            break;
565        case Y_AXIS:
566            _cutPlane->SetNormal(0, 1, 0);
567            _cutPlane->SetOrigin(0,
568                                 bounds[2] + (bounds[3]-bounds[2]) * ratio,
569                                 0);
570            break;
571        case Z_AXIS:
572            _cutPlane->SetNormal(0, 0, 1);
573            _cutPlane->SetOrigin(0,
574                                 0,
575                                 bounds[4] + (bounds[5]-bounds[4]) * ratio);
576            break;
577        default:
578            ERROR("Invalid Axis");
579            return;
580        }
581    } else {
582        int dims[3];
583        if (_pointSplatter != NULL) {
584            _pointSplatter->GetSampleDimensions(dims);
585        } else {
586            vtkImageData *imageData = vtkImageData::SafeDownCast(_dataSet->getVtkDataSet());
587            if (imageData == NULL) {
588                ERROR("Not a volume data set");
589                return;
590            }
591            imageData->GetDimensions(dims);
592        }
593        int voi[6];
594
595        switch (axis) {
596        case X_AXIS:
597            voi[0] = voi[1] = (int)((dims[0]-1) * ratio);
598            voi[2] = 0;
599            voi[3] = dims[1]-1;
600            voi[4] = 0;
601            voi[5] = dims[2]-1;
602            break;
603        case Y_AXIS:
604            voi[0] = 0;
605            voi[1] = dims[0]-1;
606            voi[2] = voi[3] = (int)((dims[1]-1) * ratio);
607            voi[4] = 0;
608            voi[5] = dims[2]-1;
609            break;
610        case Z_AXIS:
611            voi[0] = 0;
612            voi[1] = dims[0]-1;
613            voi[2] = 0;
614            voi[3] = dims[1]-1;
615            voi[4] = voi[5] = (int)((dims[2]-1) * ratio);
616            break;
617        default:
618            ERROR("Invalid Axis");
619            return;
620        }
621
622        _volumeSlicer->SetVOI(voi);
623    }
624
625    if (_dsMapper != NULL)
626        _dsMapper->Update();
627    if (_contourMapper != NULL)
628        _contourMapper->Update();
629}
630
631/**
632 * \brief Called when the color map has been edited
633 */
634void HeightMap::updateColorMap()
635{
636    setColorMap(_colorMap);
637}
638
639/**
640 * \brief Associate a colormap lookup table with the DataSet
641 */
642void HeightMap::setColorMap(ColorMap *cmap)
643{
644    if (cmap == NULL)
645        return;
646
647    _colorMap = cmap;
648 
649    if (_lut == NULL) {
650        _lut = vtkSmartPointer<vtkLookupTable>::New();
651        if (_dsMapper != NULL) {
652            _dsMapper->UseLookupTableScalarRangeOn();
653            _dsMapper->SetLookupTable(_lut);
654        }
655    }
656
657    _lut->DeepCopy(cmap->getLookupTable());
658    _lut->SetRange(_dataRange);
659}
660
661void HeightMap::updateRanges(Renderer *renderer)
662{
663    VtkGraphicsObject::updateRanges(renderer);
664
665    if (_lut != NULL) {
666        _lut->SetRange(_dataRange);
667    }
668
669    if (_contours.empty() && _numContours > 0) {
670        // Need to recompute isovalues
671        update();
672    }
673}
674
675/**
676 * \brief Specify number of evenly spaced contour lines to render
677 *
678 * Will override any existing contours
679 */
680void HeightMap::setNumContours(int numContours)
681{
682    _contours.clear();
683    _numContours = numContours;
684
685    update();
686}
687
688/**
689 * \brief Specify an explicit list of contour isovalues
690 *
691 * Will override any existing contours
692 */
693void HeightMap::setContourList(const std::vector<double>& contours)
694{
695    _contours = contours;
696    _numContours = (int)_contours.size();
697
698    update();
699}
700
701/**
702 * \brief Get the number of contours
703 */
704int HeightMap::getNumContours() const
705{
706    return _numContours;
707}
708
709/**
710 * \brief Get the contour list (may be empty if number of contours
711 * was specified in place of a list)
712 */
713const std::vector<double>& HeightMap::getContourList() const
714{
715    return _contours;
716}
717
718/**
719 * \brief Turn on/off lighting of this object
720 */
721void HeightMap::setLighting(bool state)
722{
723    if (_dsActor != NULL)
724        _dsActor->GetProperty()->SetLighting((state ? 1 : 0));
725}
726
727/**
728 * \brief Turn on/off rendering of mesh edges
729 */
730void HeightMap::setEdgeVisibility(bool state)
731{
732    if (_dsActor != NULL) {
733        _dsActor->GetProperty()->SetEdgeVisibility((state ? 1 : 0));
734    }
735}
736
737/**
738 * \brief Turn on/off rendering of colormaped surface
739 */
740void HeightMap::setContourSurfaceVisibility(bool state)
741{
742    if (_dsActor != NULL) {
743        _dsActor->SetVisibility((state ? 1 : 0));
744    }
745}
746
747/**
748 * \brief Turn on/off rendering of contour isolines
749 */
750void HeightMap::setContourLineVisibility(bool state)
751{
752    if (_contourActor != NULL) {
753        _contourActor->SetVisibility((state ? 1 : 0));
754    }
755}
756
757/**
758 * \brief Set RGB color of polygon edges
759 */
760void HeightMap::setEdgeColor(float color[3])
761{
762    _edgeColor[0] = color[0];
763    _edgeColor[1] = color[1];
764    _edgeColor[2] = color[2];
765    if (_dsActor != NULL)
766        _dsActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
767}
768
769/**
770 * \brief Set pixel width of polygon edges (may be a no-op)
771 */
772void HeightMap::setEdgeWidth(float edgeWidth)
773{
774    _edgeWidth = edgeWidth;
775    if (_dsActor != NULL)
776        _dsActor->GetProperty()->SetLineWidth(_edgeWidth);
777}
778
779/**
780 * \brief Set RGB color of contour isolines
781 */
782void HeightMap::setContourEdgeColor(float color[3])
783{
784    _contourEdgeColor[0] = color[0];
785    _contourEdgeColor[1] = color[1];
786    _contourEdgeColor[2] = color[2];
787    if (_contourActor != NULL) {
788        _contourActor->GetProperty()->SetColor(_contourEdgeColor[0],
789                                               _contourEdgeColor[1],
790                                               _contourEdgeColor[2]);
791        _contourActor->GetProperty()->SetEdgeColor(_contourEdgeColor[0],
792                                                   _contourEdgeColor[1],
793                                                   _contourEdgeColor[2]);
794    }
795}
796
797/**
798 * \brief Set pixel width of contour isolines (may be a no-op)
799 */
800void HeightMap::setContourEdgeWidth(float edgeWidth)
801{
802    _contourEdgeWidth = edgeWidth;
803    if (_contourActor != NULL)
804        _contourActor->GetProperty()->SetLineWidth(_contourEdgeWidth);
805}
806
807/**
808 * \brief Set a group of world coordinate planes to clip rendering
809 *
810 * Passing NULL for planes will remove all cliping planes
811 */
812void HeightMap::setClippingPlanes(vtkPlaneCollection *planes)
813{
814    if (_dsMapper != NULL) {
815        _dsMapper->SetClippingPlanes(planes);
816    }
817    if (_contourMapper != NULL) {
818        _contourMapper->SetClippingPlanes(planes);
819    }
820}
821
Note: See TracBrowser for help on using the repository browser.