source: vtkvis/trunk/LIC.cpp @ 6716

Last change on this file since 6716 was 5835, checked in by ldelgass, 9 years ago

Require VTK >= 6.0.0

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