source: vtkvis/trunk/Streamlines.cpp @ 6226

Last change on this file since 6226 was 5835, checked in by ldelgass, 9 years ago

Require VTK >= 6.0.0

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