source: vtkvis/trunk/HeightMap.cpp @ 6372

Last change on this file since 6372 was 6162, checked in by ldelgass, 9 years ago

merge r5843:5844 from 1.8 branch

  • Property svn:eol-style set to native
File size: 41.6 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2004-2012  HUBzero Foundation, LLC
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cassert>
9
10#include <vtkDataSet.h>
11#include <vtkPointData.h>
12#include <vtkCellData.h>
13#include <vtkCellDataToPointData.h>
14#include <vtkPolyDataMapper.h>
15#include <vtkUnstructuredGrid.h>
16#include <vtkProperty.h>
17#include <vtkImageData.h>
18#include <vtkLookupTable.h>
19#include <vtkTransform.h>
20#include <vtkDelaunay2D.h>
21#include <vtkDelaunay3D.h>
22#include <vtkGaussianSplatter.h>
23#include <vtkExtractVOI.h>
24#include <vtkDataSetSurfaceFilter.h>
25#include <vtkContourFilter.h>
26#include <vtkStripper.h>
27#include <vtkWarpScalar.h>
28#include <vtkPropAssembly.h>
29#include <vtkCutter.h>
30#include <vtkPlane.h>
31
32#include "HeightMap.h"
33#include "Renderer.h"
34#include "Trace.h"
35
36using namespace VtkVis;
37
38#define TRANSLATE_TO_ZORIGIN
39
40HeightMap::HeightMap(int numContours, double heightScale) :
41    GraphicsObject(),
42    _numContours(numContours),
43    _contourColorMap(false),
44    _contourEdgeWidth(1.0),
45    _warpScale(heightScale),
46    _sliceAxis(Z_AXIS),
47    _pipelineInitialized(false),
48    _cloudStyle(CLOUD_MESH),
49    _colorMap(NULL),
50    _colorMode(COLOR_BY_SCALAR),
51    _colorFieldType(DataSet::POINT_DATA),
52    _renderer(NULL)
53{
54    _contourEdgeColor[0] = 0.0f;
55    _contourEdgeColor[1] = 0.0f;
56    _contourEdgeColor[2] = 0.0f;
57    _colorFieldRange[0] = DBL_MAX;
58    _colorFieldRange[1] = -DBL_MAX;
59}
60
61HeightMap::HeightMap(const std::vector<double>& contours, double heightScale) :
62    GraphicsObject(),
63    _numContours(contours.size()),
64    _contours(contours),
65    _contourColorMap(false),
66    _contourEdgeWidth(1.0),
67    _warpScale(heightScale),
68    _sliceAxis(Z_AXIS),
69    _pipelineInitialized(false),
70    _cloudStyle(CLOUD_MESH),
71    _colorMap(NULL),
72    _colorMode(COLOR_BY_SCALAR),
73    _colorFieldType(DataSet::POINT_DATA),
74    _renderer(NULL)
75{
76    _contourEdgeColor[0] = 0.0f;
77    _contourEdgeColor[1] = 0.0f;
78    _contourEdgeColor[2] = 0.0f;
79    _colorFieldRange[0] = DBL_MAX;
80    _colorFieldRange[1] = -DBL_MAX;
81}
82
83HeightMap::~HeightMap()
84{
85#ifdef WANT_TRACE
86    if (_dataSet != NULL)
87        TRACE("Deleting HeightMap for %s", _dataSet->getName().c_str());
88    else
89        TRACE("Deleting HeightMap with NULL DataSet");
90#endif
91}
92
93void HeightMap::setDataSet(DataSet *dataSet,
94                           Renderer *renderer)
95{
96    if (_dataSet != dataSet) {
97        _dataSet = dataSet;
98        _renderer = renderer;
99
100        if (_dataSet != NULL) {
101            TRACE("DataSet name: '%s' type: %s",
102                  _dataSet->getName().c_str(),
103                  _dataSet->getVtkType());
104
105            if (renderer->getUseCumulativeRange()) {
106                renderer->getCumulativeDataRange(_dataRange,
107                                                 _dataSet->getActiveScalarsName(),
108                                                 1);
109                const char *activeVectors = _dataSet->getActiveVectorsName();
110                if (activeVectors != NULL) {
111                    renderer->getCumulativeDataRange(_vectorMagnitudeRange,
112                                                     activeVectors,
113                                                     3);
114                    for (int i = 0; i < 3; i++) {
115                        renderer->getCumulativeDataRange(_vectorComponentRange[i],
116                                                         activeVectors,
117                                                         3, i);
118                    }
119                }
120            } else {
121                _dataSet->getScalarRange(_dataRange);
122                _dataSet->getVectorRange(_vectorMagnitudeRange);
123                for (int i = 0; i < 3; i++) {
124                    _dataSet->getVectorRange(_vectorComponentRange[i], i);
125                }
126            }
127        }
128
129        update();
130    }
131}
132
133/**
134 * Compute a data scaling factor to make maximum height (at default _warpScale)
135 * equivalent to largest dimension of data set
136 */
137void HeightMap::computeDataScale()
138{
139    if (_dataSet == NULL)
140        return;
141
142    double bounds[6];
143    double boundsRange = 0.0;
144    _dataSet->getBounds(bounds);
145    for (int i = 0; i < 6; i+=2) {
146        double r = bounds[i+1] - bounds[i];
147        if (r > boundsRange)
148            boundsRange = r;
149    }
150    double datRange = _dataRange[1] - _dataRange[0];
151    if (datRange < 1.0e-17) {
152        _dataScale = 1.0;
153    } else {
154        _dataScale = boundsRange / datRange;
155    }
156    TRACE("Bounds range: %g data range: %g scaling: %g",
157          boundsRange, datRange, _dataScale);
158}
159
160/**
161 * \brief Create and initialize VTK Props to render the colormapped dataset
162 */
163void HeightMap::initProp()
164{
165    if (_dsActor == NULL) {
166        _dsActor = vtkSmartPointer<vtkActor>::New();
167        _dsActor->GetProperty()->SetOpacity(_opacity);
168        _dsActor->GetProperty()->SetColor(_color[0],
169                                          _color[1],
170                                          _color[2]);
171        _dsActor->GetProperty()->SetEdgeColor(_edgeColor[0],
172                                              _edgeColor[1],
173                                              _edgeColor[2]);
174        _dsActor->GetProperty()->SetLineWidth(_edgeWidth);
175        _dsActor->GetProperty()->EdgeVisibilityOff();
176        _dsActor->GetProperty()->SetAmbient(0.2);
177        _dsActor->GetProperty()->LightingOn();
178    }
179    if (_contourActor == NULL) {
180        _contourActor = vtkSmartPointer<vtkActor>::New();
181        _contourActor->GetProperty()->SetOpacity(_opacity);
182        _contourActor->GetProperty()->SetColor(_contourEdgeColor[0],
183                                               _contourEdgeColor[1],
184                                               _contourEdgeColor[2]);
185        _contourActor->GetProperty()->SetEdgeColor(_contourEdgeColor[0],
186                                                   _contourEdgeColor[1],
187                                                   _contourEdgeColor[2]);
188        _contourActor->GetProperty()->SetLineWidth(_contourEdgeWidth);
189        _contourActor->GetProperty()->EdgeVisibilityOff();
190        _contourActor->GetProperty()->SetAmbient(1.0);
191        _contourActor->GetProperty()->LightingOff();
192    }
193    // Contour actor is added in update() if contours are produced
194    if (_prop == NULL) {
195        _prop = vtkSmartPointer<vtkAssembly>::New();
196        getAssembly()->AddPart(_dsActor);
197    } else {
198        getAssembly()->RemovePart(_contourActor);
199    }
200}
201
202/**
203 * \brief Internal method to set up pipeline after a state change
204 */
205void HeightMap::update()
206{
207    if (_dataSet == NULL)
208        return;
209
210    TRACE("DataSet: %s", _dataSet->getName().c_str());
211
212    vtkDataSet *ds = _dataSet->getVtkDataSet();
213
214    computeDataScale();
215
216    // Contour filter to generate isolines
217    if (_contourFilter == NULL) {
218        _contourFilter = vtkSmartPointer<vtkContourFilter>::New();
219    }
220
221    // Mapper, actor to render color-mapped data set
222    if (_mapper == NULL) {
223        _mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
224        // Map scalars through lookup table regardless of type
225        _mapper->SetColorModeToMapScalars();
226    }
227
228    if (!_pipelineInitialized) {
229        vtkAlgorithmOutput *dsOutput = NULL;
230        vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
231
232        if (ds->GetPointData() == NULL ||
233            ds->GetPointData()->GetScalars() == NULL) {
234            TRACE("No scalar point data in dataset %s", _dataSet->getName().c_str());
235            if (ds->GetCellData() != NULL &&
236                ds->GetCellData()->GetScalars() != NULL) {
237                cellToPtData =
238                    vtkSmartPointer<vtkCellDataToPointData>::New();
239                cellToPtData->SetInputData(ds);
240                cellToPtData->PassCellDataOn();
241                dsOutput = cellToPtData->GetOutputPort();
242            } else {
243                USER_ERROR("No scalar field was found in the data set.");
244                return;
245            }
246        } else {
247            dsOutput = _dataSet->getProducerPort();
248        }
249
250        _splatter = NULL;
251
252        if (_dataSet->isCloud()) {
253            // DataSet is a point cloud
254            PrincipalPlane plane;
255            double offset;
256            if (_dataSet->is2D(&plane, &offset)) {
257                if (_cloudStyle == CLOUD_MESH) {
258                    // Result of Delaunay2D is a PolyData
259                    vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
260                    if (plane == PLANE_ZY) {
261                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
262                        trans->RotateWXYZ(90, 0, 1, 0);
263                        if (offset != 0.0) {
264                            trans->Translate(-offset, 0, 0);
265                        }
266                        mesher->SetTransform(trans);
267                        _sliceAxis = X_AXIS;
268                    } else if (plane == PLANE_XZ) {
269                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
270                        trans->RotateWXYZ(-90, 1, 0, 0);
271                        if (offset != 0.0) {
272                            trans->Translate(0, -offset, 0);
273                        }
274                        mesher->SetTransform(trans);
275                        _sliceAxis = Y_AXIS;
276                    } else if (offset != 0.0) {
277                        // XY with Z offset
278                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
279                        trans->Translate(0, 0, -offset);
280                        mesher->SetTransform(trans);
281                    }
282                    mesher->SetInputConnection(dsOutput);
283                    vtkAlgorithmOutput *warpOutput = initWarp(mesher->GetOutputPort());
284                    _mapper->SetInputConnection(warpOutput);
285                    _contourFilter->SetInputConnection(warpOutput);
286                } else {
287                    // _cloudStyle == CLOUD_SPLAT
288                    if (_splatter == NULL)
289                        _splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
290                    if (_volumeSlicer == NULL)
291                        _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
292                    _splatter->SetInputConnection(dsOutput);
293                    int dims[3];
294                    _splatter->GetSampleDimensions(dims);
295                    TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
296                    if (plane == PLANE_ZY) {
297                        dims[0] = 3;
298                        _volumeSlicer->SetVOI(1, 1, 0, dims[1]-1, 0, dims[1]-1);
299                        _sliceAxis = X_AXIS;
300                    } else if (plane == PLANE_XZ) {
301                        dims[1] = 3;
302                        _volumeSlicer->SetVOI(0, dims[0]-1, 1, 1, 0, dims[2]-1);
303                        _sliceAxis = Y_AXIS;
304                    } else {
305                        dims[2] = 3;
306                        _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
307                    }
308                    _splatter->SetSampleDimensions(dims);
309                    double bounds[6];
310                    _splatter->Update();
311                    _splatter->GetModelBounds(bounds);
312                    TRACE("Model bounds: %g %g %g %g %g %g",
313                          bounds[0], bounds[1],
314                          bounds[2], bounds[3],
315                          bounds[4], bounds[5]);
316                    _volumeSlicer->SetInputConnection(_splatter->GetOutputPort());
317                    _volumeSlicer->SetSampleRate(1, 1, 1);
318                    vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
319                    gf->UseStripsOn();
320                    gf->SetInputConnection(_volumeSlicer->GetOutputPort());
321                    vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
322                    _mapper->SetInputConnection(warpOutput);
323                    _contourFilter->SetInputConnection(warpOutput);
324                }
325            } else {
326                // 3D point cloud
327                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
328                gf->UseStripsOn();
329                if (_cloudStyle == CLOUD_MESH) {
330                     // Result of Delaunay3D mesher is unstructured grid
331                    vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
332                    mesher->SetInputConnection(dsOutput);
333                    // Run the mesher
334                    mesher->Update();
335                    // Get bounds of resulting grid
336                    double bounds[6];
337                    mesher->GetOutput()->GetBounds(bounds);
338                    // Sample a plane within the grid bounding box
339                    vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
340                    cutter->SetInputConnection(mesher->GetOutputPort());
341                    if (_cutPlane == NULL) {
342                        _cutPlane = vtkSmartPointer<vtkPlane>::New();
343                    }
344                    _cutPlane->SetNormal(0, 0, 1);
345                    _cutPlane->SetOrigin(0,
346                                         0,
347                                         bounds[4] + (bounds[5]-bounds[4])/2.);
348                    cutter->SetCutFunction(_cutPlane);
349                    gf->SetInputConnection(cutter->GetOutputPort());
350                } else {
351                    // _cloudStyle == CLOUD_SPLAT
352                    if (_splatter == NULL)
353                        _splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
354                    _splatter->SetInputConnection(dsOutput);
355                    int dims[3];
356                    _splatter->GetSampleDimensions(dims);
357                    TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
358                    dims[2] = 3;
359                    _splatter->SetSampleDimensions(dims);
360                    double bounds[6];
361                    _splatter->Update();
362                    _splatter->GetModelBounds(bounds);
363                    TRACE("Model bounds: %g %g %g %g %g %g",
364                          bounds[0], bounds[1],
365                          bounds[2], bounds[3],
366                          bounds[4], bounds[5]);
367                    if (_volumeSlicer == NULL)
368                        _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
369                    _volumeSlicer->SetInputConnection(_splatter->GetOutputPort());
370                    _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
371                    _volumeSlicer->SetSampleRate(1, 1, 1);
372                    gf->SetInputConnection(_volumeSlicer->GetOutputPort());
373                }
374                vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
375                _mapper->SetInputConnection(warpOutput);
376                _contourFilter->SetInputConnection(warpOutput);
377            }
378        } else if (vtkPolyData::SafeDownCast(ds) != NULL) {
379            // DataSet is a vtkPolyData with lines and/or polygons
380            vtkAlgorithmOutput *warpOutput = initWarp(dsOutput);
381            _mapper->SetInputConnection(warpOutput);
382            _contourFilter->SetInputConnection(warpOutput);
383        } else {
384            // DataSet is NOT a vtkPolyData
385            // Can be: image/volume/uniform grid, structured grid, unstructured grid, rectilinear grid
386            vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
387            gf->UseStripsOn();
388            vtkImageData *imageData = vtkImageData::SafeDownCast(ds);
389            if (!_dataSet->is2D() && imageData != NULL) {
390                // 3D image/volume/uniform grid
391                if (_volumeSlicer == NULL)
392                    _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
393                int dims[3];
394                imageData->GetDimensions(dims);
395                TRACE("Image data dimensions: %d %d %d", dims[0], dims[1], dims[2]);
396                _volumeSlicer->SetInputConnection(dsOutput);
397                _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, (dims[2]-1)/2, (dims[2]-1)/2);
398                _volumeSlicer->SetSampleRate(1, 1, 1);
399                gf->SetInputConnection(_volumeSlicer->GetOutputPort());
400            } else if (!_dataSet->is2D() && imageData == NULL) {
401                // 3D structured grid, unstructured grid, or rectilinear grid
402                double bounds[6];
403                ds->GetBounds(bounds);
404                // Sample a plane within the grid bounding box
405                vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
406                cutter->SetInputConnection(dsOutput);
407                if (_cutPlane == NULL) {
408                    _cutPlane = vtkSmartPointer<vtkPlane>::New();
409                }
410                _cutPlane->SetNormal(0, 0, 1);
411                _cutPlane->SetOrigin(0,
412                                     0,
413                                     bounds[4] + (bounds[5]-bounds[4])/2.);
414                cutter->SetCutFunction(_cutPlane);
415                gf->SetInputConnection(cutter->GetOutputPort());
416            } else {
417                // 2D data
418                gf->SetInputConnection(dsOutput);
419            }
420            vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
421            _mapper->SetInputConnection(warpOutput);
422            _contourFilter->SetInputConnection(warpOutput);
423        }
424    }
425
426    _contourFilter->ComputeNormalsOff();
427    _contourFilter->ComputeGradientsOff();
428
429    // Speed up multiple contour computation at cost of extra memory use
430    if (_numContours > 1) {
431        _contourFilter->UseScalarTreeOn();
432    } else {
433        _contourFilter->UseScalarTreeOff();
434    }
435
436    _contourFilter->SetNumberOfContours(_numContours);
437
438    if (_contours.empty()) {
439        // Evenly spaced isovalues
440        if (_colorFieldRange[0] <= _colorFieldRange[1]) {
441            _contourFilter->GenerateValues(_numContours, _colorFieldRange[0], _colorFieldRange[1]);
442        } else {
443            _contourFilter->GenerateValues(_numContours, _dataRange[0], _dataRange[1]);
444        }
445    } else {
446        // User-supplied isovalues
447        for (int i = 0; i < _numContours; i++) {
448            _contourFilter->SetValue(i, _contours[i]);
449        }
450    }
451
452    initProp();
453
454    if (_contourMapper == NULL) {
455        _contourMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
456        _contourMapper->ScalarVisibilityOff();
457        _contourMapper->SetResolveCoincidentTopologyToPolygonOffset();
458        vtkSmartPointer<vtkStripper> stripper = vtkSmartPointer<vtkStripper>::New();
459        stripper->SetInputConnection(_contourFilter->GetOutputPort());
460        _contourMapper->SetInputConnection(stripper->GetOutputPort());
461        _contourActor->SetMapper(_contourMapper);
462        //_contourActor->InterpolateScalarsBeforeMappingOn();
463    }
464
465    setInterpolateBeforeMapping(true);
466
467    if (_lut == NULL) {
468        setColorMap(ColorMap::getDefault());
469        setColorMode(_colorMode);
470    } else if (!_pipelineInitialized) {
471        double *rangePtr = _colorFieldRange;
472        if (_colorFieldRange[0] > _colorFieldRange[1]) {
473            rangePtr = NULL;
474        }
475        setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
476    }
477
478    _pipelineInitialized = true;
479
480    TRACE("DataSet active scalars: %s", _dataSet->getActiveScalarsName());
481    // Ensure updated dataScale is applied
482    if (_warp != NULL) {
483        _warp->SetScaleFactor(_warpScale * _dataScale);
484    }
485#ifdef TRANSLATE_TO_ZORIGIN
486    double pos[3];
487    pos[0] = 0;
488    pos[1] = 0;
489    pos[2] = - _warpScale * _dataScale * _dataRange[0];
490    setPosition(pos);
491    TRACE("Z_POS: %g", pos[2]);
492#endif
493    _dsActor->SetMapper(_mapper);
494
495    _mapper->Update();
496    _contourMapper->Update();
497
498    // Only add contour actor to assembly if contours were
499    // produced, in order to prevent messing up assembly bounds
500    double bounds[6];
501    _contourActor->GetBounds(bounds);
502    if (bounds[0] <= bounds[1]) {
503        getAssembly()->AddPart(_contourActor);
504    }
505}
506
507void HeightMap::setCloudStyle(CloudStyle style)
508{
509    if (style != _cloudStyle) {
510        _cloudStyle = style;
511        if (_dataSet != NULL) {
512            _pipelineInitialized = false;
513            update();
514        }
515    }
516}
517
518void HeightMap::setInterpolateBeforeMapping(bool state)
519{
520    if (_mapper != NULL) {
521        _mapper->SetInterpolateScalarsBeforeMapping((state ? 1 : 0));
522    }
523}
524
525void HeightMap::setAspect(double aspect)
526{
527    double bounds[6];
528    vtkDataSet *ds = _dataSet->getVtkDataSet();
529    ds->GetBounds(bounds);
530    double size[3];
531    size[0] = bounds[1] - bounds[0];
532    size[1] = bounds[3] - bounds[2];
533    size[2] = bounds[5] - bounds[4];
534    double scale[3];
535    scale[0] = scale[1] = scale[2] = 1.;
536
537    if (aspect == 1.0) {
538        // Square
539        switch (_sliceAxis) {
540        case X_AXIS: {
541            if (size[1] > size[2] && size[2] > 0.0) {
542                scale[2] = size[1] / size[2];
543            } else if (size[2] > size[1] && size[1] > 0.0) {
544                scale[1] = size[2] / size[1];
545            }
546        }
547            break;
548        case Y_AXIS: {
549            if (size[0] > size[2] && size[2] > 0.0) {
550                scale[2] = size[0] / size[2];
551            } else if (size[2] > size[0] && size[0] > 0.0) {
552                scale[0] = size[2] / size[0];
553            }
554        }
555            break;
556        case Z_AXIS: {
557            if (size[0] > size[1] && size[1] > 0.0) {
558                scale[1] = size[0] / size[1];
559            } else if (size[1] > size[0] && size[0] > 0.0) {
560                scale[0] = size[1] / size[0];
561            }
562        }
563            break;
564        }
565    } else if (aspect != 0.0) {
566        switch (_sliceAxis) {
567        case X_AXIS: {
568            if (aspect > 1.0) {
569                if (size[2] > size[1] && size[1] > 0.0) {
570                    scale[1] = (size[2] / aspect) / size[1];
571                } else if (size[2] > 0.0) {
572                    scale[2] = (size[1] * aspect) / size[2];
573                }
574            } else {
575                if (size[1] > size[2] && size[2] > 0.0) {
576                    scale[2] = (size[1] * aspect) / size[2];
577                } else if (size[1] > 0.0) {
578                    scale[1] = (size[2] / aspect) / size[1];
579                }
580            }
581        }
582            break;
583        case Y_AXIS: {
584            if (aspect > 1.0) {
585                if (size[0] > size[2] && size[2] > 0.0) {
586                    scale[2] = (size[0] / aspect) / size[2];
587                } else if (size[0] > 0.0) {
588                    scale[0] = (size[2] * aspect) / size[0];
589                }
590            } else {
591                if (size[2] > size[0] && size[0] > 0.0) {
592                    scale[0] = (size[2] * aspect) / size[0];
593                } else if (size[2] > 0.0) {
594                    scale[2] = (size[0] / aspect) / size[2];
595                }
596            }
597        }
598            break;
599        case Z_AXIS: {
600            if (aspect > 1.0) {
601                if (size[0] > size[1] && size[1] > 0.0) {
602                    scale[1] = (size[0] / aspect) / size[1];
603                } else if (size[0] > 0.0) {
604                    scale[0] = (size[1] * aspect) / size[0];
605                }
606            } else {
607                if (size[1] > size[0] && size[0] > 0.0) {
608                    scale[0] = (size[1] * aspect) / size[0];
609                } else if (size[1] > 0.0) {
610                    scale[1] = (size[0] / aspect) / size[1];
611                }
612            }
613        }
614            break;
615        }
616    }
617
618    TRACE("%s dims %g,%g,%g", _dataSet->getName().c_str(),
619          size[0], size[1], size[2]);
620    TRACE("Setting scale to %g,%g,%g", scale[0], scale[1], scale[2]);
621    setScale(scale);
622}
623
624vtkAlgorithmOutput *HeightMap::initWarp(vtkAlgorithmOutput *input)
625{
626    TRACE("Warp scale: %g", _warpScale);
627    if (_warpScale == 0.0) {
628        _warp = NULL;
629        _normalsGenerator = NULL;
630        return input;
631    } else if (input == NULL) {
632        ERROR("NULL input");
633        return input;
634    } else {
635        if (_warp == NULL)
636            _warp = vtkSmartPointer<vtkWarpScalar>::New();
637        if (_normalsGenerator == NULL)
638            _normalsGenerator = vtkSmartPointer<vtkPolyDataNormals>::New();
639        switch (_sliceAxis) {
640        case X_AXIS:
641            _warp->SetNormal(1, 0, 0);
642            break;
643        case Y_AXIS:
644            _warp->SetNormal(0, 1, 0);
645            break;
646        default:
647            _warp->SetNormal(0, 0, 1);
648        }
649        _warp->UseNormalOn();
650        _warp->SetScaleFactor(_warpScale * _dataScale);
651        _warp->SetInputConnection(input);
652        _normalsGenerator->SetInputConnection(_warp->GetOutputPort());
653        _normalsGenerator->SetFeatureAngle(90.);
654        _normalsGenerator->AutoOrientNormalsOff();
655        _normalsGenerator->ComputePointNormalsOn();
656        return _normalsGenerator->GetOutputPort();
657    }
658}
659
660/**
661 * \brief Controls relative scaling of height of mountain plot
662 */
663void HeightMap::setHeightScale(double scale)
664{
665    if (_warpScale == scale)
666        return;
667
668    _warpScale = scale;
669    if (_warp == NULL && scale != 0.0) {
670        vtkAlgorithmOutput *warpOutput = initWarp(_mapper->GetInputConnection(0, 0));
671        _mapper->SetInputConnection(warpOutput);
672        _contourFilter->SetInputConnection(warpOutput);
673    } else if (scale == 0.0) {
674        vtkAlgorithmOutput *warpInput = _warp->GetInputConnection(0, 0);
675        _mapper->SetInputConnection(warpInput);
676        _contourFilter->SetInputConnection(warpInput);
677        _warp = NULL;
678    } else {
679        _warp->SetScaleFactor(_warpScale * _dataScale);
680    }
681#ifdef TRANSLATE_TO_ZORIGIN
682    double pos[3];
683    pos[0] = 0;
684    pos[1] = 0;
685    pos[2] = - _warpScale * _dataScale * _dataRange[0];
686    setPosition(pos);
687    TRACE("Z_POS: %g", pos[2]);
688#endif
689    if (_mapper != NULL)
690        _mapper->Update();
691    if (_contourMapper != NULL)
692        _contourMapper->Update();
693}
694
695/**
696 * \brief Select a 2D slice plane from a 3D DataSet
697 *
698 * \param[in] axis Axis of slice plane
699 * \param[in] ratio Position [0,1] of slice plane along axis
700 */
701void HeightMap::selectVolumeSlice(Axis axis, double ratio)
702{
703    if (_dataSet->is2D()) {
704        WARN("DataSet not 3D, returning");
705        return;
706    }
707
708    if (_volumeSlicer == NULL &&
709        _cutPlane == NULL) {
710        WARN("Called before update() or DataSet is not a volume");
711        return;
712    }
713
714    _sliceAxis = axis;
715    if (_warp != NULL) {
716        switch (axis) {
717        case X_AXIS:
718            _warp->SetNormal(1, 0, 0);
719            break;
720        case Y_AXIS:
721            _warp->SetNormal(0, 1, 0);
722            break;
723        case Z_AXIS:
724            _warp->SetNormal(0, 0, 1);
725            break;
726        default:
727            ERROR("Invalid Axis");
728            return;
729        }
730    }
731
732    if (_cutPlane != NULL) {
733        double bounds[6];
734        _dataSet->getBounds(bounds);
735        switch (axis) {
736        case X_AXIS:
737            _cutPlane->SetNormal(1, 0, 0);
738            _cutPlane->SetOrigin(bounds[0] + (bounds[1]-bounds[0]) * ratio,
739                                 0,
740                                 0);
741            break;
742        case Y_AXIS:
743            _cutPlane->SetNormal(0, 1, 0);
744            _cutPlane->SetOrigin(0,
745                                 bounds[2] + (bounds[3]-bounds[2]) * ratio,
746                                 0);
747            break;
748        case Z_AXIS:
749            _cutPlane->SetNormal(0, 0, 1);
750            _cutPlane->SetOrigin(0,
751                                 0,
752                                 bounds[4] + (bounds[5]-bounds[4]) * ratio);
753            break;
754        default:
755            ERROR("Invalid Axis");
756            return;
757        }
758    } else {
759        int dims[3];
760        if (_splatter != NULL) {
761            _splatter->GetSampleDimensions(dims);
762        } else {
763            vtkImageData *imageData = vtkImageData::SafeDownCast(_dataSet->getVtkDataSet());
764            if (imageData == NULL) {
765                ERROR("Not a volume data set");
766                return;
767            }
768            imageData->GetDimensions(dims);
769        }
770        int voi[6];
771
772        switch (axis) {
773        case X_AXIS:
774            voi[0] = voi[1] = (int)((dims[0]-1) * ratio);
775            voi[2] = 0;
776            voi[3] = dims[1]-1;
777            voi[4] = 0;
778            voi[5] = dims[2]-1;
779            break;
780        case Y_AXIS:
781            voi[0] = 0;
782            voi[1] = dims[0]-1;
783            voi[2] = voi[3] = (int)((dims[1]-1) * ratio);
784            voi[4] = 0;
785            voi[5] = dims[2]-1;
786            break;
787        case Z_AXIS:
788            voi[0] = 0;
789            voi[1] = dims[0]-1;
790            voi[2] = 0;
791            voi[3] = dims[1]-1;
792            voi[4] = voi[5] = (int)((dims[2]-1) * ratio);
793            break;
794        default:
795            ERROR("Invalid Axis");
796            return;
797        }
798
799        _volumeSlicer->SetVOI(voi);
800    }
801
802    if (_mapper != NULL)
803        _mapper->Update();
804    if (_contourMapper != NULL)
805        _contourMapper->Update();
806}
807
808void HeightMap::updateRanges(Renderer *renderer)
809{
810    TRACE("Enter");
811
812    if (_dataSet == NULL) {
813        ERROR("called before setDataSet");
814        return;
815    }
816
817    if (renderer->getUseCumulativeRange()) {
818        renderer->getCumulativeDataRange(_dataRange,
819                                         _dataSet->getActiveScalarsName(),
820                                         1);
821        const char *activeVectors = _dataSet->getActiveVectorsName();
822        if (activeVectors != NULL) {
823            renderer->getCumulativeDataRange(_vectorMagnitudeRange,
824                                             activeVectors,
825                                             3);
826            for (int i = 0; i < 3; i++) {
827                renderer->getCumulativeDataRange(_vectorComponentRange[i],
828                                                 activeVectors,
829                                                 3, i);
830            }
831        }
832    } else {
833        _dataSet->getScalarRange(_dataRange);
834        _dataSet->getVectorRange(_vectorMagnitudeRange);
835        for (int i = 0; i < 3; i++) {
836            _dataSet->getVectorRange(_vectorComponentRange[i], i);
837        }
838    }
839 
840    // Need to update color map ranges
841    double *rangePtr = _colorFieldRange;
842    if (_colorFieldRange[0] > _colorFieldRange[1]) {
843        rangePtr = NULL;
844    }
845    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
846
847    if ((_contours.empty() && _numContours > 0) ||
848        _warpScale > 0.0) {
849        // Contour isovalues and/or warp need to be recomputed
850        update();
851    } else {
852        computeDataScale();
853    }
854
855    TRACE("Leave");
856}
857
858void HeightMap::setContourLineColorMapEnabled(bool mode)
859{
860    _contourColorMap = mode;
861
862    double *rangePtr = _colorFieldRange;
863    if (_colorFieldRange[0] > _colorFieldRange[1]) {
864        rangePtr = NULL;
865    }
866    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
867}
868
869void HeightMap::setColorMode(ColorMode mode)
870{
871    _colorMode = mode;
872    if (_dataSet == NULL)
873        return;
874
875    switch (mode) {
876    case COLOR_BY_SCALAR:
877        setColorMode(mode,
878                     _dataSet->getActiveScalarsType(),
879                     _dataSet->getActiveScalarsName());
880        break;
881    case COLOR_BY_VECTOR_MAGNITUDE:
882        setColorMode(mode,
883                     _dataSet->getActiveVectorsType(),
884                     _dataSet->getActiveVectorsName());
885        break;
886    case COLOR_BY_VECTOR_X:
887        setColorMode(mode,
888                     _dataSet->getActiveVectorsType(),
889                     _dataSet->getActiveVectorsName());
890        break;
891    case COLOR_BY_VECTOR_Y:
892        setColorMode(mode,
893                     _dataSet->getActiveVectorsType(),
894                     _dataSet->getActiveVectorsName());
895        break;
896    case COLOR_BY_VECTOR_Z:
897        setColorMode(mode,
898                     _dataSet->getActiveVectorsType(),
899                     _dataSet->getActiveVectorsName());
900        break;
901    case COLOR_CONSTANT:
902    default:
903        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
904        break;
905    }
906}
907
908void HeightMap::setColorMode(ColorMode mode,
909                             const char *name, double range[2])
910{
911    if (_dataSet == NULL)
912        return;
913    DataSet::DataAttributeType type = DataSet::POINT_DATA;
914    int numComponents = 1;
915    if (name != NULL && strlen(name) > 0 &&
916        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
917        ERROR("Field not found: %s", name);
918        return;
919    }
920    setColorMode(mode, type, name, range);
921}
922
923void HeightMap::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
924                             const char *name, double range[2])
925{
926    _colorMode = mode;
927    _colorFieldType = type;
928    if (name == NULL)
929        _colorFieldName.clear();
930    else
931        _colorFieldName = name;
932    if (range == NULL) {
933        _colorFieldRange[0] = DBL_MAX;
934        _colorFieldRange[1] = -DBL_MAX;
935    } else {
936        memcpy(_colorFieldRange, range, sizeof(double)*2);
937    }
938
939    if (_dataSet == NULL || _mapper == NULL)
940        return;
941
942    switch (type) {
943    case DataSet::POINT_DATA:
944        _mapper->SetScalarModeToUsePointFieldData();
945        _contourMapper->SetScalarModeToUsePointFieldData();
946        break;
947    case DataSet::CELL_DATA:
948        _mapper->SetScalarModeToUseCellFieldData();
949        _contourMapper->SetScalarModeToUseCellFieldData();
950        break;
951    default:
952        ERROR("Unsupported DataAttributeType: %d", type);
953        return;
954    }
955
956    if (_splatter != NULL) {
957        if (name != NULL && strlen(name) > 0) {
958            _splatter->SetInputArrayToProcess(0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_POINTS, name);
959        }
960        _mapper->SelectColorArray("SplatterValues");
961        _contourMapper->SelectColorArray("SplatterValues");
962    } else if (name != NULL && strlen(name) > 0) {
963        _mapper->SelectColorArray(name);
964        _contourMapper->SelectColorArray(name);
965    } else {
966        _mapper->SetScalarModeToDefault();
967        _contourMapper->SetScalarModeToDefault();
968    }
969
970    if (_lut != NULL) {
971        if (range != NULL) {
972            _lut->SetRange(range);
973            TRACE("range %g, %g", range[0], range[1]);
974        } else if (name != NULL && strlen(name) > 0) {
975            double r[2];
976            int comp = -1;
977            if (mode == COLOR_BY_VECTOR_X)
978                comp = 0;
979            else if (mode == COLOR_BY_VECTOR_Y)
980                comp = 1;
981            else if (mode == COLOR_BY_VECTOR_Z)
982                comp = 2;
983
984            if (_renderer->getUseCumulativeRange()) {
985                int numComponents;
986                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
987                    ERROR("Field not found: %s, type: %d", name, type);
988                    return;
989                } else if (numComponents < comp+1) {
990                    ERROR("Request for component %d in field with %d components",
991                          comp, numComponents);
992                    return;
993                }
994                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
995            } else {
996                _dataSet->getDataRange(r, name, type, comp);
997            }
998            _lut->SetRange(r);
999            TRACE("%s range %g, %g", name, r[0], r[1]);
1000        }
1001        _lut->Modified();
1002    }
1003
1004    if (_contourColorMap) {
1005        _contourMapper->ScalarVisibilityOn();
1006    } else {
1007        _contourMapper->ScalarVisibilityOff();
1008    }
1009
1010    switch (mode) {
1011    case COLOR_BY_SCALAR:
1012        _mapper->ScalarVisibilityOn();
1013        break;
1014    case COLOR_BY_VECTOR_MAGNITUDE:
1015        _mapper->ScalarVisibilityOn();
1016        if (_lut != NULL) {
1017            _lut->SetVectorModeToMagnitude();
1018        }
1019        break;
1020    case COLOR_BY_VECTOR_X:
1021        _mapper->ScalarVisibilityOn();
1022        if (_lut != NULL) {
1023            _lut->SetVectorModeToComponent();
1024            _lut->SetVectorComponent(0);
1025        }
1026        break;
1027    case COLOR_BY_VECTOR_Y:
1028        _mapper->ScalarVisibilityOn();
1029        if (_lut != NULL) {
1030            _lut->SetVectorModeToComponent();
1031            _lut->SetVectorComponent(1);
1032        }
1033        break;
1034    case COLOR_BY_VECTOR_Z:
1035        _mapper->ScalarVisibilityOn();
1036        if (_lut != NULL) {
1037            _lut->SetVectorModeToComponent();
1038            _lut->SetVectorComponent(2);
1039        }
1040        break;
1041    case COLOR_CONSTANT:
1042    default:
1043        _mapper->ScalarVisibilityOff();
1044        break;
1045    }
1046}
1047
1048/**
1049 * \brief Called when the color map has been edited
1050 */
1051void HeightMap::updateColorMap()
1052{
1053    setColorMap(_colorMap);
1054}
1055
1056/**
1057 * \brief Associate a colormap lookup table with the DataSet
1058 */
1059void HeightMap::setColorMap(ColorMap *cmap)
1060{
1061    if (cmap == NULL)
1062        return;
1063
1064    _colorMap = cmap;
1065
1066    if (_lut == NULL) {
1067        _lut = vtkSmartPointer<vtkLookupTable>::New();
1068        if (_mapper != NULL) {
1069            _mapper->UseLookupTableScalarRangeOn();
1070            _mapper->SetLookupTable(_lut);
1071        }
1072        if (_contourMapper != NULL) {
1073            _contourMapper->UseLookupTableScalarRangeOn();
1074            _contourMapper->SetLookupTable(_lut);
1075        }
1076        _lut->DeepCopy(cmap->getLookupTable());
1077        switch (_colorMode) {
1078        case COLOR_CONSTANT:
1079        case COLOR_BY_SCALAR:
1080            _lut->SetRange(_dataRange);
1081            break;
1082        case COLOR_BY_VECTOR_MAGNITUDE:
1083            _lut->SetRange(_vectorMagnitudeRange);
1084            break;
1085        case COLOR_BY_VECTOR_X:
1086            _lut->SetRange(_vectorComponentRange[0]);
1087            break;
1088        case COLOR_BY_VECTOR_Y:
1089            _lut->SetRange(_vectorComponentRange[1]);
1090            break;
1091        case COLOR_BY_VECTOR_Z:
1092            _lut->SetRange(_vectorComponentRange[2]);
1093            break;
1094        default:
1095            break;
1096        }
1097    } else {
1098        double range[2];
1099        _lut->GetTableRange(range);
1100        _lut->DeepCopy(cmap->getLookupTable());
1101        _lut->SetRange(range);
1102        _lut->Modified();
1103    }
1104
1105    switch (_colorMode) {
1106    case COLOR_BY_VECTOR_MAGNITUDE:
1107        _lut->SetVectorModeToMagnitude();
1108        break;
1109    case COLOR_BY_VECTOR_X:
1110        _lut->SetVectorModeToComponent();
1111        _lut->SetVectorComponent(0);
1112        break;
1113    case COLOR_BY_VECTOR_Y:
1114        _lut->SetVectorModeToComponent();
1115        _lut->SetVectorComponent(1);
1116        break;
1117    case COLOR_BY_VECTOR_Z:
1118        _lut->SetVectorModeToComponent();
1119        _lut->SetVectorComponent(2);
1120        break;
1121    default:
1122         break;
1123    }
1124}
1125
1126/**
1127 * \brief Specify number of evenly spaced contour lines to render
1128 *
1129 * Will override any existing contours
1130 */
1131void HeightMap::setNumContours(int numContours)
1132{
1133    _contours.clear();
1134    _numContours = numContours;
1135
1136    update();
1137}
1138
1139/**
1140 * \brief Specify an explicit list of contour isovalues
1141 *
1142 * Will override any existing contours
1143 */
1144void HeightMap::setContourList(const std::vector<double>& contours)
1145{
1146    _contours = contours;
1147    _numContours = (int)_contours.size();
1148
1149    update();
1150}
1151
1152/**
1153 * \brief Get the number of contours
1154 */
1155int HeightMap::getNumContours() const
1156{
1157    return _numContours;
1158}
1159
1160/**
1161 * \brief Get the contour list (may be empty if number of contours
1162 * was specified in place of a list)
1163 */
1164const std::vector<double>& HeightMap::getContourList() const
1165{
1166    return _contours;
1167}
1168
1169/**
1170 * \brief Turn on/off lighting of this object
1171 */
1172void HeightMap::setLighting(bool state)
1173{
1174    _lighting = state;
1175    if (_dsActor != NULL)
1176        _dsActor->GetProperty()->SetLighting((state ? 1 : 0));
1177}
1178
1179/**
1180 * \brief Turn on/off rendering of mesh edges
1181 */
1182void HeightMap::setEdgeVisibility(bool state)
1183{
1184    if (_dsActor != NULL) {
1185        _dsActor->GetProperty()->SetEdgeVisibility((state ? 1 : 0));
1186    }
1187}
1188
1189/**
1190 * \brief Turn on/off rendering of colormaped surface
1191 */
1192void HeightMap::setContourSurfaceVisibility(bool state)
1193{
1194    if (_dsActor != NULL) {
1195        _dsActor->SetVisibility((state ? 1 : 0));
1196    }
1197}
1198
1199/**
1200 * \brief Turn on/off rendering of contour isolines
1201 */
1202void HeightMap::setContourLineVisibility(bool state)
1203{
1204    if (_contourActor != NULL) {
1205        _contourActor->SetVisibility((state ? 1 : 0));
1206    }
1207}
1208
1209/**
1210 * \brief Set RGB color of mesh
1211 */
1212void HeightMap::setColor(float color[3])
1213{
1214    _color[0] = color[0];
1215    _color[1] = color[1];
1216    _color[2] = color[2];
1217    if (_dsActor != NULL)
1218        _dsActor->GetProperty()->SetColor(_color[0], _color[1], _color[2]);
1219}
1220
1221/**
1222 * \brief Set RGB color of polygon edges
1223 */
1224void HeightMap::setEdgeColor(float color[3])
1225{
1226    _edgeColor[0] = color[0];
1227    _edgeColor[1] = color[1];
1228    _edgeColor[2] = color[2];
1229    if (_dsActor != NULL)
1230        _dsActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
1231}
1232
1233/**
1234 * \brief Set pixel width of polygon edges (may be a no-op)
1235 */
1236void HeightMap::setEdgeWidth(float edgeWidth)
1237{
1238    _edgeWidth = edgeWidth;
1239    if (_dsActor != NULL)
1240        _dsActor->GetProperty()->SetLineWidth(_edgeWidth);
1241}
1242
1243/**
1244 * \brief Set RGB color of contour isolines
1245 */
1246void HeightMap::setContourEdgeColor(float color[3])
1247{
1248    _contourEdgeColor[0] = color[0];
1249    _contourEdgeColor[1] = color[1];
1250    _contourEdgeColor[2] = color[2];
1251    if (_contourActor != NULL) {
1252        _contourActor->GetProperty()->SetColor(_contourEdgeColor[0],
1253                                               _contourEdgeColor[1],
1254                                               _contourEdgeColor[2]);
1255        _contourActor->GetProperty()->SetEdgeColor(_contourEdgeColor[0],
1256                                                   _contourEdgeColor[1],
1257                                                   _contourEdgeColor[2]);
1258    }
1259}
1260
1261/**
1262 * \brief Set pixel width of contour isolines (may be a no-op)
1263 */
1264void HeightMap::setContourEdgeWidth(float edgeWidth)
1265{
1266    _contourEdgeWidth = edgeWidth;
1267    if (_contourActor != NULL)
1268        _contourActor->GetProperty()->SetLineWidth(_contourEdgeWidth);
1269}
1270
1271/**
1272 * \brief Set a group of world coordinate planes to clip rendering
1273 *
1274 * Passing NULL for planes will remove all cliping planes
1275 */
1276void HeightMap::setClippingPlanes(vtkPlaneCollection *planes)
1277{
1278    if (_mapper != NULL) {
1279        _mapper->SetClippingPlanes(planes);
1280    }
1281    if (_contourMapper != NULL) {
1282        _contourMapper->SetClippingPlanes(planes);
1283    }
1284}
1285
Note: See TracBrowser for help on using the repository browser.