source: branches/nanovis2/packages/vizservers/vtkvis/RpHeightMap.cpp @ 3305

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

sync with trunk

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