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

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

Apply a scaling heuristic to streamline tubes and ribbons based on cell sizes.
Protocol scale values are now relative. Also cache cell size info in
DataSet? class.

  • Property svn:eol-style set to native
File size: 33.9 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2011, Purdue Research Foundation
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cstdlib>
9#include <ctime>
10#include <cfloat>
11#include <cmath>
12
13#include <vtkMath.h>
14#include <vtkActor.h>
15#include <vtkProperty.h>
16#include <vtkPoints.h>
17#include <vtkCellArray.h>
18#include <vtkPolyLine.h>
19#include <vtkRegularPolygonSource.h>
20#include <vtkPointData.h>
21#include <vtkCellData.h>
22#include <vtkCellDataToPointData.h>
23#include <vtkPolyData.h>
24#include <vtkTubeFilter.h>
25#include <vtkRibbonFilter.h>
26#include <vtkTransform.h>
27#include <vtkTransformPolyDataFilter.h>
28
29#include "RpStreamlines.h"
30#include "Trace.h"
31
32using namespace Rappture::VtkVis;
33
34Streamlines::Streamlines() :
35    VtkGraphicsObject(),
36    _lineType(LINES),
37    _colorMode(COLOR_BY_VECTOR_MAGNITUDE),
38    _colorMap(NULL),
39    _seedVisible(true),
40    _dataScale(1)
41{
42    _faceCulling = true;
43    _color[0] = 1.0f;
44    _color[1] = 1.0f;
45    _color[2] = 1.0f;
46    _seedColor[0] = 1.0f;
47    _seedColor[1] = 1.0f;
48    _seedColor[2] = 1.0f;
49    vtkMath::RandomSeed((int)time(NULL));
50    srand((unsigned int)time(NULL));
51}
52
53Streamlines::~Streamlines()
54{
55}
56
57void Streamlines::setDataSet(DataSet *dataSet,
58                             bool useCumulative,
59                             double scalarRange[2],
60                             double vectorMagnitudeRange[2],
61                             double vectorComponentRange[3][2])
62{
63    if (_dataSet != dataSet) {
64        _dataSet = dataSet;
65
66        if (useCumulative) {
67            _dataRange[0] = scalarRange[0];
68            _dataRange[1] = scalarRange[1];
69            _vectorMagnitudeRange[0] = vectorMagnitudeRange[0];
70            _vectorMagnitudeRange[1] = vectorMagnitudeRange[1];
71            for (int i = 0; i < 3; i++) {
72                _vectorComponentRange[i][0] = vectorComponentRange[i][0];
73                _vectorComponentRange[i][1] = vectorComponentRange[i][1];
74            }
75        } else {
76            _dataSet->getScalarRange(_dataRange);
77            _dataSet->getVectorRange(_vectorMagnitudeRange);
78            for (int i = 0; i < 3; i++) {
79                _dataSet->getVectorRange(_vectorComponentRange[i], i);
80            }
81        }
82
83        update();
84    }
85}
86
87/**
88 * \brief Create and initialize a VTK Prop to render Streamlines
89 */
90void Streamlines::initProp()
91{
92    if (_linesActor == NULL) {
93        _linesActor = vtkSmartPointer<vtkActor>::New();
94        _linesActor->GetProperty()->SetColor(_color[0], _color[1], _color[2]);
95        _linesActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
96        _linesActor->GetProperty()->SetLineWidth(_edgeWidth);
97        _linesActor->GetProperty()->SetOpacity(_opacity);
98        _linesActor->GetProperty()->SetAmbient(.2);
99        if (!_lighting)
100            _linesActor->GetProperty()->LightingOff();
101        switch (_lineType) {
102        case LINES:
103            setCulling(_linesActor->GetProperty(), false);
104            _linesActor->GetProperty()->SetRepresentationToWireframe();
105            _linesActor->GetProperty()->EdgeVisibilityOff();
106            break;
107        case TUBES:
108            if (_faceCulling && _opacity == 1.0)
109                setCulling(_linesActor->GetProperty(), true);
110            _linesActor->GetProperty()->SetRepresentationToSurface();
111            _linesActor->GetProperty()->EdgeVisibilityOff();
112            break;
113        case RIBBONS:
114            setCulling(_linesActor->GetProperty(), false);
115            _linesActor->GetProperty()->SetRepresentationToSurface();
116            _linesActor->GetProperty()->EdgeVisibilityOff();
117            break;
118        default:
119            ;
120        }
121    }
122    if (_seedActor == NULL) {
123        _seedActor = vtkSmartPointer<vtkActor>::New();
124        _seedActor->GetProperty()->SetColor(_seedColor[0], _seedColor[1], _seedColor[2]);
125        _seedActor->GetProperty()->SetLineWidth(1);
126        _seedActor->GetProperty()->SetPointSize(2);
127        _seedActor->GetProperty()->SetOpacity(_opacity);
128        _seedActor->GetProperty()->SetRepresentationToWireframe();
129        _seedActor->GetProperty()->LightingOff();
130    }
131    if (_prop == NULL) {
132        _prop = vtkSmartPointer<vtkAssembly>::New();
133        getAssembly()->AddPart(_linesActor);
134        getAssembly()->AddPart(_seedActor);
135    }
136}
137
138/**
139 * \brief Get a pseudo-random number in range [min,max]
140 */
141double Streamlines::getRandomNum(double min, double max)
142{
143#if 1
144    return vtkMath::Random(min, max);
145#else
146    int r = rand();
147    return (min + ((double)r / RAND_MAX) * (max - min));
148#endif
149}
150
151/**
152 * \brief Get a random 3D point within an AABB
153 *
154 * \param[out] pt The random point
155 * \param[in] bounds The bounds of the AABB
156 */
157void Streamlines::getRandomPoint(double pt[3], const double bounds[6])
158{
159    pt[0] = getRandomNum(bounds[0], bounds[1]);
160    pt[1] = getRandomNum(bounds[2], bounds[3]);
161    pt[2] = getRandomNum(bounds[4], bounds[5]);
162}
163
164/**
165 * \brief Get a random point within a triangle (including edges)
166 *
167 * \param[out] pt The random point
168 * \param[in] v1 Triangle vertex 1
169 * \param[in] v2 Triangle vertex 2
170 * \param[in] v3 Triangle vertex 3
171 */
172void Streamlines::getRandomPointInTriangle(double pt[3],
173                                           const double v1[3],
174                                           const double v2[3],
175                                           const double v3[3])
176{
177    // Choose random barycentric coordinates
178    double bary[3];
179    bary[0] = getRandomNum(0, 1);
180    bary[1] = getRandomNum(0, 1);
181    if (bary[0] + bary[1] > 1.0) {
182        bary[0] = 1.0 - bary[0];
183        bary[1] = 1.0 - bary[1];
184    }
185    bary[2] = 1.0 - bary[0] - bary[1];
186
187    TRACE("bary %g %g %g", bary[0], bary[1], bary[2]);
188    // Convert to cartesian coords
189    for (int i = 0; i < 3; i++) {
190        pt[i] = v1[i] * bary[0] + v2[i] * bary[1] + v3[i] * bary[2];
191    }
192}
193
194/**
195 * \brief Get a random point on a line segment (including endpoints)
196 */
197void Streamlines::getRandomPointOnLineSegment(double pt[3],
198                                              const double endpt[3],
199                                              const double endpt2[3])
200{
201    double ratio = getRandomNum(0, 1);
202    pt[0] = endpt[0] + ratio * (endpt2[0] - endpt[0]);
203    pt[1] = endpt[1] + ratio * (endpt2[1] - endpt[1]);
204    pt[2] = endpt[2] + ratio * (endpt2[2] - endpt[2]);
205}
206
207/**
208 * \brief Get a random point within a vtkDataSet's mesh
209 *
210 * Note: This currently doesn't give a uniform distribution of
211 * points in space and can generate points outside the mesh
212 */
213void Streamlines::getRandomCellPt(double pt[3], vtkDataSet *ds)
214{
215    int numCells = (int)ds->GetNumberOfCells();
216    // XXX: Not uniform distribution (shouldn't use mod, and assumes
217    // all cells are equal area/volume)
218    int cell = rand() % numCells;
219    double bounds[6];
220    ds->GetCellBounds(cell, bounds);
221    // Note: point is inside AABB of cell, but may be outside the cell
222    getRandomPoint(pt, bounds);
223}
224
225/**
226 * \brief Internal method to set up pipeline after a state change
227 */
228void Streamlines::update()
229{
230    if (_dataSet == NULL) {
231        return;
232    }
233
234    vtkDataSet *ds = _dataSet->getVtkDataSet();
235
236    double bounds[6];
237    _dataSet->getBounds(bounds);
238    double maxBound = 0.0;
239    if (bounds[1] - bounds[0] > maxBound) {
240        maxBound = bounds[1] - bounds[0];
241    }
242    if (bounds[3] - bounds[2] > maxBound) {
243        maxBound = bounds[3] - bounds[2];
244    }
245    if (bounds[5] - bounds[4] > maxBound) {
246        maxBound = bounds[5] - bounds[4];
247    }
248
249    double cellSizeRange[2];
250    double avgSize;
251    _dataSet->getCellSizeRange(cellSizeRange, &avgSize);
252    _dataScale = avgSize / 8.;
253
254    vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
255
256    if (ds->GetPointData() == NULL ||
257        ds->GetPointData()->GetVectors() == NULL) {
258        TRACE("No vector point data found in DataSet %s", _dataSet->getName().c_str());
259        if (ds->GetCellData() == NULL ||
260            ds->GetCellData()->GetVectors() == NULL) {
261            ERROR("No vectors found in DataSet %s", _dataSet->getName().c_str());
262        } else {
263            cellToPtData =
264                vtkSmartPointer<vtkCellDataToPointData>::New();
265            cellToPtData->SetInput(ds);
266            //cellToPtData->PassCellDataOn();
267            cellToPtData->Update();
268            ds = cellToPtData->GetOutput();
269        }
270    }
271
272    if (_streamTracer == NULL) {
273        _streamTracer = vtkSmartPointer<vtkStreamTracer>::New();
274    }
275
276    _streamTracer->SetInput(ds);
277    _streamTracer->SetMaximumPropagation(maxBound);
278
279    if (_pdMapper == NULL) {
280        _pdMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
281        _pdMapper->SetResolveCoincidentTopologyToPolygonOffset();
282        _pdMapper->ScalarVisibilityOn();
283    }
284    if (_seedMapper == NULL) {
285        _seedMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
286        _seedMapper->SetResolveCoincidentTopologyToPolygonOffset();
287    }
288
289    // Set up seed source object
290    setSeedToRandomPoints(200);
291
292    switch (_lineType) {
293    case LINES: {
294        _streamTracer->SetComputeVorticity(false);
295        _pdMapper->SetInputConnection(_streamTracer->GetOutputPort());
296    }
297        break;
298    case TUBES: {
299        _streamTracer->SetComputeVorticity(true);
300        _lineFilter = vtkSmartPointer<vtkTubeFilter>::New();
301        vtkTubeFilter *tubeFilter = vtkTubeFilter::SafeDownCast(_lineFilter);
302        tubeFilter->SetNumberOfSides(5);
303        _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
304        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
305    }
306        break;
307    case RIBBONS: {
308        _streamTracer->SetComputeVorticity(true);
309        _lineFilter = vtkSmartPointer<vtkRibbonFilter>::New();
310        _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
311        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
312    }
313        break;
314    default:
315        ERROR("Unknown LineType: %d", _lineType);
316    }
317
318#if defined(DEBUG) && defined(WANT_TRACE)
319    _streamTracer->Update();
320    vtkPolyData *pd = _streamTracer->GetOutput();
321    TRACE("Verts: %d Lines: %d Polys: %d Strips: %d",
322                  pd->GetNumberOfVerts(),
323                  pd->GetNumberOfLines(),
324                  pd->GetNumberOfPolys(),
325                  pd->GetNumberOfStrips());
326#endif
327
328    initProp();
329
330    _seedActor->SetMapper(_seedMapper);
331
332    if (_lut == NULL) {
333        setColorMap(ColorMap::getDefault());
334    }
335
336    setColorMode(_colorMode);
337
338    _linesActor->SetMapper(_pdMapper);
339    _pdMapper->Update();
340    _seedMapper->Update();
341}
342
343/**
344 * \brief Use randomly distributed seed points
345 *
346 * Note: The current implementation doesn't give a uniform
347 * distribution of points, and points outside the mesh bounds
348 * may be generated
349 *
350 * \param[in] numPoints Number of random seed points to generate
351 */
352void Streamlines::setSeedToRandomPoints(int numPoints)
353{
354    if (_streamTracer != NULL) {
355        // Set up seed source object
356        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
357        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
358        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
359
360        for (int i = 0; i < numPoints; i++) {
361            double pt[3];
362            getRandomCellPt(pt, _dataSet->getVtkDataSet());
363            TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
364            pts->InsertNextPoint(pt);
365            cells->InsertNextCell(1);
366            cells->InsertCellPoint(i);
367        }
368
369        seed->SetPoints(pts);
370        seed->SetVerts(cells);
371
372        TRACE("Seed points: %d", seed->GetNumberOfPoints());
373        vtkSmartPointer<vtkDataSet> oldSeed;
374        if (_streamTracer->GetSource() != NULL) {
375            oldSeed = _streamTracer->GetSource();
376        }
377
378        _streamTracer->SetSource(seed);
379        if (oldSeed != NULL) {
380            oldSeed->SetPipelineInformation(NULL);
381        }
382
383        _seedMapper->SetInput(seed);
384    }
385}
386
387/**
388 * \brief Use seed points along a line
389 *
390 * \param[in] start Starting point of rake line
391 * \param[in] end End point of rake line
392 * \param[in] numPoints Number of points along line to generate
393 */
394void Streamlines::setSeedToRake(double start[3], double end[3], int numPoints)
395{
396    if (numPoints < 2)
397        return;
398    if (_streamTracer != NULL) {
399        // Set up seed source object
400        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
401        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
402        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
403        vtkSmartPointer<vtkPolyLine> polyline = vtkSmartPointer<vtkPolyLine>::New();
404
405        double dir[3];
406        for (int i = 0; i < 3; i++) {
407            dir[i] = end[i] - start[i];
408        }
409
410        polyline->GetPointIds()->SetNumberOfIds(numPoints);
411        for (int i = 0; i < numPoints; i++) {
412            double pt[3];
413            for (int ii = 0; ii < 3; ii++) {
414                pt[ii] = start[ii] + dir[ii] * ((double)i / (numPoints-1));
415            }
416            TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
417            pts->InsertNextPoint(pt);
418            polyline->GetPointIds()->SetId(i, i);
419        }
420
421        cells->InsertNextCell(polyline);
422        seed->SetPoints(pts);
423        seed->SetLines(cells);
424
425        TRACE("Seed points: %d", seed->GetNumberOfPoints());
426        vtkSmartPointer<vtkDataSet> oldSeed;
427        if (_streamTracer->GetSource() != NULL) {
428            oldSeed = _streamTracer->GetSource();
429        }
430
431        _streamTracer->SetSource(seed);
432        if (oldSeed != NULL) {
433            oldSeed->SetPipelineInformation(NULL);
434        }
435
436        _seedMapper->SetInput(seed);
437    }
438}
439
440/**
441 * \brief Create seed points inside a disk with an optional hole
442 *
443 * \param[in] center Center point of disk
444 * \param[in] normal Normal vector to orient disk
445 * \param[in] radius Radius of disk
446 * \param[in] innerRadius Radius of hole at center of disk
447 * \param[in] numPoints Number of random points to generate
448 */
449void Streamlines::setSeedToDisk(double center[3],
450                                double normal[3],
451                                double radius,
452                                double innerRadius,
453                                int numPoints)
454{
455    if (_streamTracer != NULL) {
456        // Set up seed source object
457        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
458        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
459        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
460
461        // The following code is based on vtkRegularPolygonSource::RequestData
462
463        double px[3];
464        double py[3];
465        double axis[3] = {1., 0., 0.};
466
467        if (vtkMath::Normalize(normal) == 0.0) {
468            normal[0] = 0.0;
469            normal[1] = 0.0;
470            normal[2] = 1.0;
471        }
472
473        // Find axis in plane (orthogonal to normal)
474        bool done = false;
475        vtkMath::Cross(normal, axis, px);
476        if (vtkMath::Normalize(px) > 1.0e-3) {
477            done = true;
478        }
479        if (!done) {
480            axis[0] = 0.0;
481            axis[1] = 1.0;
482            axis[2] = 0.0;
483            vtkMath::Cross(normal, axis, px);
484            if (vtkMath::Normalize(px) > 1.0e-3) {
485                done = true;
486            }
487        }
488        if (!done) {
489            axis[0] = 0.0;
490            axis[1] = 0.0;
491            axis[2] = 1.0;
492            vtkMath::Cross(normal, axis, px);
493            vtkMath::Normalize(px);
494        }
495        // Create third orthogonal basis vector
496        vtkMath::Cross(px, normal, py);
497
498        double minSquared = (innerRadius*innerRadius)/(radius*radius);
499        for (int j = 0; j < numPoints; j++) {
500            // Get random sweep angle and radius
501            double angle = getRandomNum(0, 2.0 * vtkMath::DoublePi());
502            // Need sqrt to get uniform distribution
503            double r = sqrt(getRandomNum(minSquared, 1)) * radius;
504            double pt[3];
505            for (int i = 0; i < 3; i++) {
506                pt[i] = center[i] + r * (px[i] * cos(angle) + py[i] * sin(angle));
507            }
508            TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
509            pts->InsertNextPoint(pt);
510            cells->InsertNextCell(1);
511            cells->InsertCellPoint(j);
512        }
513
514        seed->SetPoints(pts);
515        seed->SetVerts(cells);
516
517        TRACE("Seed points: %d", seed->GetNumberOfPoints());
518        vtkSmartPointer<vtkDataSet> oldSeed;
519        if (_streamTracer->GetSource() != NULL) {
520            oldSeed = _streamTracer->GetSource();
521        }
522
523        _streamTracer->SetSource(seed);
524        if (oldSeed != NULL) {
525            oldSeed->SetPipelineInformation(NULL);
526        }
527
528        _seedMapper->SetInput(seed);
529    }
530}
531
532/**
533 * \brief Use seed points from an n-sided polygon
534 *
535 * \param[in] center Center point of polygon
536 * \param[in] normal Normal vector to orient polygon
537 * \param[in] angle Angle in degrees to rotate about normal
538 * \param[in] radius Radius of circumscribing circle
539 * \param[in] numSides Number of polygon sides (and points) to generate
540 */
541void Streamlines::setSeedToPolygon(double center[3],
542                                   double normal[3],
543                                   double angle,
544                                   double radius,
545                                   int numSides)
546{
547    if (_streamTracer != NULL) {
548        // Set up seed source object
549        vtkSmartPointer<vtkRegularPolygonSource> seed = vtkSmartPointer<vtkRegularPolygonSource>::New();
550
551        seed->SetCenter(center);
552        seed->SetNormal(normal);
553        seed->SetRadius(radius);
554        seed->SetNumberOfSides(numSides);
555        seed->GeneratePolygonOn();
556
557        if (angle != 0.0) {
558            vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
559            trans->RotateWXYZ(angle, normal);
560            vtkSmartPointer<vtkTransformPolyDataFilter> transFilt =
561                vtkSmartPointer<vtkTransformPolyDataFilter>::New();
562            transFilt->SetInputConnection(seed->GetOutputPort());
563            transFilt->SetTransform(trans);
564        }
565
566        TRACE("Seed points: %d", numSides);
567        vtkSmartPointer<vtkDataSet> oldSeed;
568        if (_streamTracer->GetSource() != NULL) {
569            oldSeed = _streamTracer->GetSource();
570        }
571
572        if (angle != 0.0) {
573            vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
574            trans->Translate(+center[0], +center[1], +center[2]);
575            trans->RotateWXYZ(angle, normal);
576            trans->Translate(-center[0], -center[1], -center[2]);
577            vtkSmartPointer<vtkTransformPolyDataFilter> transFilt =
578                vtkSmartPointer<vtkTransformPolyDataFilter>::New();
579            transFilt->SetInputConnection(seed->GetOutputPort());
580            transFilt->SetTransform(trans);
581            _streamTracer->SetSourceConnection(transFilt->GetOutputPort());
582            _seedMapper->SetInputConnection(transFilt->GetOutputPort());
583        } else {
584            _streamTracer->SetSourceConnection(seed->GetOutputPort());
585            _seedMapper->SetInputConnection(seed->GetOutputPort());
586        }
587
588        if (oldSeed != NULL) {
589            oldSeed->SetPipelineInformation(NULL);
590        }
591    }
592}
593
594/**
595 * \brief Use seed points from an n-sided polygon
596 *
597 * \param[in] center Center point of polygon
598 * \param[in] normal Normal vector to orient polygon
599 * \param[in] angle Angle in degrees to rotate about normal
600 * \param[in] radius Radius of circumscribing circle
601 * \param[in] numSides Number of polygon sides (and points) to generate
602 * \param[in] numPoints Number of random points to generate
603 */
604void Streamlines::setSeedToFilledPolygon(double center[3],
605                                         double normal[3],
606                                         double angle,
607                                         double radius,
608                                         int numSides,
609                                         int numPoints)
610{
611    if (_streamTracer != NULL) {
612         // Set up seed source object
613        vtkSmartPointer<vtkPolyData> seed = vtkSmartPointer<vtkPolyData>::New();
614        vtkSmartPointer<vtkPoints> pts = vtkSmartPointer<vtkPoints>::New();
615        vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
616
617        // The following code is based on vtkRegularPolygonSource::RequestData
618
619        double px[3];
620        double py[3];
621        double axis[3] = {1., 0., 0.};
622
623        if (vtkMath::Normalize(normal) == 0.0) {
624            normal[0] = 0.0;
625            normal[1] = 0.0;
626            normal[2] = 1.0;
627        }
628
629        // Find axis in plane (orthogonal to normal)
630        bool done = false;
631        vtkMath::Cross(normal, axis, px);
632        if (vtkMath::Normalize(px) > 1.0e-3) {
633            done = true;
634        }
635        if (!done) {
636            axis[0] = 0.0;
637            axis[1] = 1.0;
638            axis[2] = 0.0;
639            vtkMath::Cross(normal, axis, px);
640            if (vtkMath::Normalize(px) > 1.0e-3) {
641                done = true;
642            }
643        }
644        if (!done) {
645            axis[0] = 0.0;
646            axis[1] = 0.0;
647            axis[2] = 1.0;
648            vtkMath::Cross(normal, axis, px);
649            vtkMath::Normalize(px);
650        }
651        // Create third orthogonal basis vector
652        vtkMath::Cross(px, normal, py);
653
654        double verts[numSides][3];
655        double sliceTheta = 2.0 * vtkMath::DoublePi() / (double)numSides;
656        angle = vtkMath::RadiansFromDegrees(angle);
657        for (int j = 0; j < numSides; j++) {
658            for (int i = 0; i < 3; i++) {
659                double theta = sliceTheta * (double)j - angle;
660                verts[j][i] = center[i] + radius * (px[i] * cos(theta) +
661                                                    py[i] * sin(theta));
662            }
663            TRACE("Vert %d: %g %g %g", j, verts[j][0], verts[j][1], verts[j][2]);
664        }
665
666        // Note: this gives a uniform distribution because the polygon is regular and
667        // the triangular sections have equal area
668        if (numSides == 3) {
669            for (int j = 0; j < numPoints; j++) {
670                double pt[3];
671                getRandomPointInTriangle(pt, verts[0], verts[1], verts[2]);
672                TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
673                pts->InsertNextPoint(pt);
674                cells->InsertNextCell(1);
675                cells->InsertCellPoint(j);
676            }
677        } else {
678            for (int j = 0; j < numPoints; j++) {
679                // Get random triangle section
680                int tri = rand() % numSides;
681                double pt[3];
682                getRandomPointInTriangle(pt, center, verts[tri], verts[(tri+1) % numSides]);
683                TRACE("Seed pt: %g %g %g", pt[0], pt[1], pt[2]);
684                pts->InsertNextPoint(pt);
685                cells->InsertNextCell(1);
686                cells->InsertCellPoint(j);
687            }
688        }
689
690        seed->SetPoints(pts);
691        seed->SetVerts(cells);
692
693        TRACE("Seed points: %d", seed->GetNumberOfPoints());
694        vtkSmartPointer<vtkDataSet> oldSeed;
695        if (_streamTracer->GetSource() != NULL) {
696            oldSeed = _streamTracer->GetSource();
697        }
698
699        _streamTracer->SetSource(seed);
700        if (oldSeed != NULL) {
701            oldSeed->SetPipelineInformation(NULL);
702        }
703
704        _seedMapper->SetInput(seed);
705    }
706}
707
708/**
709 * \brief Set maximum length of stream lines in world coordinates
710 */
711void Streamlines::setMaxPropagation(double length)
712{
713    if (_streamTracer != NULL) {
714        _streamTracer->SetMaximumPropagation(length);
715    }
716}
717
718/**
719 * \brief Set streamline type to polylines
720 */
721void Streamlines::setLineTypeToLines()
722{
723    _lineType = LINES;
724    if (_streamTracer != NULL &&
725        _pdMapper != NULL) {
726        _streamTracer->SetComputeVorticity(false);
727        _pdMapper->SetInputConnection(_streamTracer->GetOutputPort());
728        _lineFilter = NULL;
729        setCulling(_linesActor->GetProperty(), false);
730        _linesActor->GetProperty()->SetRepresentationToWireframe();
731        _linesActor->GetProperty()->LightingOff();
732    }
733}
734
735/**
736 * \brief Set streamline type to 3D tubes
737 *
738 * \param[in] numSides Number of sides (>=3) for tubes
739 * \param[in] radius World coordinate minimum tube radius
740 */
741void Streamlines::setLineTypeToTubes(int numSides, double radius)
742{
743    _lineType = TUBES;
744    if (_streamTracer != NULL) {
745        _streamTracer->SetComputeVorticity(true);
746        if (vtkTubeFilter::SafeDownCast(_lineFilter) == NULL) {
747            _lineFilter = vtkSmartPointer<vtkTubeFilter>::New();
748            _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
749        }
750        vtkTubeFilter *tubeFilter = vtkTubeFilter::SafeDownCast(_lineFilter);
751        if (numSides < 3)
752            numSides = 3;
753        tubeFilter->SetNumberOfSides(numSides);
754        tubeFilter->SetRadius(_dataScale * radius);
755        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
756        if (_faceCulling && _opacity == 1.0)
757            setCulling(_linesActor->GetProperty(), true);
758        _linesActor->GetProperty()->SetRepresentationToSurface();
759        _linesActor->GetProperty()->LightingOn();
760     }
761}
762
763/**
764 * \brief Set streamline type to 3D ribbons
765 *
766 * \param[in] width Minimum half-width of ribbons
767 * \param[in] angle Default ribbon angle in degrees from normal
768 */
769void Streamlines::setLineTypeToRibbons(double width, double angle)
770{
771    _lineType = RIBBONS;
772    if (_streamTracer != NULL) {
773        _streamTracer->SetComputeVorticity(true);
774        if (vtkRibbonFilter::SafeDownCast(_lineFilter) == NULL) {
775            _lineFilter = vtkSmartPointer<vtkRibbonFilter>::New();
776            _lineFilter->SetInputConnection(_streamTracer->GetOutputPort());
777        }
778        vtkRibbonFilter *ribbonFilter = vtkRibbonFilter::SafeDownCast(_lineFilter);
779        ribbonFilter->SetWidth(_dataScale * width);
780        ribbonFilter->SetAngle(angle);
781        ribbonFilter->UseDefaultNormalOn();
782        _pdMapper->SetInputConnection(_lineFilter->GetOutputPort());
783        setCulling(_linesActor->GetProperty(), false);
784        _linesActor->GetProperty()->SetRepresentationToSurface();
785        _linesActor->GetProperty()->LightingOn();
786    }
787}
788
789void Streamlines::updateRanges(bool useCumulative,
790                               double scalarRange[2],
791                               double vectorMagnitudeRange[2],
792                               double vectorComponentRange[3][2])
793{
794    if (useCumulative) {
795        _dataRange[0] = scalarRange[0];
796        _dataRange[1] = scalarRange[1];
797        _vectorMagnitudeRange[0] = vectorMagnitudeRange[0];
798        _vectorMagnitudeRange[1] = vectorMagnitudeRange[1];
799        for (int i = 0; i < 3; i++) {
800            _vectorComponentRange[i][0] = vectorComponentRange[i][0];
801            _vectorComponentRange[i][1] = vectorComponentRange[i][1];
802        }
803    } else {
804        _dataSet->getScalarRange(_dataRange);
805        _dataSet->getVectorRange(_vectorMagnitudeRange);
806        for (int i = 0; i < 3; i++) {
807            _dataSet->getVectorRange(_vectorComponentRange[i], i);
808        }
809    }
810
811    // Need to update color map ranges and/or active vector field
812    setColorMode(_colorMode);
813}
814
815void Streamlines::setColorMode(ColorMode mode)
816{
817    _colorMode = mode;
818    if (_dataSet == NULL || _pdMapper == NULL)
819        return;
820
821    vtkDataSet *ds = _dataSet->getVtkDataSet();
822
823    switch (mode) {
824    case COLOR_BY_SCALAR: {
825        _pdMapper->ScalarVisibilityOn();
826        _pdMapper->SetScalarModeToDefault();
827        if (_lut != NULL) {
828            _lut->SetRange(_dataRange);
829        }
830    }
831        break;
832    case COLOR_BY_VECTOR_MAGNITUDE: {
833        _pdMapper->ScalarVisibilityOn();
834        _pdMapper->SetScalarModeToUsePointFieldData();
835        if (ds->GetPointData() != NULL &&
836            ds->GetPointData()->GetVectors() != NULL) {
837            _pdMapper->SelectColorArray(ds->GetPointData()->GetVectors()->GetName());
838        }
839        if (_lut != NULL) {
840            _lut->SetRange(_vectorMagnitudeRange);
841            _lut->SetVectorModeToMagnitude();
842        }
843    }
844        break;
845    case COLOR_BY_VECTOR_X:
846        _pdMapper->ScalarVisibilityOn();
847        _pdMapper->SetScalarModeToUsePointFieldData();
848        if (ds->GetPointData() != NULL &&
849            ds->GetPointData()->GetVectors() != NULL) {
850            _pdMapper->SelectColorArray(ds->GetPointData()->GetVectors()->GetName());
851        }
852        if (_lut != NULL) {
853            _lut->SetRange(_vectorComponentRange[0]);
854            _lut->SetVectorModeToComponent();
855            _lut->SetVectorComponent(0);
856        }
857        break;
858    case COLOR_BY_VECTOR_Y:
859        _pdMapper->ScalarVisibilityOn();
860        _pdMapper->SetScalarModeToUsePointFieldData();
861        if (ds->GetPointData() != NULL &&
862            ds->GetPointData()->GetVectors() != NULL) {
863            _pdMapper->SelectColorArray(ds->GetPointData()->GetVectors()->GetName());
864        }
865        if (_lut != NULL) {
866            _lut->SetRange(_vectorComponentRange[1]);
867            _lut->SetVectorModeToComponent();
868            _lut->SetVectorComponent(1);
869        }
870        break;
871    case COLOR_BY_VECTOR_Z:
872        _pdMapper->ScalarVisibilityOn();
873        _pdMapper->SetScalarModeToUsePointFieldData();
874        if (ds->GetPointData() != NULL &&
875            ds->GetPointData()->GetVectors() != NULL) {
876            _pdMapper->SelectColorArray(ds->GetPointData()->GetVectors()->GetName());
877        }
878        if (_lut != NULL) {
879            _lut->SetRange(_vectorComponentRange[2]);
880            _lut->SetVectorModeToComponent();
881            _lut->SetVectorComponent(2);
882        }
883        break;
884    case COLOR_CONSTANT:
885    default:
886        _pdMapper->ScalarVisibilityOff();
887        break;
888    }
889}
890
891/**
892 * \brief Called when the color map has been edited
893 */
894void Streamlines::updateColorMap()
895{
896    setColorMap(_colorMap);
897}
898
899/**
900 * \brief Associate a colormap lookup table with the DataSet
901 */
902void Streamlines::setColorMap(ColorMap *cmap)
903{
904    if (cmap == NULL)
905        return;
906
907    _colorMap = cmap;
908 
909    if (_lut == NULL) {
910        _lut = vtkSmartPointer<vtkLookupTable>::New();
911        if (_pdMapper != NULL) {
912            _pdMapper->UseLookupTableScalarRangeOn();
913            _pdMapper->SetLookupTable(_lut);
914        }
915    }
916
917    _lut->DeepCopy(cmap->getLookupTable());
918
919    switch (_colorMode) {
920    case COLOR_CONSTANT:
921    case COLOR_BY_SCALAR:
922        _lut->SetRange(_dataRange);
923        break;
924    case COLOR_BY_VECTOR_MAGNITUDE:
925        _lut->SetVectorModeToMagnitude();
926        _lut->SetRange(_vectorMagnitudeRange);
927        break;
928    case COLOR_BY_VECTOR_X:
929        _lut->SetVectorModeToComponent();
930        _lut->SetVectorComponent(0);
931        _lut->SetRange(_vectorComponentRange[0]);
932        break;
933    case COLOR_BY_VECTOR_Y:
934        _lut->SetVectorModeToComponent();
935        _lut->SetVectorComponent(1);
936        _lut->SetRange(_vectorComponentRange[1]);
937        break;
938    case COLOR_BY_VECTOR_Z:
939        _lut->SetVectorModeToComponent();
940        _lut->SetVectorComponent(2);
941        _lut->SetRange(_vectorComponentRange[2]);
942        break;
943    default:
944         break;
945    }
946}
947
948/**
949 * \brief Turn on/off lighting of this object
950 */
951void Streamlines::setLighting(bool state)
952{
953    _lighting = state;
954    if (_linesActor != NULL)
955        _linesActor->GetProperty()->SetLighting((state ? 1 : 0));
956}
957
958/**
959 * \brief Set opacity of this object
960 */
961void Streamlines::setOpacity(double opacity)
962{
963    _opacity = opacity;
964    if (_linesActor != NULL) {
965        _linesActor->GetProperty()->SetOpacity(_opacity);
966        if (_opacity < 1.0)
967            setCulling(_linesActor->GetProperty(), false);
968        else if (_faceCulling && _lineType == TUBES)
969            setCulling(_linesActor->GetProperty(), true);
970    }
971    if (_seedActor != NULL) {
972        _seedActor->GetProperty()->SetOpacity(_opacity);
973    }
974}
975
976/**
977 * \brief Turn on/off rendering of this Streamlines
978 */
979void Streamlines::setVisibility(bool state)
980{
981    if (_linesActor != NULL) {
982        _linesActor->SetVisibility((state ? 1 : 0));
983    }
984    if (_seedActor != NULL) {
985        if (!state ||
986            (state && _seedVisible)) {
987            _seedActor->SetVisibility((state ? 1 : 0));
988        }
989    }
990}
991
992/**
993 * \brief Turn on/off rendering of the seed geometry
994 */
995void Streamlines::setSeedVisibility(bool state)
996{
997    _seedVisible = state;
998    if (_seedActor != NULL) {
999        _seedActor->SetVisibility((state ? 1 : 0));
1000    }
1001}
1002
1003/**
1004 * \brief Get visibility state of the Streamlines
1005 *
1006 * \return Are the Streamlines visible?
1007 */
1008bool Streamlines::getVisibility()
1009{
1010    if (_linesActor == NULL) {
1011        return false;
1012    } else {
1013        return (_linesActor->GetVisibility() != 0);
1014    }
1015}
1016
1017/**
1018 * \brief Turn on/off rendering of edges
1019 */
1020void Streamlines::setEdgeVisibility(bool state)
1021{
1022    if (_linesActor != NULL) {
1023        _linesActor->GetProperty()->SetEdgeVisibility((state ? 1 : 0));
1024    }
1025}
1026
1027/**
1028 * \brief Set RGB color of stream lines
1029 */
1030void Streamlines::setColor(float color[3])
1031{
1032    _color[0] = color[0];
1033    _color[1] = color[1];
1034    _color[2] = color[2];
1035    if (_linesActor != NULL)
1036        _linesActor->GetProperty()->SetColor(_color[0], _color[1], _color[2]);
1037}
1038
1039/**
1040 * \brief Set RGB color of stream line edges
1041 */
1042void Streamlines::setEdgeColor(float color[3])
1043{
1044    _edgeColor[0] = color[0];
1045    _edgeColor[1] = color[1];
1046    _edgeColor[2] = color[2];
1047    if (_linesActor != NULL)
1048        _linesActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
1049}
1050
1051/**
1052 * \brief Set RGB color of seed geometry
1053 */
1054void Streamlines::setSeedColor(float color[3])
1055{
1056    _seedColor[0] = color[0];
1057    _seedColor[1] = color[1];
1058    _seedColor[2] = color[2];
1059    if (_seedActor != NULL)
1060        _seedActor->GetProperty()->SetColor(_seedColor[0], _seedColor[1], _seedColor[2]);
1061}
1062
1063/**
1064 * \brief Set pixel width of stream lines (may be a no-op)
1065 */
1066void Streamlines::setEdgeWidth(float edgeWidth)
1067{
1068    _edgeWidth = edgeWidth;
1069    if (_linesActor != NULL)
1070        _linesActor->GetProperty()->SetLineWidth(_edgeWidth);
1071}
1072
1073/**
1074 * \brief Set a group of world coordinate planes to clip rendering
1075 *
1076 * Passing NULL for planes will remove all cliping planes
1077 */
1078void Streamlines::setClippingPlanes(vtkPlaneCollection *planes)
1079{
1080    if (_pdMapper != NULL) {
1081        _pdMapper->SetClippingPlanes(planes);
1082    }
1083    if (_seedMapper != NULL) {
1084        _seedMapper->SetClippingPlanes(planes);
1085    }
1086}
Note: See TracBrowser for help on using the repository browser.