source: trunk/packages/vizservers/vtkvis/RpStreamlines.cpp @ 2513

Last change on this file since 2513 was 2513, checked in by ldelgass, 13 years ago

Streamlines max propagation is _always_ world coordinates -- fix the default
and add a bit of default length. New heuristic is sum of x,y,z bounds. Fix
colormapping when input dataset has only cell data.

  • Property svn:eol-style set to native
File size: 45.9 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2011, Purdue Research Foundation
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cstdlib>
9#include <ctime>
10#include <cfloat>
11#include <cmath>
12
13#include <vtkMath.h>
14#include <vtkActor.h>
15#include <vtkProperty.h>
16#include <vtkPoints.h>
17#include <vtkCellArray.h>
18#include <vtkPolyLine.h>
19#include <vtkRegularPolygonSource.h>
20#include <vtkPointData.h>
21#include <vtkCellData.h>
22#include <vtkCellDataToPointData.h>
23#include <vtkPolygon.h>
24#include <vtkPolyData.h>
25#include <vtkTubeFilter.h>
26#include <vtkRibbonFilter.h>
27#include <vtkTransform.h>
28#include <vtkTransformPolyDataFilter.h>
29#include <vtkVertexGlyphFilter.h>
30
31#include "RpStreamlines.h"
32#include "Trace.h"
33
34using namespace Rappture::VtkVis;
35
36Streamlines::Streamlines() :
37    VtkGraphicsObject(),
38    _lineType(LINES),
39    _colorMode(COLOR_BY_VECTOR_MAGNITUDE),
40    _colorMap(NULL),
41    _seedVisible(true),
42    _dataScale(1)
43{
44    _faceCulling = true;
45    _color[0] = 1.0f;
46    _color[1] = 1.0f;
47    _color[2] = 1.0f;
48    _seedColor[0] = 1.0f;
49    _seedColor[1] = 1.0f;
50    _seedColor[2] = 1.0f;
51    vtkMath::RandomSeed((int)time(NULL));
52    srand((unsigned int)time(NULL));
53}
54
55Streamlines::~Streamlines()
56{
57}
58
59void Streamlines::setDataSet(DataSet *dataSet,
60                             bool useCumulative,
61                             double scalarRange[2],
62                             double vectorMagnitudeRange[2],
63                             double vectorComponentRange[3][2])
64{
65    if (_dataSet != dataSet) {
66        _dataSet = dataSet;
67
68        if (useCumulative) {
69            _dataRange[0] = scalarRange[0];
70            _dataRange[1] = scalarRange[1];
71            _vectorMagnitudeRange[0] = vectorMagnitudeRange[0];
72            _vectorMagnitudeRange[1] = vectorMagnitudeRange[1];
73            for (int i = 0; i < 3; i++) {
74                _vectorComponentRange[i][0] = vectorComponentRange[i][0];
75                _vectorComponentRange[i][1] = vectorComponentRange[i][1];
76            }
77        } else {
78            _dataSet->getScalarRange(_dataRange);
79            _dataSet->getVectorRange(_vectorMagnitudeRange);
80            for (int i = 0; i < 3; i++) {
81                _dataSet->getVectorRange(_vectorComponentRange[i], i);
82            }
83        }
84
85        update();
86    }
87}
88
89/**
90 * \brief Create and initialize a VTK Prop to render Streamlines
91 */
92void Streamlines::initProp()
93{
94    if (_linesActor == NULL) {
95        _linesActor = vtkSmartPointer<vtkActor>::New();
96        _linesActor->GetProperty()->SetColor(_color[0], _color[1], _color[2]);
97        _linesActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
98        _linesActor->GetProperty()->SetLineWidth(_edgeWidth);
99        _linesActor->GetProperty()->SetOpacity(_opacity);
100        _linesActor->GetProperty()->SetAmbient(.2);
101        if (!_lighting)
102            _linesActor->GetProperty()->LightingOff();
103        switch (_lineType) {
104        case LINES:
105            setCulling(_linesActor->GetProperty(), false);
106            _linesActor->GetProperty()->SetRepresentationToWireframe();
107            _linesActor->GetProperty()->EdgeVisibilityOff();
108            break;
109        case TUBES:
110            if (_faceCulling && _opacity == 1.0)
111                setCulling(_linesActor->GetProperty(), true);
112            _linesActor->GetProperty()->SetRepresentationToSurface();
113            _linesActor->GetProperty()->EdgeVisibilityOff();
114            break;
115        case RIBBONS:
116            setCulling(_linesActor->GetProperty(), false);
117            _linesActor->GetProperty()->SetRepresentationToSurface();
118            _linesActor->GetProperty()->EdgeVisibilityOff();
119            break;
120        default:
121            ;
122        }
123    }
124    if (_seedActor == NULL) {
125        _seedActor = vtkSmartPointer<vtkActor>::New();
126        _seedActor->GetProperty()->SetColor(_seedColor[0], _seedColor[1], _seedColor[2]);
127        _seedActor->GetProperty()->SetEdgeColor(_seedColor[0], _seedColor[1], _seedColor[2]);
128        _seedActor->GetProperty()->SetLineWidth(1);
129        _seedActor->GetProperty()->SetPointSize(2);
130        _seedActor->GetProperty()->SetOpacity(_opacity);
131        _seedActor->GetProperty()->SetRepresentationToWireframe();
132        _seedActor->GetProperty()->LightingOff();
133        setSeedVisibility(_seedVisible);
134    }
135    if (_prop == NULL) {
136        _prop = vtkSmartPointer<vtkAssembly>::New();
137        getAssembly()->AddPart(_linesActor);
138        getAssembly()->AddPart(_seedActor);
139    }
140}
141
142/**
143 * \brief Get a pseudo-random number in range [min,max]
144 */
145double Streamlines::getRandomNum(double min, double max)
146{
147#if 1
148    return vtkMath::Random(min, max);
149#else
150    int r = rand();
151    return (min + ((double)r / RAND_MAX) * (max - min));
152#endif
153}
154
155/**
156 * \brief Get a random 3D point within an AABB
157 *
158 * \param[out] pt The random point
159 * \param[in] bounds The bounds of the AABB
160 */
161void Streamlines::getRandomPoint(double pt[3], const double bounds[6])
162{
163    pt[0] = getRandomNum(bounds[0], bounds[1]);
164    pt[1] = getRandomNum(bounds[2], bounds[3]);
165    pt[2] = getRandomNum(bounds[4], bounds[5]);
166}
167
168/**
169 * \brief Get a random point within a triangle (including edges)
170 *
171 * \param[out] pt The random point
172 * \param[in] v1 Triangle vertex 1
173 * \param[in] v2 Triangle vertex 2
174 * \param[in] v3 Triangle vertex 3
175 */
176void Streamlines::getRandomPointInTriangle(double pt[3],
177                                           const double v0[3],
178                                           const double v1[3],
179                                           const double v2[3])
180{
181    // Choose random barycentric coordinates
182    double bary[3];
183    bary[0] = getRandomNum(0, 1);
184    bary[1] = getRandomNum(0, 1);
185    if (bary[0] + bary[1] > 1.0) {
186        bary[0] = 1.0 - bary[0];
187        bary[1] = 1.0 - bary[1];
188    }
189    bary[2] = 1.0 - bary[0] - bary[1];
190
191    TRACE("bary %g %g %g", bary[0], bary[1], bary[2]);
192    // Convert to cartesian coords
193    for (int i = 0; i < 3; i++) {
194        pt[i] = v0[i] * bary[0] + v1[i] * bary[1] + v2[i] * bary[2];
195    }
196}
197
198void Streamlines::getRandomPointInTetrahedron(double pt[3],
199                                              const double v0[3],
200                                              const double v1[3],
201                                              const double v2[3],
202                                              const double v3[3])
203{
204    // Choose random barycentric coordinates
205    double bary[4];
206    bary[0] = getRandomNum(0, 1);
207    bary[1] = getRandomNum(0, 1);
208    bary[2] = getRandomNum(0, 1);
209    if (bary[0] + bary[1] > 1.0) {
210        bary[0] = 1.0 - bary[0];
211        bary[1] = 1.0 - bary[1];
212    }
213    if (bary[1] + bary[2] > 1.0) {
214        double tmp = bary[2];
215        bary[2] = 1.0 - bary[0] - bary[1];
216        bary[1] = 1.0 - tmp;
217    } else if (bary[0] + bary[1] + bary[2] > 1.0) {
218        double tmp = bary[2];
219        bary[2] = bary[0] + bary[1] + bary[2] - 1.0;
220        bary[0] = 1.0 - bary[1] - tmp;
221    }
222    bary[3] = 1.0 - bary[0] - bary[1] - bary[2];
223    TRACE("bary %g %g %g %g", bary[0], bary[1], bary[2], bary[3]);
224    // Convert to cartesian coords
225    for (int i = 0; i < 3; i++) {
226#if 0
227        pt[i] = (v0[i] - v3[i]) * bary[0] +
228                (v1[i] - v3[i]) * bary[1] +
229                (v2[i] - v3[i]) * bary[2] + v3[i];
230#else
231        pt[i] = v0[i] * bary[0] + v1[i] * bary[1] +
232                v2[i] * bary[2] + v3[i] * bary[3];
233#endif
234    }
235}
236
237/**
238 * \brief Get a random point on a line segment (including endpoints)
239 */
240void Streamlines::getRandomPointOnLineSegment(double pt[3],
241                                              const double endpt[3],
242                                              const double endpt2[3])
243{
244    double ratio = getRandomNum(0, 1);
245    pt[0] = endpt[0] + ratio * (endpt2[0] - endpt[0]);
246    pt[1] = endpt[1] + ratio * (endpt2[1] - endpt[1]);
247    pt[2] = endpt[2] + ratio * (endpt2[2] - endpt[2]);
248}
249
250/**
251 * \brief Get a random point within a vtkDataSet's mesh
252 *
253 * Note: This currently doesn't always give a uniform distribution
254 * of points in space and can generate points outside the mesh for
255 * unusual cell types
256 */
257void Streamlines::getRandomCellPt(double pt[3], vtkDataSet *ds)
258{
259    int numCells = (int)ds->GetNumberOfCells();
260    // XXX: Not uniform distribution (shouldn't use mod, and assumes
261    // all cells are equal area/volume)
262    int cell = rand() % numCells;
263    int type = ds->GetCellType(cell);
264    if (type == VTK_VERTEX) {
265        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
266        ds->GetCellPoints(cell, ptIds);
267        assert(ptIds->GetNumberOfIds() == 1);
268        ds->GetPoint(ptIds->GetId(0), pt);
269    } else if (type == VTK_POLY_VERTEX) {
270        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
271        ds->GetCellPoints(cell, ptIds);
272        assert(ptIds->GetNumberOfIds() >= 1);
273        int id = rand() % ptIds->GetNumberOfIds();
274        ds->GetPoint(ptIds->GetId(id), pt);
275    } else if (type == VTK_LINE) {
276        double v[2][3];
277        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
278        ds->GetCellPoints(cell, ptIds);
279        assert(ptIds->GetNumberOfIds() == 2);
280        for (int i = 0; i < 2; i++) {
281            ds->GetPoint(ptIds->GetId(i), v[i]);
282        }
283        getRandomPointOnLineSegment(pt, v[0], v[1]);
284    } else if (type == VTK_POLY_LINE) {
285        double v[2][3];
286        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
287        ds->GetCellPoints(cell, ptIds);
288        assert(ptIds->GetNumberOfIds() >= 2);
289        int id = rand() % (ptIds->GetNumberOfIds()-1);
290        for (int i = 0; i < 2; i++) {
291            ds->GetPoint(ptIds->GetId(id+i), v[i]);
292        }
293        getRandomPointOnLineSegment(pt, v[0], v[1]);
294    } else if (type == VTK_TRIANGLE) {
295        double v[3][3];
296        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
297        ds->GetCellPoints(cell, ptIds);
298        assert(ptIds->GetNumberOfIds() == 3);
299        for (int i = 0; i < 3; i++) {
300            ds->GetPoint(ptIds->GetId(i), v[i]);
301        }
302        getRandomPointInTriangle(pt, v[0], v[1], v[2]);
303    } else if (type == VTK_TRIANGLE_STRIP) {
304        double v[3][3];
305        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
306        ds->GetCellPoints(cell, ptIds);
307        assert(ptIds->GetNumberOfIds() >= 3);
308        int id = rand() % (ptIds->GetNumberOfIds()-2);
309        for (int i = 0; i < 3; i++) {
310            ds->GetPoint(ptIds->GetId(id+i), v[i]);
311        }
312        getRandomPointInTriangle(pt, v[0], v[1], v[2]);
313    } else if (type == VTK_POLYGON) {
314        vtkPolygon *poly = vtkPolygon::SafeDownCast(ds->GetCell(cell));
315        assert (poly != NULL);
316        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
317        poly->Triangulate(ptIds);
318        assert(ptIds->GetNumberOfIds() >= 3 && ptIds->GetNumberOfIds() % 3 == 0);
319        int tri = rand() % (ptIds->GetNumberOfIds()/3);
320        double v[3][3];
321        for (int i = 0; i < 3; i++) {
322            ds->GetPoint(ptIds->GetId(i + tri * 3), v[i]);
323        }
324        getRandomPointInTriangle(pt, v[0], v[1], v[2]);
325    } else if (type == VTK_QUAD) {
326        double v[4][3];
327        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
328        ds->GetCellPoints(cell, ptIds);
329        assert(ptIds->GetNumberOfIds() == 4);
330        for (int i = 0; i < 4; i++) {
331            ds->GetPoint(ptIds->GetId(i), v[i]);
332        }
333        int tri = rand() & 0x1;
334        if (tri) {
335            getRandomPointInTriangle(pt, v[0], v[1], v[2]);
336        } else {
337            getRandomPointInTriangle(pt, v[0], v[2], v[3]);
338        }
339    } else if (type == VTK_TETRA) {
340        double v[4][3];
341        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
342        ds->GetCellPoints(cell, ptIds);
343        assert(ptIds->GetNumberOfIds() == 4);
344        for (int i = 0; i < 4; i++) {
345            ds->GetPoint(ptIds->GetId(i), v[i]);
346        }
347        getRandomPointInTetrahedron(pt, v[0], v[1], v[2], v[3]);
348    } else if (type == VTK_HEXAHEDRON) {
349        double v[8][3];
350        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
351        ds->GetCellPoints(cell, ptIds);
352        assert(ptIds->GetNumberOfIds() == 8);
353        for (int i = 0; i < 8; i++) {
354            ds->GetPoint(ptIds->GetId(i), v[i]);
355        }
356        int tetra = rand() % 5;
357        switch (tetra) {
358        case 0:
359            getRandomPointInTetrahedron(pt, v[0], v[1], v[2], v[5]);
360            break;
361        case 1:
362            getRandomPointInTetrahedron(pt, v[0], v[2], v[7], v[3]);
363            break;
364        case 2:
365            getRandomPointInTetrahedron(pt, v[0], v[5], v[7], v[4]);
366            break;
367        case 3:
368            getRandomPointInTetrahedron(pt, v[5], v[2], v[7], v[6]);
369            break;
370        case 4:
371        default:
372            getRandomPointInTetrahedron(pt, v[0], v[2], v[7], v[5]);
373            break;
374        }
375    } else if (type == VTK_WEDGE) {
376        double v[6][3];
377        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
378        ds->GetCellPoints(cell, ptIds);
379        assert(ptIds->GetNumberOfIds() == 6);
380        for (int i = 0; i < 6; i++) {
381            ds->GetPoint(ptIds->GetId(i), v[i]);
382        }
383        double vv[3][3];
384        getRandomPointOnLineSegment(vv[0], v[0], v[3]);
385        getRandomPointOnLineSegment(vv[1], v[1], v[4]);
386        getRandomPointOnLineSegment(vv[2], v[2], v[5]);
387        getRandomPointInTriangle(pt, vv[0], vv[1], vv[2]);
388    } else if (type == VTK_PYRAMID) {
389        double v[5][3];
390        vtkSmartPointer<vtkIdList> ptIds = vtkSmartPointer<vtkIdList>::New();
391        ds->GetCellPoints(cell, ptIds);
392        assert(ptIds->GetNumberOfIds() == 5);
393        for (int i = 0; i < 5; i++) {
394            ds->GetPoint(ptIds->GetId(i), v[i]);
395        }
396        int tetra = rand() & 0x1;
397        if (tetra) {
398            getRandomPointInTetrahedron(pt, v[0], v[1], v[2], v[4]);
399        } else {
400            getRandomPointInTetrahedron(pt, v[0], v[2], v[3], v[4]);
401        }
402    } else {
403        double bounds[6];
404        ds->GetCellBounds(cell, bounds);
405        // Note: For pixel/voxel cells, this is exact.  However, if the cell is
406        // not an axis aligned box, the point may be outside the cell
407        getRandomPoint(pt, bounds);
408    }
409}
410
411/**
412 * \brief Internal method to set up pipeline after a state change
413 */
414void Streamlines::update()
415{
416    if (_dataSet == NULL) {
417        return;
418    }
419
420    vtkDataSet *ds = _dataSet->getVtkDataSet();
421
422    double bounds[6];
423    _dataSet->getBounds(bounds);
424    double xLen = bounds[1] - bounds[0];
425    double yLen = bounds[3] - bounds[2];
426    double zLen = bounds[5] - bounds[4];
427    double maxBound = 0.0;
428    if (xLen > maxBound) {
429        maxBound = xLen;
430    }
431    if (yLen > maxBound) {
432        maxBound = yLen;
433    }
434    if (zLen > maxBound) {
435        maxBound = zLen;
436    }
437
438    double cellSizeRange[2];
439    double avgSize;
440    _dataSet->getCellSizeRange(cellSizeRange, &avgSize);
441    _dataScale = avgSize / 8.;
442
443    vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
444
445    if (ds->GetPointData() == NULL ||
446        ds->GetPointData()->GetVectors() == NULL) {
447        TRACE("No vector point data found in DataSet %s", _dataSet->getName().c_str());
448        if (ds->GetCellData() == NULL ||
449            ds->GetCellData()->GetVectors() == NULL) {
450            ERROR("No vectors found in DataSet %s", _dataSet->getName().c_str());
451        } else {
452            cellToPtData =
453                vtkSmartPointer<vtkCellDataToPointData>::New();
454            cellToPtData->SetInput(ds);
455            //cellToPtData->PassCellDataOn();
456            cellToPtData->Update();
457            ds = cellToPtData->GetOutput();
458        }
459    }
460
461    if (_streamTracer == NULL) {
462        _streamTracer = vtkSmartPointer<vtkStreamTracer>::New();
463    }
464
465    _streamTracer->SetInput(ds);
466    _streamTracer->SetMaximumPropagation(xLen + yLen + zLen);
467    _streamTracer->SetIntegratorTypeToRungeKutta45();
468
469    TRACE("Setting streamlines max length to %g",
470          _streamTracer->GetMaximumPropagation());
471
472    if (_pdMapper == NULL) {
473        _pdMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
474        _pdMapper->SetResolveCoincidentTopologyToPolygonOffset();
475        _pdMapper->ScalarVisibilityOn();
476    }
477    if (_seedMapper == NULL) {
478        _seedMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
479        _seedMapper->SetResolveCoincidentTopologyToPolygonOffset();
480        _seedMapper->ScalarVisibilityOff();
481    }
482
483    // Set up seed source object
484    setSeedToFilledMesh(200);
485
486    switch (_lineType) {
487    case LINES: {
488        _streamTracer->SetComputeVorticity(false);
489        _pdMapper->SetInputConnection(_streamTracer->GetOutputPort());
490    }
491        break;
492    case TUBES: {
493        _streamTracer->SetComputeVorticity(true);
494        _lineFilter = vtkSmartPointer<vtkTubeFilter>::New();
495        vtkTubeFilter *tubeFilter = vtkTubeFilter::SafeDownCast(_lineFilter);
496        tubeFilter->SetNumberOfSides(5);
497        _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
498        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
499    }
500        break;
501    case RIBBONS: {
502        _streamTracer->SetComputeVorticity(true);
503        _lineFilter = vtkSmartPointer<vtkRibbonFilter>::New();
504        _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
505        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
506    }
507        break;
508    default:
509        ERROR("Unknown LineType: %d", _lineType);
510    }
511
512#if defined(DEBUG) && defined(WANT_TRACE)
513    _streamTracer->Update();
514    vtkPolyData *pd = _streamTracer->GetOutput();
515    TRACE("Verts: %d Lines: %d Polys: %d Strips: %d",
516                  pd->GetNumberOfVerts(),
517                  pd->GetNumberOfLines(),
518                  pd->GetNumberOfPolys(),
519                  pd->GetNumberOfStrips());
520#endif
521
522    initProp();
523
524    _seedActor->SetMapper(_seedMapper);
525
526    if (_lut == NULL) {
527        setColorMap(ColorMap::getDefault());
528    }
529
530    setColorMode(_colorMode);
531
532    _linesActor->SetMapper(_pdMapper);
533    _pdMapper->Update();
534    _seedMapper->Update();
535}
536
537/**
538 * \brief Use points of the DataSet associated with this
539 * Streamlines as seeds
540 */
541void Streamlines::setSeedToMeshPoints()
542{
543    setSeedToMeshPoints(_dataSet->getVtkDataSet());
544}
545
546/**
547 * \brief Use seed points randomly distributed within the cells
548 * of the DataSet associated with this Streamlines
549 *
550 * Note: The current implementation doesn't give a uniform
551 * distribution of points, and points outside the mesh bounds
552 * may be generated
553 *
554 * \param[in] numPoints Number of random seed points to generate
555 */
556void Streamlines::setSeedToFilledMesh(int numPoints)
557{
558    setSeedToFilledMesh(_dataSet->getVtkDataSet(), numPoints);
559}
560
561/**
562 * \brief Use points of a supplied vtkDataSet as seeds
563 *
564 * \param[in] seed vtkDataSet with points to use as seeds
565 */
566void Streamlines::setSeedToMeshPoints(vtkDataSet *seed)
567{
568    if (_streamTracer != NULL) {
569        TRACE("Seed points: %d", seed->GetNumberOfPoints());
570        vtkSmartPointer<vtkDataSet> oldSeed;
571        if (_streamTracer->GetSource() != NULL) {
572            oldSeed = _streamTracer->GetSource();
573        }
574
575        _streamTracer->SetSource(seed);
576        if (oldSeed != NULL) {
577            oldSeed->SetPipelineInformation(NULL);
578        }
579
580        if (vtkPolyData::SafeDownCast(seed) != NULL) {
581            _seedMapper->SetInput(vtkPolyData::SafeDownCast(seed));
582        } else {
583            vtkSmartPointer<vtkVertexGlyphFilter> vertFilter = vtkSmartPointer<vtkVertexGlyphFilter>::New();
584            vertFilter->SetInput(seed);
585            _seedMapper->SetInputConnection(vertFilter->GetOutputPort());
586        }
587    }
588}
589
590/**
591 * \brief Use seed points randomly distributed within the cells
592 * of a supplied vtkDataSet
593 *
594 * Note: The current implementation doesn't give a uniform
595 * distribution of points, and points outside the mesh bounds
596 * may be generated
597 *
598 * \param[in] ds vtkDataSet containing cells
599 * \param[in] numPoints Number of random seed points to generate
600 */
601void Streamlines::setSeedToFilledMesh(vtkDataSet *ds, int numPoints)
602{
603    if (_streamTracer != NULL) {
604        // Set up seed source object
605        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
606        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
607        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
608
609        if (ds->GetNumberOfCells() < 1) {
610            ERROR("No cells in mesh");
611        }
612
613        for (int i = 0; i < numPoints; i++) {
614            double pt[3];
615            getRandomCellPt(pt, ds);
616            //TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
617            pts->InsertNextPoint(pt);
618            cells->InsertNextCell(1);
619            cells->InsertCellPoint(i);
620        }
621
622        seed->SetPoints(pts);
623        seed->SetVerts(cells);
624
625        TRACE("Seed points: %d", seed->GetNumberOfPoints());
626        vtkSmartPointer<vtkDataSet> oldSeed;
627        if (_streamTracer->GetSource() != NULL) {
628            oldSeed = _streamTracer->GetSource();
629        }
630
631        _streamTracer->SetSource(seed);
632        if (oldSeed != NULL) {
633            oldSeed->SetPipelineInformation(NULL);
634        }
635
636        _seedMapper->SetInput(seed);
637    }
638}
639
640/**
641 * \brief Use seed points along a line
642 *
643 * \param[in] start Starting point of rake line
644 * \param[in] end End point of rake line
645 * \param[in] numPoints Number of points along line to generate
646 */
647void Streamlines::setSeedToRake(double start[3], double end[3], int numPoints)
648{
649    if (numPoints < 2)
650        return;
651    if (_streamTracer != NULL) {
652        // Set up seed source object
653        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
654        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
655        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
656        vtkSmartPointer<vtkPolyLine> polyline = vtkSmartPointer<vtkPolyLine>::New();
657
658        double dir[3];
659        for (int i = 0; i < 3; i++) {
660            dir[i] = end[i] - start[i];
661        }
662
663        polyline->GetPointIds()->SetNumberOfIds(numPoints);
664        for (int i = 0; i < numPoints; i++) {
665            double pt[3];
666            for (int ii = 0; ii < 3; ii++) {
667                pt[ii] = start[ii] + dir[ii] * ((double)i / (numPoints-1));
668            }
669            //TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
670            pts->InsertNextPoint(pt);
671            polyline->GetPointIds()->SetId(i, i);
672        }
673
674        cells->InsertNextCell(polyline);
675        seed->SetPoints(pts);
676        seed->SetLines(cells);
677
678        TRACE("Seed points: %d", seed->GetNumberOfPoints());
679        vtkSmartPointer<vtkDataSet> oldSeed;
680        if (_streamTracer->GetSource() != NULL) {
681            oldSeed = _streamTracer->GetSource();
682        }
683
684        _streamTracer->SetSource(seed);
685        if (oldSeed != NULL) {
686            oldSeed->SetPipelineInformation(NULL);
687        }
688
689        _seedMapper->SetInput(seed);
690    }
691}
692
693/**
694 * \brief Create seed points inside a disk with an optional hole
695 *
696 * \param[in] center Center point of disk
697 * \param[in] normal Normal vector to orient disk
698 * \param[in] radius Radius of disk
699 * \param[in] innerRadius Radius of hole at center of disk
700 * \param[in] numPoints Number of random points to generate
701 */
702void Streamlines::setSeedToDisk(double center[3],
703                                double normal[3],
704                                double radius,
705                                double innerRadius,
706                                int numPoints)
707{
708    if (_streamTracer != NULL) {
709        // Set up seed source object
710        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
711        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
712        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
713
714        // The following code is based on vtkRegularPolygonSource::RequestData
715
716        double px[3];
717        double py[3];
718        double axis[3] = {1., 0., 0.};
719
720        if (vtkMath::Normalize(normal) == 0.0) {
721            normal[0] = 0.0;
722            normal[1] = 0.0;
723            normal[2] = 1.0;
724        }
725
726        // Find axis in plane (orthogonal to normal)
727        bool done = false;
728        vtkMath::Cross(normal, axis, px);
729        if (vtkMath::Normalize(px) > 1.0e-3) {
730            done = true;
731        }
732        if (!done) {
733            axis[0] = 0.0;
734            axis[1] = 1.0;
735            axis[2] = 0.0;
736            vtkMath::Cross(normal, axis, px);
737            if (vtkMath::Normalize(px) > 1.0e-3) {
738                done = true;
739            }
740        }
741        if (!done) {
742            axis[0] = 0.0;
743            axis[1] = 0.0;
744            axis[2] = 1.0;
745            vtkMath::Cross(normal, axis, px);
746            vtkMath::Normalize(px);
747        }
748        // Create third orthogonal basis vector
749        vtkMath::Cross(px, normal, py);
750
751        double minSquared = (innerRadius*innerRadius)/(radius*radius);
752        for (int j = 0; j < numPoints; j++) {
753            // Get random sweep angle and radius
754            double angle = getRandomNum(0, 2.0 * vtkMath::DoublePi());
755            // Need sqrt to get uniform distribution
756            double r = sqrt(getRandomNum(minSquared, 1)) * radius;
757            double pt[3];
758            for (int i = 0; i < 3; i++) {
759                pt[i] = center[i] + r * (px[i] * cos(angle) + py[i] * sin(angle));
760            }
761            //TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
762            pts->InsertNextPoint(pt);
763            cells->InsertNextCell(1);
764            cells->InsertCellPoint(j);
765        }
766
767        seed->SetPoints(pts);
768        seed->SetVerts(cells);
769
770        TRACE("Seed points: %d", seed->GetNumberOfPoints());
771        vtkSmartPointer<vtkDataSet> oldSeed;
772        if (_streamTracer->GetSource() != NULL) {
773            oldSeed = _streamTracer->GetSource();
774        }
775
776        _streamTracer->SetSource(seed);
777        if (oldSeed != NULL) {
778            oldSeed->SetPipelineInformation(NULL);
779        }
780
781        _seedMapper->SetInput(seed);
782    }
783}
784
785/**
786 * \brief Use seed points from an n-sided polygon
787 *
788 * \param[in] center Center point of polygon
789 * \param[in] normal Normal vector to orient polygon
790 * \param[in] angle Angle in degrees to rotate about normal
791 * \param[in] radius Radius of circumscribing circle
792 * \param[in] numSides Number of polygon sides (and points) to generate
793 */
794void Streamlines::setSeedToPolygon(double center[3],
795                                   double normal[3],
796                                   double angle,
797                                   double radius,
798                                   int numSides)
799{
800    if (_streamTracer != NULL) {
801        // Set up seed source object
802        vtkSmartPointer<vtkRegularPolygonSource> seed = vtkSmartPointer<vtkRegularPolygonSource>::New();
803
804        seed->SetCenter(center);
805        seed->SetNormal(normal);
806        seed->SetRadius(radius);
807        seed->SetNumberOfSides(numSides);
808        seed->GeneratePolygonOn();
809
810        if (angle != 0.0) {
811            vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
812            trans->RotateWXYZ(angle, normal);
813            vtkSmartPointer<vtkTransformPolyDataFilter> transFilt =
814                vtkSmartPointer<vtkTransformPolyDataFilter>::New();
815            transFilt->SetInputConnection(seed->GetOutputPort());
816            transFilt->SetTransform(trans);
817        }
818
819        TRACE("Seed points: %d", numSides);
820        vtkSmartPointer<vtkDataSet> oldSeed;
821        if (_streamTracer->GetSource() != NULL) {
822            oldSeed = _streamTracer->GetSource();
823        }
824
825        if (angle != 0.0) {
826            vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
827            trans->Translate(+center[0], +center[1], +center[2]);
828            trans->RotateWXYZ(angle, normal);
829            trans->Translate(-center[0], -center[1], -center[2]);
830            vtkSmartPointer<vtkTransformPolyDataFilter> transFilt =
831                vtkSmartPointer<vtkTransformPolyDataFilter>::New();
832            transFilt->SetInputConnection(seed->GetOutputPort());
833            transFilt->SetTransform(trans);
834            _streamTracer->SetSourceConnection(transFilt->GetOutputPort());
835            _seedMapper->SetInputConnection(transFilt->GetOutputPort());
836        } else {
837            _streamTracer->SetSourceConnection(seed->GetOutputPort());
838            _seedMapper->SetInputConnection(seed->GetOutputPort());
839        }
840
841        if (oldSeed != NULL) {
842            oldSeed->SetPipelineInformation(NULL);
843        }
844    }
845}
846
847/**
848 * \brief Use seed points from an n-sided polygon
849 *
850 * \param[in] center Center point of polygon
851 * \param[in] normal Normal vector to orient polygon
852 * \param[in] angle Angle in degrees to rotate about normal
853 * \param[in] radius Radius of circumscribing circle
854 * \param[in] numSides Number of polygon sides (and points) to generate
855 * \param[in] numPoints Number of random points to generate
856 */
857void Streamlines::setSeedToFilledPolygon(double center[3],
858                                         double normal[3],
859                                         double angle,
860                                         double radius,
861                                         int numSides,
862                                         int numPoints)
863{
864    if (_streamTracer != NULL) {
865         // Set up seed source object
866        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
867        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
868        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
869
870        // The following code is based on vtkRegularPolygonSource::RequestData
871
872        double px[3];
873        double py[3];
874        double axis[3] = {1., 0., 0.};
875
876        if (vtkMath::Normalize(normal) == 0.0) {
877            normal[0] = 0.0;
878            normal[1] = 0.0;
879            normal[2] = 1.0;
880        }
881
882        // Find axis in plane (orthogonal to normal)
883        bool done = false;
884        vtkMath::Cross(normal, axis, px);
885        if (vtkMath::Normalize(px) > 1.0e-3) {
886            done = true;
887        }
888        if (!done) {
889            axis[0] = 0.0;
890            axis[1] = 1.0;
891            axis[2] = 0.0;
892            vtkMath::Cross(normal, axis, px);
893            if (vtkMath::Normalize(px) > 1.0e-3) {
894                done = true;
895            }
896        }
897        if (!done) {
898            axis[0] = 0.0;
899            axis[1] = 0.0;
900            axis[2] = 1.0;
901            vtkMath::Cross(normal, axis, px);
902            vtkMath::Normalize(px);
903        }
904        // Create third orthogonal basis vector
905        vtkMath::Cross(px, normal, py);
906
907        double verts[numSides][3];
908        double sliceTheta = 2.0 * vtkMath::DoublePi() / (double)numSides;
909        angle = vtkMath::RadiansFromDegrees(angle);
910        for (int j = 0; j < numSides; j++) {
911            for (int i = 0; i < 3; i++) {
912                double theta = sliceTheta * (double)j - angle;
913                verts[j][i] = center[i] + radius * (px[i] * cos(theta) +
914                                                    py[i] * sin(theta));
915            }
916            //TRACE("Vert %d: %g %g %g", j, verts[j][0], verts[j][1], verts[j][2]);
917        }
918
919        // Note: this gives a uniform distribution because the polygon is regular and
920        // the triangular sections have equal area
921        if (numSides == 3) {
922            for (int j = 0; j < numPoints; j++) {
923                double pt[3];
924                getRandomPointInTriangle(pt, verts[0], verts[1], verts[2]);
925                //TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
926                pts->InsertNextPoint(pt);
927                cells->InsertNextCell(1);
928                cells->InsertCellPoint(j);
929            }
930        } else {
931            for (int j = 0; j < numPoints; j++) {
932                // Get random triangle section
933                int tri = rand() % numSides;
934                double pt[3];
935                getRandomPointInTriangle(pt, center, verts[tri], verts[(tri+1) % numSides]);
936                //TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
937                pts->InsertNextPoint(pt);
938                cells->InsertNextCell(1);
939                cells->InsertCellPoint(j);
940            }
941        }
942
943        seed->SetPoints(pts);
944        seed->SetVerts(cells);
945
946        TRACE("Seed points: %d", seed->GetNumberOfPoints());
947        vtkSmartPointer<vtkDataSet> oldSeed;
948        if (_streamTracer->GetSource() != NULL) {
949            oldSeed = _streamTracer->GetSource();
950        }
951
952        _streamTracer->SetSource(seed);
953        if (oldSeed != NULL) {
954            oldSeed->SetPipelineInformation(NULL);
955        }
956
957        _seedMapper->SetInput(seed);
958    }
959}
960
961/**
962 * \brief Set the integration method used
963 */
964void Streamlines::setIntegrator(IntegratorType integrator)
965{
966    if (_streamTracer != NULL) {
967        switch (integrator) {
968        case RUNGE_KUTTA2:
969            _streamTracer->SetIntegratorTypeToRungeKutta2();
970            break;
971        case RUNGE_KUTTA4:
972            _streamTracer->SetIntegratorTypeToRungeKutta4();
973            break;
974        case RUNGE_KUTTA45:
975            _streamTracer->SetIntegratorTypeToRungeKutta45();
976            break;
977        default:
978            ;
979        }
980    }
981}
982
983/**
984 * \brief Set the direction of integration
985 */
986void Streamlines::setIntegrationDirection(IntegrationDirection dir)
987{
988    if (_streamTracer != NULL) {
989        switch (dir) {
990        case FORWARD:
991            _streamTracer->SetIntegrationDirectionToForward();
992            break;
993        case BACKWARD:
994            _streamTracer->SetIntegrationDirectionToBackward();
995            break;
996        case BOTH:
997            _streamTracer->SetIntegrationDirectionToBoth();
998            break;
999        default:
1000            ;
1001        }
1002    }
1003}
1004
1005/**
1006 * \brief Set the step size units.  Length units are world
1007 * coordinates, and cell units means steps are from cell to
1008 * cell.  Default is cell units.
1009 *
1010 * Note: calling this function will not convert existing
1011 * initial, minimum or maximum step value settings to the
1012 * new units, so this function should be called before
1013 * setting step values.
1014 */
1015void Streamlines::setIntegrationStepUnit(StepUnit unit)
1016{
1017    if (_streamTracer != NULL) {
1018        switch (unit) {
1019        case LENGTH_UNIT:
1020            _streamTracer->SetIntegrationStepUnit(vtkStreamTracer::LENGTH_UNIT);
1021            break;
1022        case CELL_LENGTH_UNIT:
1023            _streamTracer->SetIntegrationStepUnit(vtkStreamTracer::CELL_LENGTH_UNIT);
1024            break;
1025        default:
1026            ;
1027        }
1028    }
1029}
1030
1031/**
1032 * \brief Set initial step size for adaptive step integrator in
1033 * step units (see setIntegrationStepUnit).  For non-adaptive
1034 * integrators, this is the fixed step size.
1035 */
1036void Streamlines::setInitialIntegrationStep(double step)
1037{
1038    if (_streamTracer != NULL) {
1039        _streamTracer->SetInitialIntegrationStep(step);
1040    }
1041}
1042
1043/**
1044 * \brief Set minimum step for adaptive step integrator in
1045 * step units (see setIntegrationStepUnit)
1046 */
1047void Streamlines::setMinimumIntegrationStep(double step)
1048{
1049    if (_streamTracer != NULL) {
1050        _streamTracer->SetMinimumIntegrationStep(step);
1051    }
1052}
1053
1054/**
1055 * \brief Set maximum step for adaptive step integrator in
1056 * step units (see setIntegrationStepUnit)
1057 */
1058void Streamlines::setMaximumIntegrationStep(double step)
1059{
1060    if (_streamTracer != NULL) {
1061        _streamTracer->SetMaximumIntegrationStep(step);
1062    }
1063}
1064
1065/**
1066 * \brief Set maximum error tolerance
1067 */
1068void Streamlines::setMaximumError(double error)
1069{
1070    if (_streamTracer != NULL) {
1071        _streamTracer->SetMaximumError(error);
1072    }
1073}
1074
1075/**
1076 * \brief Set maximum length of stream lines in world
1077 * coordinates
1078 */
1079void Streamlines::setMaxPropagation(double length)
1080{
1081    if (_streamTracer != NULL) {
1082        _streamTracer->SetMaximumPropagation(length);
1083    }
1084}
1085
1086/**
1087 * \brief Set maximum number of integration steps
1088 */
1089void Streamlines::setMaxNumberOfSteps(int steps)
1090{
1091    if (_streamTracer != NULL) {
1092        _streamTracer->SetMaximumNumberOfSteps(steps);
1093    }
1094}
1095
1096/**
1097 * \brief Set the minimum speed before integration stops
1098 */
1099void Streamlines::setTerminalSpeed(double speed)
1100{
1101    if (_streamTracer != NULL) {
1102        _streamTracer->SetTerminalSpeed(speed);
1103    }
1104}
1105
1106/**
1107 * \brief Set streamline type to polylines
1108 */
1109void Streamlines::setLineTypeToLines()
1110{
1111    _lineType = LINES;
1112    if (_streamTracer != NULL &&
1113        _pdMapper != NULL) {
1114        _streamTracer->SetComputeVorticity(false);
1115        _pdMapper->SetInputConnection(_streamTracer->GetOutputPort());
1116        _lineFilter = NULL;
1117        setCulling(_linesActor->GetProperty(), false);
1118        _linesActor->GetProperty()->SetRepresentationToWireframe();
1119        _linesActor->GetProperty()->LightingOff();
1120    }
1121}
1122
1123/**
1124 * \brief Set streamline type to 3D tubes
1125 *
1126 * \param[in] numSides Number of sides (>=3) for tubes
1127 * \param[in] radius World coordinate minimum tube radius
1128 */
1129void Streamlines::setLineTypeToTubes(int numSides, double radius)
1130{
1131    _lineType = TUBES;
1132    if (_streamTracer != NULL) {
1133        _streamTracer->SetComputeVorticity(true);
1134        if (vtkTubeFilter::SafeDownCast(_lineFilter) == NULL) {
1135            _lineFilter = vtkSmartPointer<vtkTubeFilter>::New();
1136            _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
1137        }
1138        vtkTubeFilter *tubeFilter = vtkTubeFilter::SafeDownCast(_lineFilter);
1139        if (numSides < 3)
1140            numSides = 3;
1141        tubeFilter->SetNumberOfSides(numSides);
1142        tubeFilter->SetRadius(_dataScale * radius);
1143        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
1144        if (_faceCulling && _opacity == 1.0)
1145            setCulling(_linesActor->GetProperty(), true);
1146        _linesActor->GetProperty()->SetRepresentationToSurface();
1147        _linesActor->GetProperty()->LightingOn();
1148     }
1149}
1150
1151/**
1152 * \brief Set streamline type to 3D ribbons
1153 *
1154 * \param[in] width Minimum half-width of ribbons
1155 * \param[in] angle Default ribbon angle in degrees from normal
1156 */
1157void Streamlines::setLineTypeToRibbons(double width, double angle)
1158{
1159    _lineType = RIBBONS;
1160    if (_streamTracer != NULL) {
1161        _streamTracer->SetComputeVorticity(true);
1162        if (vtkRibbonFilter::SafeDownCast(_lineFilter) == NULL) {
1163            _lineFilter = vtkSmartPointer<vtkRibbonFilter>::New();
1164            _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
1165        }
1166        vtkRibbonFilter *ribbonFilter = vtkRibbonFilter::SafeDownCast(_lineFilter);
1167        ribbonFilter->SetWidth(_dataScale * width);
1168        ribbonFilter->SetAngle(angle);
1169        ribbonFilter->UseDefaultNormalOn();
1170        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
1171        setCulling(_linesActor->GetProperty(), false);
1172        _linesActor->GetProperty()->SetRepresentationToSurface();
1173        _linesActor->GetProperty()->LightingOn();
1174    }
1175}
1176
1177void Streamlines::updateRanges(bool useCumulative,
1178                               double scalarRange[2],
1179                               double vectorMagnitudeRange[2],
1180                               double vectorComponentRange[3][2])
1181{
1182    if (useCumulative) {
1183        _dataRange[0] = scalarRange[0];
1184        _dataRange[1] = scalarRange[1];
1185        _vectorMagnitudeRange[0] = vectorMagnitudeRange[0];
1186        _vectorMagnitudeRange[1] = vectorMagnitudeRange[1];
1187        for (int i = 0; i < 3; i++) {
1188            _vectorComponentRange[i][0] = vectorComponentRange[i][0];
1189            _vectorComponentRange[i][1] = vectorComponentRange[i][1];
1190        }
1191    } else {
1192        _dataSet->getScalarRange(_dataRange);
1193        _dataSet->getVectorRange(_vectorMagnitudeRange);
1194        for (int i = 0; i < 3; i++) {
1195            _dataSet->getVectorRange(_vectorComponentRange[i], i);
1196        }
1197    }
1198
1199    // Need to update color map ranges and/or active vector field
1200    setColorMode(_colorMode);
1201}
1202
1203void Streamlines::setColorMode(ColorMode mode)
1204{
1205    _colorMode = mode;
1206    if (_dataSet == NULL || _pdMapper == NULL)
1207        return;
1208
1209    switch (mode) {
1210    case COLOR_BY_SCALAR: {
1211        _pdMapper->ScalarVisibilityOn();
1212        _pdMapper->SetScalarModeToDefault();
1213        if (_lut != NULL) {
1214            _lut->SetRange(_dataRange);
1215        }
1216    }
1217        break;
1218    case COLOR_BY_VECTOR_MAGNITUDE: {
1219        _pdMapper->ScalarVisibilityOn();
1220        _pdMapper->SetScalarModeToUsePointFieldData();
1221        _pdMapper->SelectColorArray(_dataSet->getActiveVectorsName());
1222        if (_lut != NULL) {
1223            _lut->SetRange(_vectorMagnitudeRange);
1224            _lut->SetVectorModeToMagnitude();
1225        }
1226    }
1227        break;
1228    case COLOR_BY_VECTOR_X:
1229        _pdMapper->ScalarVisibilityOn();
1230        _pdMapper->SetScalarModeToUsePointFieldData();
1231        _pdMapper->SelectColorArray(_dataSet->getActiveVectorsName());
1232        if (_lut != NULL) {
1233            _lut->SetRange(_vectorComponentRange[0]);
1234            _lut->SetVectorModeToComponent();
1235            _lut->SetVectorComponent(0);
1236        }
1237        break;
1238    case COLOR_BY_VECTOR_Y:
1239        _pdMapper->ScalarVisibilityOn();
1240        _pdMapper->SetScalarModeToUsePointFieldData();
1241        _pdMapper->SelectColorArray(_dataSet->getActiveVectorsName());
1242        if (_lut != NULL) {
1243            _lut->SetRange(_vectorComponentRange[1]);
1244            _lut->SetVectorModeToComponent();
1245            _lut->SetVectorComponent(1);
1246        }
1247        break;
1248    case COLOR_BY_VECTOR_Z:
1249        _pdMapper->ScalarVisibilityOn();
1250        _pdMapper->SetScalarModeToUsePointFieldData();
1251        _pdMapper->SelectColorArray(_dataSet->getActiveVectorsName());
1252        if (_lut != NULL) {
1253            _lut->SetRange(_vectorComponentRange[2]);
1254            _lut->SetVectorModeToComponent();
1255            _lut->SetVectorComponent(2);
1256        }
1257        break;
1258    case COLOR_CONSTANT:
1259    default:
1260        _pdMapper->ScalarVisibilityOff();
1261        break;
1262    }
1263}
1264
1265/**
1266 * \brief Called when the color map has been edited
1267 */
1268void Streamlines::updateColorMap()
1269{
1270    setColorMap(_colorMap);
1271}
1272
1273/**
1274 * \brief Associate a colormap lookup table with the DataSet
1275 */
1276void Streamlines::setColorMap(ColorMap *cmap)
1277{
1278    if (cmap == NULL)
1279        return;
1280
1281    _colorMap = cmap;
1282 
1283    if (_lut == NULL) {
1284        _lut = vtkSmartPointer<vtkLookupTable>::New();
1285        if (_pdMapper != NULL) {
1286            _pdMapper->UseLookupTableScalarRangeOn();
1287            _pdMapper->SetLookupTable(_lut);
1288        }
1289    }
1290
1291    _lut->DeepCopy(cmap->getLookupTable());
1292
1293    switch (_colorMode) {
1294    case COLOR_CONSTANT:
1295    case COLOR_BY_SCALAR:
1296        _lut->SetRange(_dataRange);
1297        break;
1298    case COLOR_BY_VECTOR_MAGNITUDE:
1299        _lut->SetVectorModeToMagnitude();
1300        _lut->SetRange(_vectorMagnitudeRange);
1301        break;
1302    case COLOR_BY_VECTOR_X:
1303        _lut->SetVectorModeToComponent();
1304        _lut->SetVectorComponent(0);
1305        _lut->SetRange(_vectorComponentRange[0]);
1306        break;
1307    case COLOR_BY_VECTOR_Y:
1308        _lut->SetVectorModeToComponent();
1309        _lut->SetVectorComponent(1);
1310        _lut->SetRange(_vectorComponentRange[1]);
1311        break;
1312    case COLOR_BY_VECTOR_Z:
1313        _lut->SetVectorModeToComponent();
1314        _lut->SetVectorComponent(2);
1315        _lut->SetRange(_vectorComponentRange[2]);
1316        break;
1317    default:
1318         break;
1319    }
1320}
1321
1322/**
1323 * \brief Turn on/off lighting of this object
1324 */
1325void Streamlines::setLighting(bool state)
1326{
1327    _lighting = state;
1328    if (_linesActor != NULL)
1329        _linesActor->GetProperty()->SetLighting((state ? 1 : 0));
1330}
1331
1332/**
1333 * \brief Set opacity of this object
1334 */
1335void Streamlines::setOpacity(double opacity)
1336{
1337    _opacity = opacity;
1338    if (_linesActor != NULL) {
1339        _linesActor->GetProperty()->SetOpacity(_opacity);
1340        if (_opacity < 1.0)
1341            setCulling(_linesActor->GetProperty(), false);
1342        else if (_faceCulling && _lineType == TUBES)
1343            setCulling(_linesActor->GetProperty(), true);
1344    }
1345    if (_seedActor != NULL) {
1346        _seedActor->GetProperty()->SetOpacity(_opacity);
1347    }
1348}
1349
1350/**
1351 * \brief Turn on/off rendering of this Streamlines
1352 */
1353void Streamlines::setVisibility(bool state)
1354{
1355    if (_linesActor != NULL) {
1356        _linesActor->SetVisibility((state ? 1 : 0));
1357    }
1358    if (_seedActor != NULL) {
1359        if (!state ||
1360            (state && _seedVisible)) {
1361            _seedActor->SetVisibility((state ? 1 : 0));
1362        }
1363    }
1364}
1365
1366/**
1367 * \brief Turn on/off rendering of the seed geometry
1368 */
1369void Streamlines::setSeedVisibility(bool state)
1370{
1371    _seedVisible = state;
1372    if (_seedActor != NULL) {
1373        _seedActor->SetVisibility((state ? 1 : 0));
1374    }
1375}
1376
1377/**
1378 * \brief Get visibility state of the Streamlines
1379 *
1380 * \return Are the Streamlines visible?
1381 */
1382bool Streamlines::getVisibility()
1383{
1384    if (_linesActor == NULL) {
1385        return false;
1386    } else {
1387        return (_linesActor->GetVisibility() != 0);
1388    }
1389}
1390
1391/**
1392 * \brief Turn on/off rendering of edges
1393 */
1394void Streamlines::setEdgeVisibility(bool state)
1395{
1396    if (_linesActor != NULL) {
1397        _linesActor->GetProperty()->SetEdgeVisibility((state ? 1 : 0));
1398    }
1399}
1400
1401/**
1402 * \brief Set RGB color of stream lines
1403 */
1404void Streamlines::setColor(float color[3])
1405{
1406    _color[0] = color[0];
1407    _color[1] = color[1];
1408    _color[2] = color[2];
1409    if (_linesActor != NULL)
1410        _linesActor->GetProperty()->SetColor(_color[0], _color[1], _color[2]);
1411}
1412
1413/**
1414 * \brief Set RGB color of stream line edges
1415 */
1416void Streamlines::setEdgeColor(float color[3])
1417{
1418    _edgeColor[0] = color[0];
1419    _edgeColor[1] = color[1];
1420    _edgeColor[2] = color[2];
1421    if (_linesActor != NULL)
1422        _linesActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
1423}
1424
1425/**
1426 * \brief Set RGB color of seed geometry
1427 */
1428void Streamlines::setSeedColor(float color[3])
1429{
1430    _seedColor[0] = color[0];
1431    _seedColor[1] = color[1];
1432    _seedColor[2] = color[2];
1433    if (_seedActor != NULL) {
1434        _seedActor->GetProperty()->SetColor(_seedColor[0], _seedColor[1], _seedColor[2]);
1435        _seedActor->GetProperty()->SetEdgeColor(_seedColor[0], _seedColor[1], _seedColor[2]);
1436    }
1437}
1438
1439/**
1440 * \brief Set pixel width of stream lines (may be a no-op)
1441 */
1442void Streamlines::setEdgeWidth(float edgeWidth)
1443{
1444    _edgeWidth = edgeWidth;
1445    if (_linesActor != NULL)
1446        _linesActor->GetProperty()->SetLineWidth(_edgeWidth);
1447}
1448
1449/**
1450 * \brief Set a group of world coordinate planes to clip rendering
1451 *
1452 * Passing NULL for planes will remove all cliping planes
1453 */
1454void Streamlines::setClippingPlanes(vtkPlaneCollection *planes)
1455{
1456    if (_pdMapper != NULL) {
1457        _pdMapper->SetClippingPlanes(planes);
1458    }
1459    if (_seedMapper != NULL) {
1460        _seedMapper->SetClippingPlanes(planes);
1461    }
1462}
Note: See TracBrowser for help on using the repository browser.