source: trunk/packages/vizservers/vtkvis/RpWarp.cpp @ 3132

Last change on this file since 3132 was 3124, checked in by ldelgass, 12 years ago

Add preliminary version of mesh warping by vector field (e.g. mesh deformation,
flow profile/surfaces)

  • Property svn:eol-style set to native
File size: 17.8 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2012, Purdue Research Foundation
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cassert>
9#include <cfloat>
10#include <cstring>
11
12#include <vtkDataSet.h>
13#include <vtkPointData.h>
14#include <vtkCellData.h>
15#include <vtkCellDataToPointData.h>
16#include <vtkDataSetMapper.h>
17#include <vtkUnstructuredGrid.h>
18#include <vtkProperty.h>
19#include <vtkImageData.h>
20#include <vtkLookupTable.h>
21#include <vtkTransform.h>
22#include <vtkDelaunay2D.h>
23#include <vtkDelaunay3D.h>
24#include <vtkGaussianSplatter.h>
25#include <vtkExtractVOI.h>
26#include <vtkDataSetSurfaceFilter.h>
27
28#include "RpWarp.h"
29#include "RpVtkRenderer.h"
30#include "Trace.h"
31
32#define MESH_POINT_CLOUDS
33
34using namespace Rappture::VtkVis;
35
36Warp::Warp() :
37    VtkGraphicsObject(),
38    _warpScale(1.0),
39    _colorMap(NULL),
40    _colorMode(COLOR_BY_SCALAR),
41    _colorFieldType(DataSet::POINT_DATA),
42    _renderer(NULL)
43{
44    _colorFieldRange[0] = DBL_MAX;
45    _colorFieldRange[1] = -DBL_MAX;
46}
47
48Warp::~Warp()
49{
50#ifdef WANT_TRACE
51    if (_dataSet != NULL)
52        TRACE("Deleting Warp for %s", _dataSet->getName().c_str());
53    else
54        TRACE("Deleting Warp with NULL DataSet");
55#endif
56}
57
58void Warp::setDataSet(DataSet *dataSet,
59                      Renderer *renderer)
60{
61    if (_dataSet != dataSet) {
62        _dataSet = dataSet;
63
64        _renderer = renderer;
65
66        if (renderer->getUseCumulativeRange()) {
67            renderer->getCumulativeDataRange(_dataRange,
68                                             _dataSet->getActiveScalarsName(),
69                                             1);
70            renderer->getCumulativeDataRange(_vectorMagnitudeRange,
71                                             _dataSet->getActiveVectorsName(),
72                                             3);
73            for (int i = 0; i < 3; i++) {
74                renderer->getCumulativeDataRange(_vectorComponentRange[i],
75                                                 _dataSet->getActiveVectorsName(),
76                                                 3, i);
77            }
78        } else {
79            _dataSet->getScalarRange(_dataRange);
80            _dataSet->getVectorRange(_vectorMagnitudeRange);
81            for (int i = 0; i < 3; i++) {
82                _dataSet->getVectorRange(_vectorComponentRange[i], i);
83            }
84        }
85
86        update();
87    }
88}
89
90/**
91 * \brief Internal method to set up pipeline after a state change
92 */
93void Warp::update()
94{
95    if (_dataSet == NULL)
96        return;
97
98    vtkDataSet *ds = _dataSet->getVtkDataSet();
99
100    if (_warp == NULL) {
101        _warp = vtkSmartPointer<vtkWarpVector>::New();
102    }
103
104    // Mapper, actor to render color-mapped data set
105    if (_dsMapper == NULL) {
106        _dsMapper = vtkSmartPointer<vtkDataSetMapper>::New();
107        // Map scalars through lookup table regardless of type
108        _dsMapper->SetColorModeToMapScalars();
109        //_dsMapper->InterpolateScalarsBeforeMappingOn();
110    }
111
112    vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
113    if (ds->GetPointData() == NULL ||
114        ds->GetPointData()->GetVectors() == NULL) {
115        TRACE("No vector point data in dataset %s", _dataSet->getName().c_str());
116        if (ds->GetCellData() != NULL &&
117            ds->GetCellData()->GetVectors() != NULL) {
118            cellToPtData =
119                vtkSmartPointer<vtkCellDataToPointData>::New();
120            cellToPtData->SetInput(ds);
121            //cellToPtData->PassCellDataOn();
122            cellToPtData->Update();
123            ds = cellToPtData->GetOutput();
124        } else {
125            ERROR("No vector data in dataset %s", _dataSet->getName().c_str());
126            return;
127        }
128    }
129
130    vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
131    if (pd) {
132        // DataSet is a vtkPolyData
133        if (pd->GetNumberOfLines() == 0 &&
134            pd->GetNumberOfPolys() == 0 &&
135            pd->GetNumberOfStrips() == 0) {
136            // DataSet is a point cloud
137            PrincipalPlane plane;
138            double offset;
139            if (_dataSet->is2D(&plane, &offset)) {
140#ifdef MESH_POINT_CLOUDS
141                vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
142                if (plane == PLANE_ZY) {
143                    vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
144                    trans->RotateWXYZ(90, 0, 1, 0);
145                    if (offset != 0.0) {
146                        trans->Translate(-offset, 0, 0);
147                    }
148                    mesher->SetTransform(trans);
149                } else if (plane == PLANE_XZ) {
150                    vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
151                    trans->RotateWXYZ(-90, 1, 0, 0);
152                    if (offset != 0.0) {
153                        trans->Translate(0, -offset, 0);
154                    }
155                    mesher->SetTransform(trans);
156                } else if (offset != 0.0) {
157                    // XY with Z offset
158                    vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
159                    trans->Translate(0, 0, -offset);
160                    mesher->SetTransform(trans);
161                }
162                mesher->SetInput(pd);
163                vtkAlgorithmOutput *warpOutput = initWarp(mesher->GetOutputPort());
164                _dsMapper->SetInputConnection(warpOutput);
165#else
166                vtkSmartPointer<vtkGaussianSplatter> splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
167                vtkSmartPointer<vtkExtractVOI> slicer = vtkSmartPointer<vtkExtractVOI>::New();
168                splatter->SetInput(pd);
169                int dims[3];
170                splatter->GetSampleDimensions(dims);
171                TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
172                if (plane == PLANE_ZY) {
173                    dims[0] = 3;
174                    slicer->SetVOI(1, 1, 0, dims[1]-1, 0, dims[1]-1);
175                } else if (plane == PLANE_XZ) {
176                    dims[1] = 3;
177                    slicer->SetVOI(0, dims[0]-1, 1, 1, 0, dims[2]-1);
178                } else {
179                    dims[2] = 3;
180                    slicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
181                }
182                splatter->SetSampleDimensions(dims);
183                double bounds[6];
184                splatter->Update();
185                splatter->GetModelBounds(bounds);
186                TRACE("Model bounds: %g %g %g %g %g %g",
187                      bounds[0], bounds[1],
188                      bounds[2], bounds[3],
189                      bounds[4], bounds[5]);
190                slicer->SetInputConnection(splatter->GetOutputPort());
191                slicer->SetSampleRate(1, 1, 1);
192                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
193                gf->UseStripsOn();
194                gf->SetInputConnection(slicer->GetOutputPort());
195                vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
196                _dsMapper->SetInputConnection(warpOutput);
197#endif
198            } else {
199                // Data Set is a 3D point cloud
200                vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
201                mesher->SetInput(pd);
202                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
203                gf->SetInputConnection(mesher->GetOutputPort());
204                vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
205                _dsMapper->SetInputConnection(warpOutput);
206             }
207        } else {
208            // DataSet is a vtkPolyData with lines and/or polygons
209            vtkAlgorithmOutput *warpOutput = initWarp(pd);
210            _dsMapper->SetInput(ds);
211            if (warpOutput != NULL) {
212                _dsMapper->SetInputConnection(warpOutput);
213            } else {
214                _dsMapper->SetInput(pd);
215            }
216        }
217    } else {
218        TRACE("Generating surface for data set");
219        // DataSet is NOT a vtkPolyData
220        vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
221        gf->SetInput(ds);
222        vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
223        _dsMapper->SetInputConnection(warpOutput);
224    }
225
226    if (_lut == NULL) {
227        setColorMap(ColorMap::getDefault());
228    }
229
230    setColorMode(_colorMode);
231
232    initProp();
233    getActor()->SetMapper(_dsMapper);
234    _dsMapper->Update();
235}
236
237vtkAlgorithmOutput *Warp::initWarp(vtkAlgorithmOutput *input)
238{
239    TRACE("Warp scale: %g", _warpScale);
240    if (_warpScale == 0.0) {
241        _warp = NULL;
242        return input;
243    } else if (input == NULL) {
244        ERROR("NULL input");
245        return input;
246    } else {
247        if (_warp == NULL)
248            _warp = vtkSmartPointer<vtkWarpVector>::New();
249        _warp->SetScaleFactor(_warpScale);
250        _warp->SetInputConnection(input);
251        return _warp->GetOutputPort();
252    }
253}
254
255vtkAlgorithmOutput *Warp::initWarp(vtkDataSet *ds)
256{
257    TRACE("Warp scale: %g", _warpScale);
258    if (_warpScale == 0.0) {
259        _warp = NULL;
260        return NULL;
261    } else if (ds == NULL) {
262        ERROR("NULL input");
263        return NULL;
264    } else {
265        if (_warp == NULL)
266            _warp = vtkSmartPointer<vtkWarpVector>::New();
267        _warp->SetScaleFactor(_warpScale);
268        _warp->SetInput(ds);
269        return _warp->GetOutputPort();
270    }
271}
272
273/**
274 * \brief Controls relative scaling of deformation
275 */
276void Warp::setWarpScale(double scale)
277{
278    if (_warpScale == scale)
279        return;
280
281    _warpScale = scale;
282    if (_warp == NULL) {
283        vtkAlgorithmOutput *warpOutput = initWarp(_dsMapper->GetInputConnection(0, 0));
284        _dsMapper->SetInputConnection(warpOutput);
285    } else if (scale == 0.0) {
286        vtkAlgorithmOutput *warpInput = _warp->GetInputConnection(0, 0);
287        _dsMapper->SetInputConnection(warpInput);
288        _warp = NULL;
289    } else {
290        _warp->SetScaleFactor(_warpScale);
291    }
292
293    if (_dsMapper != NULL)
294        _dsMapper->Update();
295}
296
297void Warp::updateRanges(Renderer *renderer)
298{
299    if (_dataSet == NULL) {
300        ERROR("called before setDataSet");
301        return;
302    }
303
304    if (renderer->getUseCumulativeRange()) {
305        renderer->getCumulativeDataRange(_dataRange,
306                                         _dataSet->getActiveScalarsName(),
307                                         1);
308        renderer->getCumulativeDataRange(_vectorMagnitudeRange,
309                                         _dataSet->getActiveVectorsName(),
310                                         3);
311        for (int i = 0; i < 3; i++) {
312            renderer->getCumulativeDataRange(_vectorComponentRange[i],
313                                             _dataSet->getActiveVectorsName(),
314                                             3, i);
315        }
316    } else {
317        _dataSet->getScalarRange(_dataRange);
318        _dataSet->getVectorRange(_vectorMagnitudeRange);
319        for (int i = 0; i < 3; i++) {
320            _dataSet->getVectorRange(_vectorComponentRange[i], i);
321        }
322    }
323
324    // Need to update color map ranges
325    double *rangePtr = _colorFieldRange;
326    if (_colorFieldRange[0] > _colorFieldRange[1]) {
327        rangePtr = NULL;
328    }
329    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
330}
331
332void Warp::setColorMode(ColorMode mode)
333{
334    _colorMode = mode;
335    if (_dataSet == NULL)
336        return;
337
338    switch (mode) {
339    case COLOR_BY_SCALAR:
340        setColorMode(mode,
341                     _dataSet->getActiveScalarsType(),
342                     _dataSet->getActiveScalarsName(),
343                     _dataRange);
344        break;
345    case COLOR_BY_VECTOR_MAGNITUDE:
346        setColorMode(mode,
347                     _dataSet->getActiveVectorsType(),
348                     _dataSet->getActiveVectorsName(),
349                     _vectorMagnitudeRange);
350        break;
351    case COLOR_BY_VECTOR_X:
352        setColorMode(mode,
353                     _dataSet->getActiveVectorsType(),
354                     _dataSet->getActiveVectorsName(),
355                     _vectorComponentRange[0]);
356        break;
357    case COLOR_BY_VECTOR_Y:
358        setColorMode(mode,
359                     _dataSet->getActiveVectorsType(),
360                     _dataSet->getActiveVectorsName(),
361                     _vectorComponentRange[1]);
362        break;
363    case COLOR_BY_VECTOR_Z:
364        setColorMode(mode,
365                     _dataSet->getActiveVectorsType(),
366                     _dataSet->getActiveVectorsName(),
367                     _vectorComponentRange[2]);
368        break;
369    case COLOR_CONSTANT:
370    default:
371        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
372        break;
373    }
374}
375
376void Warp::setColorMode(ColorMode mode,
377                        const char *name, double range[2])
378{
379    if (_dataSet == NULL)
380        return;
381    DataSet::DataAttributeType type = DataSet::POINT_DATA;
382    int numComponents = 1;
383    if (name != NULL && strlen(name) > 0 &&
384        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
385        ERROR("Field not found: %s", name);
386        return;
387    }
388    setColorMode(mode, type, name, range);
389}
390
391void Warp::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
392                        const char *name, double range[2])
393{
394    _colorMode = mode;
395    _colorFieldType = type;
396    if (name == NULL)
397        _colorFieldName.clear();
398    else
399        _colorFieldName = name;
400    if (range == NULL) {
401        _colorFieldRange[0] = DBL_MAX;
402        _colorFieldRange[1] = -DBL_MAX;
403    } else {
404        memcpy(_colorFieldRange, range, sizeof(double)*2);
405    }
406
407    if (_dataSet == NULL || _dsMapper == NULL)
408        return;
409
410    switch (type) {
411    case DataSet::POINT_DATA:
412        _dsMapper->SetScalarModeToUsePointFieldData();
413        break;
414    case DataSet::CELL_DATA:
415        _dsMapper->SetScalarModeToUseCellFieldData();
416        break;
417    default:
418        ERROR("Unsupported DataAttributeType: %d", type);
419        return;
420    }
421
422    if (name != NULL && strlen(name) > 0) {
423        _dsMapper->SelectColorArray(name);
424    } else {
425        _dsMapper->SetScalarModeToDefault();
426    }
427
428    if (_lut != NULL) {
429        if (range != NULL) {
430            _lut->SetRange(range);
431        } else if (name != NULL && strlen(name) > 0) {
432            double r[2];
433            int comp = -1;
434            if (mode == COLOR_BY_VECTOR_X)
435                comp = 0;
436            else if (mode == COLOR_BY_VECTOR_Y)
437                comp = 1;
438            else if (mode == COLOR_BY_VECTOR_Z)
439                comp = 2;
440
441            if (_renderer->getUseCumulativeRange()) {
442                int numComponents;
443                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
444                    ERROR("Field not found: %s, type: %d", name, type);
445                    return;
446                } else if (numComponents < comp+1) {
447                    ERROR("Request for component %d in field with %d components",
448                          comp, numComponents);
449                    return;
450                }
451                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
452            } else {
453                _dataSet->getDataRange(r, name, type, comp);
454            }
455            _lut->SetRange(r);
456        }
457    }
458
459    switch (mode) {
460    case COLOR_BY_SCALAR:
461        _dsMapper->ScalarVisibilityOn();
462        break;
463    case COLOR_BY_VECTOR_MAGNITUDE:
464        _dsMapper->ScalarVisibilityOn();
465        if (_lut != NULL) {
466            _lut->SetVectorModeToMagnitude();
467        }
468        break;
469    case COLOR_BY_VECTOR_X:
470        _dsMapper->ScalarVisibilityOn();
471        if (_lut != NULL) {
472            _lut->SetVectorModeToComponent();
473            _lut->SetVectorComponent(0);
474        }
475        break;
476    case COLOR_BY_VECTOR_Y:
477        _dsMapper->ScalarVisibilityOn();
478        if (_lut != NULL) {
479            _lut->SetVectorModeToComponent();
480            _lut->SetVectorComponent(1);
481        }
482        break;
483    case COLOR_BY_VECTOR_Z:
484        _dsMapper->ScalarVisibilityOn();
485        if (_lut != NULL) {
486            _lut->SetVectorModeToComponent();
487            _lut->SetVectorComponent(2);
488        }
489        break;
490    case COLOR_CONSTANT:
491    default:
492        _dsMapper->ScalarVisibilityOff();
493        break;
494    }
495}
496
497/**
498 * \brief Called when the color map has been edited
499 */
500void Warp::updateColorMap()
501{
502    setColorMap(_colorMap);
503}
504
505/**
506 * \brief Associate a colormap lookup table with the DataSet
507 */
508void Warp::setColorMap(ColorMap *cmap)
509{
510    if (cmap == NULL)
511        return;
512
513    _colorMap = cmap;
514 
515    if (_lut == NULL) {
516        _lut = vtkSmartPointer<vtkLookupTable>::New();
517        if (_dsMapper != NULL) {
518            _dsMapper->UseLookupTableScalarRangeOn();
519            _dsMapper->SetLookupTable(_lut);
520        }
521        _lut->DeepCopy(cmap->getLookupTable());
522        switch (_colorMode) {
523        case COLOR_CONSTANT:
524        case COLOR_BY_SCALAR:
525            _lut->SetRange(_dataRange);
526            break;
527        case COLOR_BY_VECTOR_MAGNITUDE:
528            _lut->SetRange(_vectorMagnitudeRange);
529            break;
530        case COLOR_BY_VECTOR_X:
531            _lut->SetRange(_vectorComponentRange[0]);
532            break;
533        case COLOR_BY_VECTOR_Y:
534            _lut->SetRange(_vectorComponentRange[1]);
535            break;
536        case COLOR_BY_VECTOR_Z:
537            _lut->SetRange(_vectorComponentRange[2]);
538            break;
539        default:
540            break;
541        }
542    } else {
543        double range[2];
544        _lut->GetTableRange(range);
545        _lut->DeepCopy(cmap->getLookupTable());
546        _lut->SetRange(range);
547    }
548
549    switch (_colorMode) {
550    case COLOR_BY_VECTOR_MAGNITUDE:
551        _lut->SetVectorModeToMagnitude();
552        break;
553    case COLOR_BY_VECTOR_X:
554        _lut->SetVectorModeToComponent();
555        _lut->SetVectorComponent(0);
556        break;
557    case COLOR_BY_VECTOR_Y:
558        _lut->SetVectorModeToComponent();
559        _lut->SetVectorComponent(1);
560        break;
561    case COLOR_BY_VECTOR_Z:
562        _lut->SetVectorModeToComponent();
563        _lut->SetVectorComponent(2);
564        break;
565    default:
566         break;
567    }
568}
569
570/**
571 * \brief Set a group of world coordinate planes to clip rendering
572 *
573 * Passing NULL for planes will remove all cliping planes
574 */
575void Warp::setClippingPlanes(vtkPlaneCollection *planes)
576{
577    if (_dsMapper != NULL) {
578        _dsMapper->SetClippingPlanes(planes);
579    }
580}
Note: See TracBrowser for help on using the repository browser.