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

Last change on this file since 2542 was 2542, checked in by gah, 13 years ago

update from trunk

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