source: trunk/packages/vizservers/vtkvis/Streamlines.h @ 4609

Last change on this file since 4609 was 3682, checked in by ldelgass, 11 years ago

Streamlines:

  • Set a default max. of 500 seeds for points of a mesh and make 'streamlines seed numpts' set the max when using this seed type (setting a numpts < 0 will give the previous behavior of using all points of the mesh).
  • Add point size option for seeds.
  • Add protocol for setting streamlines terminal speed.
  • NOTE: Streamlines on 2D meshes seem to only work if there is no plane offset

LIC:

  • Add support for clouds by meshing before resampling to an image.
  • Property svn:eol-style set to native
File size: 5.9 KB
RevLine 
[2317]1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
[3177]3 * Copyright (C) 2004-2012  HUBzero Foundation, LLC
[2317]4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
[3615]8#ifndef VTKVIS_STREAMLINES_H
9#define VTKVIS_STREAMLINES_H
[2317]10
11#include <vtkSmartPointer.h>
12#include <vtkProp.h>
13#include <vtkPlaneCollection.h>
14#include <vtkStreamTracer.h>
15#include <vtkPolyDataAlgorithm.h>
16#include <vtkPolyDataMapper.h>
17#include <vtkLookupTable.h>
[2328]18#include <vtkAssembly.h>
[2317]19
[2402]20#include "ColorMap.h"
[3621]21#include "GraphicsObject.h"
22#include "DataSet.h"
[2317]23
24namespace VtkVis {
25
26/**
27 * \brief Streamline visualization of vector fields
[2328]28 *
29 * The DataSet must contain vectors
[2317]30 */
[3621]31class Streamlines : public GraphicsObject {
[2317]32public:
33    enum LineType {
34        LINES,
35        TUBES,
36        RIBBONS
37    };
[2393]38    enum ColorMode {
39        COLOR_BY_SCALAR,
40        COLOR_BY_VECTOR_MAGNITUDE,
41        COLOR_BY_VECTOR_X,
42        COLOR_BY_VECTOR_Y,
43        COLOR_BY_VECTOR_Z,
44        COLOR_CONSTANT
45    };
[2510]46    enum StepUnit {
47        LENGTH_UNIT,
48        CELL_LENGTH_UNIT
49    };
50    enum IntegratorType {
51        RUNGE_KUTTA2,
52        RUNGE_KUTTA4,
53        RUNGE_KUTTA45
54    };
55    enum IntegrationDirection {
56        FORWARD,
57        BACKWARD,
58        BOTH
59    };
[2581]60    enum SeedType {
61        DATASET_MESH_POINTS,
62        DATASET_FILLED_MESH,
63        MESH_POINTS,
64        FILLED_MESH,
65        RAKE,
66        DISK,
67        POLYGON,
68        FILLED_POLYGON
69    };
[2317]70
71    Streamlines();
72    virtual ~Streamlines();
73
[2328]74    virtual const char *getClassName() const
75    {
76        return "Streamlines";
77    }
[2317]78
[2402]79    virtual void setDataSet(DataSet *dataSet,
[2612]80                            Renderer *renderer);
[2402]81
[2328]82    virtual void setLighting(bool state);
[2332]83
84    virtual void setOpacity(double opacity);
85
[2328]86    virtual void setVisibility(bool state);
[2317]87
[2328]88    virtual bool getVisibility();
[2317]89
[2393]90    virtual void setColor(float color[3]);
91
[2328]92    virtual void setEdgeVisibility(bool state);
93
94    virtual void setEdgeColor(float color[3]);
95
96    virtual void setEdgeWidth(float edgeWidth);
97
98    virtual void setClippingPlanes(vtkPlaneCollection *planes);
99
[2581]100    void setNumberOfSeedPoints(int numPoints);
101
[3682]102    void setSeedToMeshPoints(int maxPoints = 500);
[2317]103
[2506]104    void setSeedToFilledMesh(int numPoints);
105
[3682]106    void setSeedToMeshPoints(vtkDataSet *ds, int maxPoints = 500);
[2506]107
108    void setSeedToFilledMesh(vtkDataSet *ds, int numPoints);
109
[2317]110    void setSeedToRake(double start[3], double end[3], int numPoints);
111
[2349]112    void setSeedToDisk(double center[3], double normal[3],
113                       double radius, double innerRadius, int numPoints);
114
[2318]115    void setSeedToPolygon(double center[3], double normal[3],
[2349]116                          double angle, double radius,
117                          int numSides);
[2318]118
[2349]119    void setSeedToFilledPolygon(double center[3], double normal[3],
120                                double angle, double radius,
121                                int numSides, int numPoints);
122
[2510]123    void setIntegrator(IntegratorType integrator);
124
125    void setIntegrationDirection(IntegrationDirection dir);
126
127    void setIntegrationStepUnit(StepUnit unit);
128
129    void setInitialIntegrationStep(double step);
130
131    void setMinimumIntegrationStep(double step);
132
133    void setMaximumIntegrationStep(double step);
134
135    void setMaximumError(double error);
136
[2317]137    void setMaxPropagation(double length);
138
[2510]139    void setMaxNumberOfSteps(int steps);
140
141    void setTerminalSpeed(double speed);
142
[2317]143    void setLineTypeToLines();
144
145    void setLineTypeToTubes(int numSides, double radius);
146
147    void setLineTypeToRibbons(double width, double angle);
148
[2612]149    void setColorMode(ColorMode mode, DataSet::DataAttributeType type,
150                      const char *name, double range[2] = NULL);
151
152    void setColorMode(ColorMode mode,
153                      const char *name, double range[2] = NULL);
154
[2393]155    void setColorMode(ColorMode mode);
156
[2402]157    void setColorMap(ColorMap *colorMap);
[2317]158
[2402]159    /**
160     * \brief Return the ColorMap in use
161     */
162    ColorMap *getColorMap()
163    {
164        return _colorMap;
165    }
[2317]166
[2402]167    void updateColorMap();
168
[2612]169    virtual void updateRanges(Renderer *renderer);
[2402]170
[2317]171    void setSeedVisibility(bool state);
172
173    void setSeedColor(float color[3]);
174
[3682]175    void setSeedPointSize(float size);
176
[2317]177private:
[2328]178    virtual void initProp();
179    virtual void update();
[2317]180
[2349]181    static double getRandomNum(double min, double max);
[2317]182    static void getRandomPoint(double pt[3], const double bounds[6]);
[2349]183    static void getRandomPointInTriangle(double pt[3],
[2507]184                                         const double v0[3],
[2349]185                                         const double v1[3],
[2507]186                                         const double v2[3]);
[2349]187    static void getRandomPointOnLineSegment(double pt[3],
188                                            const double endpt[3],
189                                            const double endpt2[3]);
[2507]190    static void getRandomPointInTetrahedron(double pt[3],
191                                            const double v0[3],
192                                            const double v1[3],
193                                            const double v2[3],
194                                            const double v3[3]);
[2349]195    static void getRandomCellPt(double pt[3], vtkDataSet *ds);
[2317]196
197    LineType _lineType;
[2612]198    ColorMap *_colorMap;
[2393]199    ColorMode _colorMode;
[2612]200    std::string _colorFieldName;
201    DataSet::DataAttributeType _colorFieldType;
202    double _colorFieldRange[2];
[2581]203    SeedType _seedType;
[2317]204    float _seedColor[3];
205    bool _seedVisible;
[2402]206    double _vectorMagnitudeRange[2];
207    double _vectorComponentRange[3][2];
[2612]208    Renderer *_renderer;
[2454]209    double _dataScale;
[2317]210
211    vtkSmartPointer<vtkLookupTable> _lut;
[2328]212    vtkSmartPointer<vtkActor> _linesActor;
[2317]213    vtkSmartPointer<vtkActor> _seedActor;
[2581]214    vtkSmartPointer<vtkDataSet> _seedMesh;
[3680]215    vtkSmartPointer<vtkAlgorithm> _mesher;
[2317]216    vtkSmartPointer<vtkStreamTracer> _streamTracer;
217    vtkSmartPointer<vtkPolyDataAlgorithm> _lineFilter;
218    vtkSmartPointer<vtkPolyDataMapper> _pdMapper;
219    vtkSmartPointer<vtkPolyDataMapper> _seedMapper;
220};
221
222}
223
224#endif
Note: See TracBrowser for help on using the repository browser.