source: trunk/packages/vizservers/vtkvis/RpStreamlines.h @ 2618

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

Refactor vtkvis to support setting colormap fields by name/attribute type
rather than always using active scalars/vectors. Also convert common
graphics objects set methods in Renderer to template methods and separate
core and graphics object related methods to separate files.

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