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

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

Add support for setting glyph scaling and color by named field.

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