source: branches/vtkvis_threaded/RpLIC.cpp @ 2524

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

Merge 2494:2522 from trunk

  • Property svn:eol-style set to native
File size: 13.5 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 <cassert>
9
10#include <vtkActor.h>
11#include <vtkProperty.h>
12#include <vtkPointData.h>
13#include <vtkCellData.h>
14#include <vtkCellDataToPointData.h>
15#include <vtkImageData.h>
16#include <vtkPolyData.h>
17#include <vtkDataSetMapper.h>
18#include <vtkPainterPolyDataMapper.h>
19#include <vtkDataSetSurfaceFilter.h>
20#include <vtkCutter.h>
21#include <vtkImageMask.h>
22#include <vtkImageCast.h>
23
24#include "RpLIC.h"
25#include "Trace.h"
26
27using namespace Rappture::VtkVis;
28
29LIC::LIC() :
30    VtkGraphicsObject(),
31    _sliceAxis(Z_AXIS),
32    _colorMap(NULL)
33{
34}
35
36LIC::~LIC()
37{
38}
39
40void LIC::initProp()
41{
42    if (_prop == NULL) {
43        _prop = vtkSmartPointer<vtkActor>::New();
44        vtkProperty *property = getActor()->GetProperty();
45        property->SetOpacity(_opacity);
46        property->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
47        property->SetLineWidth(_edgeWidth);
48        property->EdgeVisibilityOff();
49        property->SetAmbient(.2);
50        if (!_lighting)
51            property->LightingOff();
52    }
53}
54
55void LIC::update()
56{
57    if (_dataSet == NULL)
58        return;
59
60    vtkDataSet *ds = _dataSet->getVtkDataSet();
61
62    vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
63    if (ds->GetPointData() == NULL ||
64        ds->GetPointData()->GetVectors() == NULL) {
65         WARN("No vector point data in dataset %s", _dataSet->getName().c_str());
66         if (ds->GetCellData() != NULL &&
67             ds->GetCellData()->GetVectors() != NULL) {
68             cellToPtData =
69                 vtkSmartPointer<vtkCellDataToPointData>::New();
70             cellToPtData->SetInput(ds);
71             //cellToPtData->PassCellDataOn();
72             cellToPtData->Update();
73             ds = cellToPtData->GetOutput();
74         } else {
75             ERROR("No vector cell data in dataset %s", _dataSet->getName().c_str());
76             return;
77         }
78    }
79
80    vtkImageData *imageData = vtkImageData::SafeDownCast(ds);
81    if (imageData != NULL) {
82        if (_lic == NULL) {
83            _lic = vtkSmartPointer<vtkImageDataLIC2D>::New();
84        }
85        if (!_dataSet->is2D()) {
86            // 3D image/volume/uniform grid
87            if (_volumeSlicer == NULL)
88                _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
89            int dims[3];
90            imageData->GetDimensions(dims);
91            TRACE("Image data dimensions: %d %d %d", dims[0], dims[1], dims[2]);
92            _volumeSlicer->SetInput(ds);
93            _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, (dims[2]-1)/2, (dims[2]-1)/2);
94            _volumeSlicer->SetSampleRate(1, 1, 1);
95            _lic->SetInputConnection(_volumeSlicer->GetOutputPort());
96        } else {
97            // 2D image/volume/uniform grid
98            _lic->SetInput(ds);
99        }
100        if (_mapper == NULL) {
101            _mapper = vtkSmartPointer<vtkDataSetMapper>::New();
102        }
103        _mapper->SetInputConnection(_lic->GetOutputPort());
104    } else if (vtkPolyData::SafeDownCast(ds) == NULL) {
105        // structured grid, unstructured grid, or rectilinear grid
106        if (_lic == NULL) {
107            _lic = vtkSmartPointer<vtkImageDataLIC2D>::New();
108        }
109
110        // Need to convert to vtkImageData
111        double bounds[6];
112        ds->GetBounds(bounds);
113        double xSize = bounds[1] - bounds[0];
114        double ySize = bounds[3] - bounds[2];
115        double zSize = bounds[5] - bounds[4];
116        int minDir = 2;
117        if (xSize < ySize && xSize < zSize)
118            minDir = 0;
119        if (ySize < xSize && ySize < zSize)
120            minDir = 1;
121        if (_probeFilter == NULL)
122            _probeFilter = vtkSmartPointer<vtkProbeFilter>::New();
123
124        if (!_dataSet->is2D()) {
125            // Sample a plane within the grid bounding box
126            vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
127            if (_cutPlane == NULL) {
128                _cutPlane = vtkSmartPointer<vtkPlane>::New();
129            }
130            if (minDir == 0) {
131                _cutPlane->SetNormal(1, 0, 0);
132                _cutPlane->SetOrigin(bounds[0] + (bounds[1]-bounds[0])/2.,
133                                     0,
134                                     0);
135            } else if (minDir == 1) {
136                _cutPlane->SetNormal(0, 1, 0);
137                _cutPlane->SetOrigin(0,
138                                     bounds[2] + (bounds[3]-bounds[2])/2.,
139                                     0);
140            } else {
141                _cutPlane->SetNormal(0, 0, 1);
142                _cutPlane->SetOrigin(0,
143                                     0,
144                                     bounds[4] + (bounds[5]-bounds[4])/2.);
145            }
146            cutter->SetInput(ds);
147            cutter->SetCutFunction(_cutPlane);
148            _probeFilter->SetSourceConnection(cutter->GetOutputPort());
149        } else {
150            _probeFilter->SetSource(ds);
151        }
152
153        vtkSmartPointer<vtkImageData> imageData = vtkSmartPointer<vtkImageData>::New();
154        int xdim, ydim, zdim;
155        double origin[3], spacing[3];
156        if (minDir == 0) {
157            xdim = 1;
158            ydim = 128;
159            zdim = 128;
160            origin[0] = bounds[0] + xSize/2.;
161            origin[1] = bounds[2];
162            origin[2] = bounds[4];
163            spacing[0] = 0;
164            spacing[1] = ySize/((double)(ydim-1));
165            spacing[2] = zSize/((double)(zdim-1));
166            _sliceAxis = X_AXIS;
167        } else if (minDir == 1) {
168            xdim = 128;
169            ydim = 1;
170            zdim = 128;
171            origin[0] = bounds[0];
172            origin[1] = bounds[2] + ySize/2.;
173            origin[2] = bounds[4];
174            spacing[0] = xSize/((double)(xdim-1));
175            spacing[1] = 0;
176            spacing[2] = zSize/((double)(zdim-1));
177            _sliceAxis = Y_AXIS;
178        } else {
179            xdim = 128;
180            ydim = 128;
181            zdim = 1;
182            origin[0] = bounds[0];
183            origin[1] = bounds[2];
184            origin[2] = bounds[4] + zSize/2.;
185            spacing[0] = xSize/((double)(xdim-1));
186            spacing[1] = ySize/((double)(ydim-1));
187            spacing[2] = 0;
188            _sliceAxis = Z_AXIS;
189        }
190        imageData->SetDimensions(xdim, ydim, zdim);
191        imageData->SetOrigin(origin);
192        imageData->SetSpacing(spacing);
193        _probeFilter->SetInput(imageData);
194        _lic->SetInputConnection(_probeFilter->GetOutputPort());
195
196        if (_mapper == NULL) {
197            _mapper = vtkSmartPointer<vtkDataSetMapper>::New();
198            _mapper->SetColorModeToMapScalars();
199        }
200        _lic->Update();
201        vtkSmartPointer<vtkImageCast> cast = vtkSmartPointer<vtkImageCast>::New();
202        cast->SetInputConnection(_probeFilter->GetOutputPort());
203        cast->SetInputArrayToProcess(0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_POINTS,
204                                     _probeFilter->GetValidPointMaskArrayName());
205        cast->SetOutputScalarTypeToUnsignedChar();
206        vtkSmartPointer<vtkImageMask> mask = vtkSmartPointer<vtkImageMask>::New();
207        mask->SetInputConnection(0, _lic->GetOutputPort());
208        mask->SetInputConnection(1, cast->GetOutputPort());
209        _mapper->SetInputConnection(mask->GetOutputPort());
210    } else {
211        // DataSet is a PolyData
212        vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
213        assert(pd != NULL);
214
215        // XXX: This mapper seems to conflict with offscreen rendering, so the main
216        // renderer needs to be set not to use FBOs for this to work
217        _mapper = vtkSmartPointer<vtkPainterPolyDataMapper>::New();
218        vtkPainterPolyDataMapper *ppdmapper = vtkPainterPolyDataMapper::SafeDownCast(_mapper);
219        if (_painter == NULL) {
220            _painter = vtkSmartPointer<vtkSurfaceLICPainter>::New();
221            _painter->SetDelegatePainter(ppdmapper->GetPainter());
222         }
223        ppdmapper->SetPainter(_painter);
224        ppdmapper->SetInput(pd);
225    }
226
227    if (_lut == NULL) {
228        setColorMap(ColorMap::getGrayDefault());
229    }
230
231    initProp();
232    getActor()->SetMapper(_mapper);
233
234    _mapper->Update();
235#ifdef WANT_TRACE
236    if (_lic != NULL) {
237        TRACE("FBO status: %d LIC status: %d", _lic->GetFBOSuccess(), _lic->GetLICSuccess());
238    } else if (_painter != NULL) {
239        TRACE("LIC status: %d", _painter->GetLICSuccess());
240    }
241#endif
242}
243
244/**
245 * \brief Select a 2D slice plane from a 3D DataSet
246 *
247 * \param[in] axis Axis of slice plane
248 * \param[in] ratio Position [0,1] of slice plane along axis
249 */
250void LIC::selectVolumeSlice(Axis axis, double ratio)
251{
252    if (_dataSet->is2D()) {
253        WARN("DataSet not 3D, returning");
254        return;
255    }
256
257    if (_volumeSlicer == NULL &&
258        _probeFilter == NULL) {
259        WARN("Called before update() or DataSet is not a volume");
260        return;
261    }
262
263    int dims[3];
264
265    vtkImageData *imageData = vtkImageData::SafeDownCast(_dataSet->getVtkDataSet());
266    if (imageData == NULL) {
267        if (_probeFilter != NULL) {
268            imageData = vtkImageData::SafeDownCast(_probeFilter->GetInput());
269            if (imageData == NULL) {
270                ERROR("Couldn't get probe filter input image");
271                return;
272            }
273        } else {
274            ERROR("Not a volume data set");
275            return;
276        }
277    }
278    imageData->GetDimensions(dims);
279
280    _sliceAxis = axis;
281
282    if (_probeFilter != NULL) {
283        vtkImageData *imageData = vtkImageData::SafeDownCast(_probeFilter->GetInput());
284        double bounds[6];
285        assert(vtkDataSet::SafeDownCast(_probeFilter->GetSource()) != NULL);
286        _dataSet->getBounds(bounds);
287        int dim = 128;
288
289        switch (axis) {
290        case X_AXIS:
291            imageData->SetDimensions(1, dim, dim);
292            imageData->SetOrigin(bounds[0] + (bounds[1]-bounds[0])*ratio, bounds[2], bounds[4]);
293            imageData->SetSpacing(0,
294                                  (bounds[3]-bounds[2])/((double)(dim-1)),
295                                  (bounds[5]-bounds[4])/((double)(dim-1)));
296            if (_cutPlane != NULL) {
297                _cutPlane->SetNormal(1, 0, 0);
298                _cutPlane->SetOrigin(bounds[0] + (bounds[1]-bounds[0])*ratio, 0, 0);
299            }
300            break;
301        case Y_AXIS:
302            imageData->SetDimensions(dim, 1, dim);
303            imageData->SetOrigin(bounds[0], bounds[2] + (bounds[3]-bounds[2])*ratio, bounds[4]);
304            imageData->SetSpacing((bounds[1]-bounds[0])/((double)(dim-1)),
305                                  0,
306                                  (bounds[5]-bounds[4])/((double)(dim-1)));
307            if (_cutPlane != NULL) {
308                _cutPlane->SetNormal(0, 1, 0);
309                _cutPlane->SetOrigin(0, bounds[2] + (bounds[3]-bounds[2])*ratio, 0);
310            }
311            break;
312        case Z_AXIS:
313            imageData->SetDimensions(dim, dim, 1);
314            imageData->SetOrigin(bounds[0], bounds[2], bounds[4] + (bounds[5]-bounds[4])*ratio);
315            imageData->SetSpacing((bounds[1]-bounds[0])/((double)(dim-1)),
316                                  (bounds[3]-bounds[2])/((double)(dim-1)),
317                                  0);
318            if (_cutPlane != NULL) {
319                _cutPlane->SetNormal(0, 0, 1);
320                _cutPlane->SetOrigin(0, 0, bounds[4] + (bounds[5]-bounds[4])*ratio);
321            }
322            break;
323        default:
324            ERROR("Invalid Axis");
325            return;
326        }
327    } else {
328        int voi[6];
329
330        switch (axis) {
331        case X_AXIS:
332            voi[0] = voi[1] = (int)((dims[0]-1) * ratio);
333            voi[2] = 0;
334            voi[3] = dims[1]-1;
335            voi[4] = 0;
336            voi[5] = dims[2]-1;
337            break;
338        case Y_AXIS:
339            voi[2] = voi[3] = (int)((dims[1]-1) * ratio);
340            voi[0] = 0;
341            voi[1] = dims[0]-1;
342            voi[4] = 0;
343            voi[5] = dims[2]-1;
344            break;
345        case Z_AXIS:
346            voi[4] = voi[5] = (int)((dims[2]-1) * ratio);
347            voi[0] = 0;
348            voi[1] = dims[0]-1;
349            voi[2] = 0;
350            voi[3] = dims[1]-1;
351            break;
352        default:
353            ERROR("Invalid Axis");
354            return;
355        }
356
357        _volumeSlicer->SetVOI(voi);
358    }
359
360    if (_lic != NULL)
361        _lic->Update();
362
363    if (_mapper != NULL)
364        _mapper->Update();
365}
366
367/**
368 * \brief Called when the color map has been edited
369 */
370void LIC::updateColorMap()
371{
372    setColorMap(_colorMap);
373}
374
375/**
376 * \brief Associate a colormap lookup table with the DataSet
377 */
378void LIC::setColorMap(ColorMap *cmap)
379{
380    if (cmap == NULL)
381        return;
382
383    _colorMap = cmap;
384 
385    if (_lut == NULL) {
386        _lut = vtkSmartPointer<vtkLookupTable>::New();
387        if (_mapper != NULL) {
388            _mapper->UseLookupTableScalarRangeOff();
389            _mapper->SetLookupTable(_lut);
390        }
391    }
392
393    _lut->DeepCopy(cmap->getLookupTable());
394    _lut->SetRange(_dataRange);
395}
396
397void LIC::updateRanges(bool useCumulative,
398                       double scalarRange[2],
399                       double vectorMagnitudeRange[2],
400                       double vectorComponentRange[3][2])
401{
402    if (useCumulative) {
403        _dataRange[0] = scalarRange[0];
404        _dataRange[1] = scalarRange[1];
405    } else if (_dataSet != NULL) {
406        _dataSet->getScalarRange(_dataRange);
407    }
408
409    if (_lut != NULL) {
410        _lut->SetRange(_dataRange);
411    }
412}
413
414/**
415 * \brief Set a group of world coordinate planes to clip rendering
416 *
417 * Passing NULL for planes will remove all cliping planes
418 */
419void LIC::setClippingPlanes(vtkPlaneCollection *planes)
420{
421    if (_mapper != NULL) {
422        _mapper->SetClippingPlanes(planes);
423    }
424}
Note: See TracBrowser for help on using the repository browser.