source: branches/nanovis2/packages/vizservers/vtkvis/RpVtkDataSet.cpp @ 3175

Last change on this file since 3175 was 3175, checked in by ldelgass, 12 years ago

sync with trunk

  • Property svn:eol-style set to native
File size: 30.0 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#include <cstring>
10#include <cfloat>
11#include <cmath>
12
13#include <vtkCharArray.h>
14#include <vtkDataSetReader.h>
15#include <vtkDataSetWriter.h>
16#include <vtkPolyData.h>
17#include <vtkStructuredPoints.h>
18#include <vtkStructuredGrid.h>
19#include <vtkRectilinearGrid.h>
20#include <vtkUnstructuredGrid.h>
21#include <vtkProperty.h>
22#include <vtkPointData.h>
23#include <vtkCellData.h>
24#include <vtkCell.h>
25#include <vtkLookupTable.h>
26
27#include "RpVtkDataSet.h"
28#include "Trace.h"
29
30using namespace Rappture::VtkVis;
31
32DataSet::DataSet(const std::string& name) :
33    _name(name),
34    _visible(true),
35    _opacity(1),
36    _cellSizeAverage(0)
37{
38    _cellSizeRange[0] = -1;
39    _cellSizeRange[1] = -1;
40}
41
42DataSet::~DataSet()
43{
44}
45
46/**
47 * \brief Create and initialize a VTK Prop to render the outline
48 */
49void DataSet::initProp()
50{
51    if (_prop == NULL) {
52        _prop = vtkSmartPointer<vtkActor>::New();
53        vtkProperty *property = _prop->GetProperty();
54        property->EdgeVisibilityOff();
55        property->SetOpacity(_opacity);
56        property->SetAmbient(.2);
57        property->LightingOff();
58        _prop->SetVisibility((_visible ? 1 : 0));
59    }
60}
61
62/**
63 * \brief Create and initialize a wireframe outline
64 */
65void DataSet::showOutline(bool state)
66{
67    if (state) {
68        if (_outlineFilter == NULL) {
69            _outlineFilter = vtkSmartPointer<vtkOutlineFilter>::New();
70            _outlineFilter->SetInput(_dataSet);
71        }
72        if (_outlineMapper == NULL) {
73            _outlineMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
74            _outlineMapper->SetInputConnection(_outlineFilter->GetOutputPort());
75        }
76        initProp();
77        _prop->SetMapper(_outlineMapper);
78    } else {
79        if (_prop != NULL) {
80            _prop->SetMapper(NULL);
81        }
82        if (_outlineMapper != NULL) {
83            _outlineMapper = NULL;
84        }
85        if (_outlineFilter != NULL) {
86            _outlineFilter = NULL;
87        }
88    }
89}
90
91/**
92 * \brief Set color of outline bounding box
93 */
94void DataSet::setOutlineColor(float color[3])
95{
96    if (_prop == NULL) {
97        initProp();
98    }
99    _prop->GetProperty()->SetColor(color[0], color[1], color[2]);
100}
101
102/**
103 * \brief Set opacity of DataSet outline
104 *
105 * This method is used for record-keeping and opacity of the
106 * DataSet bounds outline.  The renderer controls opacity
107 * of other related graphics objects.
108 */
109void DataSet::setOpacity(double opacity)
110{
111    _opacity = opacity;
112    if (_prop != NULL) {
113        _prop->GetProperty()->SetOpacity(opacity);
114    }
115}
116
117/**
118 * \brief Set visibility flag in DataSet
119 *
120 * This method is used for record-keeping and visibility of the
121 * DataSet bounds outline.  The renderer controls visibility
122 * of other related graphics objects.
123 */
124void DataSet::setVisibility(bool state)
125{
126    _visible = state;
127    if (_prop != NULL) {
128        _prop->SetVisibility((state ? 1 : 0));
129    }
130}
131
132/**
133 * \brief Get visibility flag in DataSet
134 *
135 * This method is used for record-keeping.  The renderer controls
136 * the visibility of related graphics objects.
137 */
138bool DataSet::getVisibility() const
139{
140    return _visible;
141}
142
143void DataSet::writeDataFile(const char *filename)
144{
145    if (_dataSet == NULL)
146        return;
147   
148    vtkSmartPointer<vtkDataSetWriter> writer = vtkSmartPointer<vtkDataSetWriter>::New();
149
150    writer->SetFileName(filename);
151    writer->SetInput(_dataSet);
152    writer->Write();
153}
154
155/**
156 * \brief Read a VTK data file
157 */
158bool DataSet::setDataFile(const char *filename)
159{
160    vtkSmartPointer<vtkDataSetReader> reader = vtkSmartPointer<vtkDataSetReader>::New();
161
162#if defined(WANT_TRACE) && defined(DEBUG)
163    reader->DebugOn();
164#endif
165
166    reader->SetFileName(filename);
167    reader->ReadAllNormalsOn();
168    //reader->ReadAllTCoordsOn();
169    reader->ReadAllScalarsOn();
170    //reader->ReadAllColorScalarsOn();
171    reader->ReadAllVectorsOn();
172    //reader->ReadAllTensorsOn();
173    reader->ReadAllFieldsOn();
174
175    return setData(reader);
176}
177
178/**
179 * \brief Read a VTK data file from a memory buffer
180 */
181bool DataSet::setData(char *data, int nbytes)
182{
183    vtkSmartPointer<vtkDataSetReader> reader = vtkSmartPointer<vtkDataSetReader>::New();
184    vtkSmartPointer<vtkCharArray> dataSetString = vtkSmartPointer<vtkCharArray>::New();
185
186#if defined(WANT_TRACE) && defined(DEBUG)
187    reader->DebugOn();
188    dataSetString->DebugOn();
189#endif
190
191    dataSetString->SetArray(data, nbytes, 1);
192    reader->SetInputArray(dataSetString);
193    reader->ReadFromInputStringOn();
194    reader->ReadAllNormalsOn();
195    //reader->ReadAllTCoordsOn();
196    reader->ReadAllScalarsOn();
197    //reader->ReadAllColorScalarsOn();
198    reader->ReadAllVectorsOn();
199    //reader->ReadAllTensorsOn();
200    reader->ReadAllFieldsOn();
201
202    return setData(reader);
203}
204
205void DataSet::print() const
206{
207    print(_dataSet);
208}
209
210void DataSet::print(vtkDataSet *ds)
211{
212    if (ds == NULL)
213        return;
214
215    TRACE("DataSet class: %s", ds->GetClassName());
216
217    TRACE("DataSet memory: %g MiB", ds->GetActualMemorySize()/1024.);
218
219    double bounds[6];
220    ds->GetBounds(bounds);
221
222    // Topology
223    TRACE("DataSet bounds: %g %g %g %g %g %g",
224          bounds[0], bounds[1],
225          bounds[2], bounds[3],
226          bounds[4], bounds[5]);
227    TRACE("Points: %d Cells: %d", ds->GetNumberOfPoints(), ds->GetNumberOfCells());
228
229    double dataRange[2];
230    if (ds->GetPointData() != NULL) {
231        TRACE("PointData arrays: %d", ds->GetPointData()->GetNumberOfArrays());
232        for (int i = 0; i < ds->GetPointData()->GetNumberOfArrays(); i++) {
233            if (ds->GetPointData()->GetArray(i) != NULL) {
234                ds->GetPointData()->GetArray(i)->GetRange(dataRange, -1);
235                TRACE("PointData[%d]: '%s' comp:%d, (%g,%g)", i,
236                      ds->GetPointData()->GetArrayName(i),
237                      ds->GetPointData()->GetAbstractArray(i)->GetNumberOfComponents(),
238                      dataRange[0], dataRange[1]);
239            } else {
240                TRACE("PointData[%d]: '%s' comp:%d", i,
241                      ds->GetPointData()->GetArrayName(i),
242                      ds->GetPointData()->GetAbstractArray(i)->GetNumberOfComponents());
243            }
244        }
245    }
246    if (ds->GetCellData() != NULL) {
247        TRACE("CellData arrays: %d", ds->GetCellData()->GetNumberOfArrays());
248        for (int i = 0; i < ds->GetCellData()->GetNumberOfArrays(); i++) {
249            if (ds->GetCellData()->GetArray(i) != NULL) {
250                ds->GetCellData()->GetArray(i)->GetRange(dataRange, -1);
251                TRACE("CellData[%d]: '%s' comp:%d, (%g,%g)", i,
252                      ds->GetCellData()->GetArrayName(i),
253                      ds->GetCellData()->GetAbstractArray(i)->GetNumberOfComponents(),
254                      dataRange[0], dataRange[1]);
255            } else {
256                TRACE("CellData[%d]: '%s' comp:%d", i,
257                      ds->GetCellData()->GetArrayName(i),
258                      ds->GetCellData()->GetAbstractArray(i)->GetNumberOfComponents());
259            }
260        }
261    }
262    if (ds->GetFieldData() != NULL) {
263        TRACE("FieldData arrays: %d", ds->GetFieldData()->GetNumberOfArrays());
264        for (int i = 0; i < ds->GetFieldData()->GetNumberOfArrays(); i++) {
265            if (ds->GetFieldData()->GetArray(i) != NULL) {
266                ds->GetFieldData()->GetArray(i)->GetRange(dataRange, -1);
267                TRACE("FieldData[%d]: '%s' comp:%d, tuples:%d (%g,%g)", i,
268                      ds->GetFieldData()->GetArrayName(i),
269                      ds->GetFieldData()->GetAbstractArray(i)->GetNumberOfComponents(),
270                      ds->GetFieldData()->GetAbstractArray(i)->GetNumberOfTuples(),
271                      dataRange[0], dataRange[1]);
272            } else {
273                TRACE("FieldData[%d]: '%s' comp:%d, tuples:%d", i,
274                      ds->GetFieldData()->GetArrayName(i),
275                      ds->GetFieldData()->GetAbstractArray(i)->GetNumberOfComponents(),
276                      ds->GetFieldData()->GetAbstractArray(i)->GetNumberOfTuples());
277            }
278        }
279    }
280}
281
282void DataSet::setDefaultArrays()
283{
284    if (_dataSet->GetPointData() != NULL &&
285        _dataSet->GetPointData()->GetScalars() == NULL &&
286        _dataSet->GetPointData()->GetNumberOfArrays() > 0) {
287        for (int i = 0; i < _dataSet->GetPointData()->GetNumberOfArrays(); i++) {
288            if (_dataSet->GetPointData()->GetArray(i) != NULL &&
289                _dataSet->GetPointData()->GetArray(i)->GetNumberOfComponents() == 1) {
290                TRACE("Setting point scalars to '%s'", _dataSet->GetPointData()->GetArrayName(i));
291                _dataSet->GetPointData()->SetActiveScalars(_dataSet->GetPointData()->GetArrayName(i));
292                break;
293            }
294        }
295    }
296    if (_dataSet->GetPointData() != NULL &&
297        _dataSet->GetPointData()->GetVectors() == NULL &&
298        _dataSet->GetPointData()->GetNumberOfArrays() > 0) {
299        for (int i = 0; i < _dataSet->GetPointData()->GetNumberOfArrays(); i++) {
300            if (_dataSet->GetPointData()->GetArray(i) != NULL &&
301                _dataSet->GetPointData()->GetArray(i)->GetNumberOfComponents() == 3) {
302                TRACE("Setting point vectors to '%s'", _dataSet->GetPointData()->GetArrayName(i));
303                _dataSet->GetPointData()->SetActiveVectors(_dataSet->GetPointData()->GetArrayName(i));
304                break;
305            }
306        }
307    }
308    if (_dataSet->GetCellData() != NULL &&
309        _dataSet->GetCellData()->GetScalars() == NULL &&
310        _dataSet->GetCellData()->GetNumberOfArrays() > 0) {
311        for (int i = 0; i < _dataSet->GetCellData()->GetNumberOfArrays(); i++) {
312            if (_dataSet->GetCellData()->GetArray(i) != NULL &&
313                _dataSet->GetCellData()->GetArray(i)->GetNumberOfComponents() == 1) {
314                TRACE("Setting cell scalars to '%s'", _dataSet->GetCellData()->GetArrayName(i));
315                _dataSet->GetCellData()->SetActiveScalars(_dataSet->GetCellData()->GetArrayName(i));
316                break;
317            }
318        }
319    }
320    if (_dataSet->GetCellData() != NULL &&
321        _dataSet->GetCellData()->GetVectors() == NULL &&
322        _dataSet->GetCellData()->GetNumberOfArrays() > 0) {
323        for (int i = 0; i < _dataSet->GetCellData()->GetNumberOfArrays(); i++) {
324            if (_dataSet->GetCellData()->GetArray(i) != NULL &&
325                _dataSet->GetCellData()->GetArray(i)->GetNumberOfComponents() == 3) {
326                TRACE("Setting cell vectors to '%s'", _dataSet->GetCellData()->GetArrayName(i));
327                _dataSet->GetCellData()->SetActiveVectors(_dataSet->GetCellData()->GetArrayName(i));
328                break;
329            }
330        }
331    }
332}
333
334/**
335 * \brief Read dataset using supplied reader
336 *
337 * Pipeline information is removed from the resulting
338 * vtkDataSet, so that the reader and its data can be
339 * released
340 */
341bool DataSet::setData(vtkDataSetReader *reader)
342{
343    TRACE("Enter");
344    // Force reading data set
345    reader->SetLookupTableName("");
346    reader->Update();
347
348    _dataSet = reader->GetOutput();
349    _dataSet->SetPipelineInformation(NULL);
350
351    if (_dataSet->GetPointData() != NULL &&
352        _dataSet->GetPointData()->GetScalars() != NULL &&
353        _dataSet->GetPointData()->GetScalars()->GetLookupTable() != NULL) {
354        ERROR("No lookup table should be specified in DataSets");
355    }
356
357    setDefaultArrays();
358
359#ifdef WANT_TRACE
360    print();
361#endif
362    TRACE("Leave");
363    return true;
364}
365
366/**
367 * \brief Set DataSet from existing vtkDataSet object
368 *
369 * Pipeline information is removed from the supplied vtkDataSet
370 */
371bool DataSet::setData(vtkDataSet *ds)
372{
373    _dataSet = ds;
374    _dataSet->SetPipelineInformation(NULL);
375
376    if (_dataSet->GetPointData() != NULL &&
377        _dataSet->GetPointData()->GetScalars() != NULL &&
378        _dataSet->GetPointData()->GetScalars()->GetLookupTable() != NULL) {
379        ERROR("No lookup table should be specified in DataSets");
380    }
381
382    setDefaultArrays();
383
384#ifdef WANT_TRACE
385    print();
386#endif
387    return true;
388}
389
390/**
391 * \brief Copy an existing vtkDataSet object
392 *
393 * Pipeline information is not copied from the supplied vtkDataSet
394 * into this DataSet, but pipeline information in the supplied
395 * vtkDataSet is not removed
396 */
397vtkDataSet *DataSet::copyData(vtkDataSet *ds)
398{
399    if (vtkPolyData::SafeDownCast(ds) != NULL) {
400        _dataSet = vtkSmartPointer<vtkPolyData>::New();
401        _dataSet->DeepCopy(vtkPolyData::SafeDownCast(ds));
402    } else if (vtkStructuredPoints::SafeDownCast(ds) != NULL) {
403        _dataSet = vtkSmartPointer<vtkStructuredPoints>::New();
404        _dataSet->DeepCopy(vtkStructuredPoints::SafeDownCast(ds));
405    } else if (vtkStructuredGrid::SafeDownCast(ds) != NULL) {
406        _dataSet = vtkSmartPointer<vtkStructuredGrid>::New();
407        _dataSet->DeepCopy(vtkStructuredGrid::SafeDownCast(ds));
408    } else if (vtkRectilinearGrid::SafeDownCast(ds) != NULL) {
409        _dataSet = vtkSmartPointer<vtkRectilinearGrid>::New();
410        _dataSet->DeepCopy(vtkRectilinearGrid::SafeDownCast(ds));
411    } else if (vtkUnstructuredGrid::SafeDownCast(ds) != NULL) {
412        _dataSet = vtkSmartPointer<vtkUnstructuredGrid>::New();
413        _dataSet->DeepCopy(vtkUnstructuredGrid::SafeDownCast(ds));
414    } else {
415        ERROR("Unknown data type");
416        return NULL;
417    }
418
419#ifdef WANT_TRACE
420    print();
421#endif   
422    return _dataSet;
423}
424
425/**
426 * \brief Does DataSet lie in the XY plane (z = 0)
427 */
428bool DataSet::isXY() const
429{
430    double bounds[6];
431    getBounds(bounds);
432    return (bounds[4] == 0. && bounds[4] == bounds[5]);
433}
434
435/**
436 * \brief Returns the dimensionality of the AABB
437 */
438int DataSet::numDimensions() const
439{
440    double bounds[6];
441    getBounds(bounds);
442    int numDims = 0;
443    if (bounds[0] != bounds[1])
444        numDims++;
445    if (bounds[2] != bounds[3])
446        numDims++;
447    if (bounds[4] != bounds[5])
448        numDims++;
449
450    return numDims;
451}
452
453/**
454 * \brief Determines if DataSet lies in a principal axis plane
455 * and if so, returns the plane normal and offset from origin
456 */
457bool DataSet::is2D(PrincipalPlane *plane, double *offset) const
458{
459    double bounds[6];
460    getBounds(bounds);
461    if (bounds[4] == bounds[5]) {
462        // Z = 0, XY plane
463        if (plane != NULL) {
464            *plane = PLANE_XY;
465        }
466        if (offset != NULL)
467            *offset = bounds[4];
468        return true;
469    } else if (bounds[0] == bounds[1]) {
470        // X = 0, ZY plane
471        if (plane != NULL) {
472            *plane = PLANE_ZY;
473        }
474        if (offset != NULL)
475            *offset = bounds[0];
476        return true;
477    } else if (bounds[2] == bounds[3]) {
478        // Y = 0, XZ plane
479        if (plane != NULL) {
480            *plane = PLANE_XZ;
481         }
482        if (offset != NULL)
483            *offset = bounds[2];
484        return true;
485    }
486    return false;
487}
488
489/**
490 * \brief Determines a principal plane with the
491 * largest two dimensions of the AABB
492 */
493PrincipalPlane DataSet::principalPlane() const
494{
495    double bounds[6];
496    getBounds(bounds);
497    double xlen = bounds[1] - bounds[0];
498    double ylen = bounds[3] - bounds[2];
499    double zlen = bounds[5] - bounds[4];
500    if (zlen <= xlen && zlen <= ylen) {
501        return PLANE_XY;
502    } else if (xlen <= ylen && xlen <= zlen) {
503        return PLANE_ZY;
504    } else {
505        return PLANE_XZ;
506    }
507}
508
509/**
510 * \brief Get the name/id of this dataset
511 */
512const std::string& DataSet::getName() const
513{
514    return _name;
515}
516
517/**
518 * \brief Get the underlying VTK DataSet object
519 */
520vtkDataSet *DataSet::getVtkDataSet()
521{
522    return _dataSet;
523}
524
525/**
526 * \brief Get the underlying VTK DataSet subclass class name
527 */
528const char *DataSet::getVtkType() const
529{
530    return _dataSet->GetClassName();
531}
532
533/**
534 * \brief Set the ative scalar array to the named field
535 */
536bool DataSet::setActiveScalars(const char *name)
537{
538    bool found = false;
539    if (_dataSet != NULL) {
540        if (_dataSet->GetPointData() != NULL) {
541            if (_dataSet->GetPointData()->SetActiveScalars(name) >= 0) {
542                TRACE("Set active point data scalars to %s", name);
543                found = true;
544            }
545        }
546        if (_dataSet->GetCellData() != NULL) {
547            if (_dataSet->GetCellData()->SetActiveScalars(name) >= 0) {
548                TRACE("Set active cell data scalars to %s", name);
549                found = true;
550            }
551        }
552    }
553#ifdef WANT_TRACE
554    if (_dataSet->GetPointData() != NULL) {
555        if (_dataSet->GetPointData()->GetScalars() != NULL) {
556            TRACE("Point data scalars: %s", _dataSet->GetPointData()->GetScalars()->GetName());
557        } else {
558            TRACE("NULL point data scalars");
559        }
560    }
561    if (_dataSet->GetCellData() != NULL) {
562        if (_dataSet->GetCellData()->GetScalars() != NULL) {
563            TRACE("Cell data scalars: %s", _dataSet->GetCellData()->GetScalars()->GetName());
564        } else {
565            TRACE("NULL cell data scalars");
566        }
567    }
568#endif
569    return found;
570}
571
572/**
573 * \brief Get the active scalar array field name
574 */
575const char *DataSet::getActiveScalarsName() const
576{
577    if (_dataSet != NULL) {
578         if (_dataSet->GetPointData() != NULL &&
579             _dataSet->GetPointData()->GetScalars() != NULL) {
580            return _dataSet->GetPointData()->GetScalars()->GetName();
581        }
582        TRACE("No point scalars");
583        if (_dataSet->GetCellData() != NULL &&
584            _dataSet->GetCellData()->GetScalars() != NULL) {
585            return _dataSet->GetCellData()->GetScalars()->GetName();
586        }
587        TRACE("No cell scalars");
588    }
589    return NULL;
590}
591
592/**
593 * \brief Get the active scalar array default attribute type
594 */
595DataSet::DataAttributeType DataSet::getActiveScalarsType() const
596{
597    if (_dataSet != NULL) {
598         if (_dataSet->GetPointData() != NULL &&
599             _dataSet->GetPointData()->GetScalars() != NULL) {
600            return POINT_DATA;
601        }
602        TRACE("No point scalars");
603        if (_dataSet->GetCellData() != NULL &&
604            _dataSet->GetCellData()->GetScalars() != NULL) {
605            return CELL_DATA;
606        }
607        TRACE("No cell scalars");
608    }
609    return POINT_DATA;
610}
611
612/**
613 * \brief Set the ative vector array to the named field
614 */
615bool DataSet::setActiveVectors(const char *name)
616{
617    bool found = false;
618    if (_dataSet != NULL) {
619        if (_dataSet->GetPointData() != NULL) {
620            if (_dataSet->GetPointData()->SetActiveVectors(name) >= 0) {
621                TRACE("Set active point data vectors to %s", name);
622                found = true;
623            }
624        }
625        if (_dataSet->GetCellData() != NULL) {
626            if (_dataSet->GetCellData()->SetActiveVectors(name) >= 0) {
627                TRACE("Set active cell data vectors to %s", name);
628                found = true;
629            }
630        }
631    }
632
633    return found;
634}
635
636/**
637 * \brief Get the active vector array default attribute type
638 */
639DataSet::DataAttributeType DataSet::getActiveVectorsType() const
640{
641    if (_dataSet != NULL) {
642         if (_dataSet->GetPointData() != NULL &&
643             _dataSet->GetPointData()->GetVectors() != NULL) {
644            return POINT_DATA;
645        }
646        TRACE("No point vectors");
647        if (_dataSet->GetCellData() != NULL &&
648            _dataSet->GetCellData()->GetVectors() != NULL) {
649            return CELL_DATA;
650        }
651        TRACE("No cell vectors");
652    }
653    return POINT_DATA;
654}
655
656/**
657 * \brief Get the active vector array field name
658 */
659const char *DataSet::getActiveVectorsName() const
660{
661    if (_dataSet != NULL) {
662        if (_dataSet->GetPointData() != NULL &&
663            _dataSet->GetPointData()->GetVectors() != NULL) {
664            return _dataSet->GetPointData()->GetVectors()->GetName();
665        }
666        TRACE("No point vectors");
667        if (_dataSet->GetCellData() != NULL &&
668            _dataSet->GetCellData()->GetVectors() != NULL) {
669            return _dataSet->GetCellData()->GetVectors()->GetName();
670        }
671        TRACE("No cell vectors");
672    }
673    return NULL;
674}
675
676bool DataSet::getFieldInfo(const char *fieldName,
677                           DataAttributeType *type,
678                           int *numComponents) const
679{
680    if (_dataSet->GetPointData() != NULL &&
681        _dataSet->GetPointData()->GetAbstractArray(fieldName) != NULL) {
682        if (type != NULL)
683            *type = POINT_DATA;
684        if (numComponents != NULL)
685            *numComponents = _dataSet->GetPointData()->GetAbstractArray(fieldName)->GetNumberOfComponents();
686        return true;
687    } else if (_dataSet->GetCellData() != NULL &&
688               _dataSet->GetCellData()->GetAbstractArray(fieldName) != NULL) {
689        if (type != NULL)
690            *type = CELL_DATA;
691        if (numComponents != NULL)
692            *numComponents = _dataSet->GetCellData()->GetAbstractArray(fieldName)->GetNumberOfComponents();
693        return true;
694    } else if (_dataSet->GetFieldData() != NULL &&
695               _dataSet->GetFieldData()->GetAbstractArray(fieldName) != NULL) {
696        if (type != NULL)
697            *type = FIELD_DATA;
698        if (numComponents != NULL)
699            *numComponents = _dataSet->GetFieldData()->GetAbstractArray(fieldName)->GetNumberOfComponents();
700        return true;
701    }
702    return false;
703}
704
705bool DataSet::getFieldInfo(const char *fieldName,
706                           DataAttributeType type,
707                           int *numComponents) const
708{
709    switch (type) {
710    case POINT_DATA:
711        if (_dataSet->GetPointData() != NULL &&
712            _dataSet->GetPointData()->GetAbstractArray(fieldName) != NULL) {
713            if (numComponents != NULL)
714                *numComponents = _dataSet->GetPointData()->GetAbstractArray(fieldName)->GetNumberOfComponents();
715            return true;
716        } else
717            return false;
718        break;
719    case CELL_DATA:
720        if (_dataSet->GetCellData() != NULL &&
721            _dataSet->GetCellData()->GetAbstractArray(fieldName) != NULL) {
722            if (numComponents != NULL)
723                *numComponents = _dataSet->GetCellData()->GetAbstractArray(fieldName)->GetNumberOfComponents();
724            return true;
725        } else
726            return false;
727        break;
728    case FIELD_DATA:
729        if (_dataSet->GetFieldData() != NULL &&
730            _dataSet->GetFieldData()->GetAbstractArray(fieldName) != NULL) {
731            if (numComponents != NULL)
732                *numComponents = _dataSet->GetFieldData()->GetAbstractArray(fieldName)->GetNumberOfComponents();
733            return true;
734        } else
735            return false;
736        break;
737    default:
738        ;
739    }
740    return false;
741}
742
743/**
744 * \brief Get the list of field names in the DataSet
745 *
746 * \param[in,out] names The field names will be appended to this list
747 * \param[in] type The DataAttributeType: e.g. POINT_DATA, CELL_DATA
748 * \param[in] numComponents Filter list by number of components, -1 means to
749 * return all fields regardless of dimension
750 */
751void DataSet::getFieldNames(std::vector<std::string>& names,
752                            DataAttributeType type, int numComponents) const
753{
754    if (_dataSet == NULL)
755        return;
756    switch (type) {
757    case POINT_DATA:
758        if (_dataSet->GetPointData() != NULL) {
759            for (int i = 0; i < _dataSet->GetPointData()->GetNumberOfArrays(); i++) {
760                if (numComponents == -1 ||
761                    (_dataSet->GetPointData()->GetAbstractArray(i) != NULL &&
762                     _dataSet->GetPointData()->GetAbstractArray(i)->GetNumberOfComponents() == numComponents)) {
763                    names.push_back(_dataSet->GetPointData()->GetArrayName(i));
764                }
765            }
766        }
767        break;
768    case CELL_DATA:
769        if (_dataSet->GetCellData() != NULL) {
770            for (int i = 0; i < _dataSet->GetCellData()->GetNumberOfArrays(); i++) {
771                if (numComponents == -1 ||
772                    (_dataSet->GetCellData()->GetAbstractArray(i) != NULL &&
773                     _dataSet->GetCellData()->GetAbstractArray(i)->GetNumberOfComponents() == numComponents)) {
774                    names.push_back(_dataSet->GetCellData()->GetArrayName(i));
775                }
776            }
777        }
778        break;
779    case FIELD_DATA:
780        if (_dataSet->GetFieldData() != NULL) {
781            for (int i = 0; i < _dataSet->GetFieldData()->GetNumberOfArrays(); i++) {
782                if (numComponents == -1 ||
783                    (_dataSet->GetFieldData()->GetAbstractArray(i) != NULL &&
784                     _dataSet->GetFieldData()->GetAbstractArray(i)->GetNumberOfComponents() == numComponents)) {
785                    names.push_back(_dataSet->GetFieldData()->GetArrayName(i));
786                }
787            }
788        }
789        break;
790    default:
791        ERROR("Unknown DataAttributeType %d", type);
792    }
793}
794
795/**
796 * \brief Get the range of scalar values in the DataSet
797 */
798void DataSet::getScalarRange(double minmax[2]) const
799{
800    if (_dataSet == NULL)
801        return;
802    if (_dataSet->GetPointData() != NULL &&
803        _dataSet->GetPointData()->GetScalars() != NULL) {
804        _dataSet->GetPointData()->GetScalars()->GetRange(minmax, -1);
805    } else if (_dataSet->GetCellData() != NULL &&
806               _dataSet->GetCellData()->GetScalars() != NULL) {
807        _dataSet->GetCellData()->GetScalars()->GetRange(minmax, -1);
808    }
809}
810
811/**
812 * \brief Get the range of a vector component in the DataSet
813 *
814 * \param[out] minmax The data range
815 * \param[in] component The field component, -1 means magnitude
816 */
817void DataSet::getVectorRange(double minmax[2], int component) const
818{
819    if (_dataSet == NULL)
820        return;
821    if (_dataSet->GetPointData() != NULL &&
822        _dataSet->GetPointData()->GetVectors() != NULL) {
823        _dataSet->GetPointData()->GetVectors()->GetRange(minmax, component);
824    } else if (_dataSet->GetCellData() != NULL &&
825               _dataSet->GetCellData()->GetVectors() != NULL) {
826        _dataSet->GetCellData()->GetVectors()->GetRange(minmax, component);
827    }
828}
829
830/**
831 * \brief Get the range of values for the named field in the DataSet
832 *
833 * \param[out] minmax The data range
834 * \param[in] fieldName The array name
835 * \param[in] type The DataAttributeType
836 * \param[in] component The field component, -1 means magnitude
837 * \return boolean indicating if field was found
838 */
839bool DataSet::getDataRange(double minmax[2], const char *fieldName,
840                           DataAttributeType type, int component) const
841{
842    if (_dataSet == NULL)
843        return false;
844    switch (type) {
845    case POINT_DATA:
846        if (_dataSet->GetPointData() != NULL &&
847            _dataSet->GetPointData()->GetArray(fieldName) != NULL) {
848            _dataSet->GetPointData()->GetArray(fieldName)->GetRange(minmax, component);
849            return true;
850        } else {
851            return false;
852        }
853        break;
854    case CELL_DATA:
855        if (_dataSet->GetCellData() != NULL &&
856            _dataSet->GetCellData()->GetArray(fieldName) != NULL) {
857            _dataSet->GetCellData()->GetArray(fieldName)->GetRange(minmax, component);
858            return true;
859        } else {
860            return false;
861        }
862        break;
863    case FIELD_DATA:
864        if (_dataSet->GetFieldData() != NULL &&
865            _dataSet->GetFieldData()->GetArray(fieldName) != NULL) {
866            _dataSet->GetFieldData()->GetArray(fieldName)->GetRange(minmax, component);
867            return true;
868        } else {
869            return false;
870        }
871        break;
872    default:
873        ERROR("Unknown DataAttributeType %d", type);
874        break;
875    }
876    return false;
877}
878
879/**
880 * \brief Get the bounds the DataSet
881 */
882void DataSet::getBounds(double bounds[6]) const
883{
884    _dataSet->GetBounds(bounds);
885}
886
887/**
888 * \brief Get the range of cell AABB diagonal lengths in the DataSet
889 */
890void DataSet::getCellSizeRange(double minmax[2], double *average)
891{
892    if (_dataSet == NULL ||
893        _dataSet->GetNumberOfCells() < 1) {
894        minmax[0] = 1;
895        minmax[1] = 1;
896        *average = 1;
897        return;
898    }
899
900    if (_cellSizeRange[0] >= 0.0) {
901        minmax[0] = _cellSizeRange[0];
902        minmax[1] = _cellSizeRange[1];
903        *average = _cellSizeAverage;
904        return;
905    }
906
907    minmax[0] = DBL_MAX;
908    minmax[1] = -DBL_MAX;
909
910    *average = 0;
911    for (int i = 0; i < _dataSet->GetNumberOfCells(); i++) {
912        double length2 = _dataSet->GetCell(i)->GetLength2();
913        if (length2 < minmax[0])
914            minmax[0] = length2;
915        if (length2 > minmax[1])
916            minmax[1] = length2;
917        *average += length2;
918    }
919    if (minmax[0] == DBL_MAX)
920        minmax[0] = 1;
921    if (minmax[1] == -DBL_MAX)
922        minmax[1] = 1;
923
924    minmax[0] = sqrt(minmax[0]);
925    minmax[1] = sqrt(minmax[1]);
926    *average = sqrt(*average/((double)_dataSet->GetNumberOfCells()));
927    _cellSizeRange[0] = minmax[0];
928    _cellSizeRange[1] = minmax[1];
929    _cellSizeAverage = *average;
930}
931
932/**
933 * \brief Get nearest data value given world coordinates x,y,z
934 *
935 * Note: no interpolation is performed on data
936 *
937 * \param[in] x World x coordinate to probe
938 * \param[in] y World y coordinate to probe
939 * \param[in] z World z coordinate to probe
940 * \param[out] value On success, contains the data value
941 * \return boolean indicating success or failure
942 */
943bool DataSet::getScalarValue(double x, double y, double z, double *value) const
944{
945    if (_dataSet == NULL)
946        return false;
947    if (_dataSet->GetPointData() == NULL ||
948        _dataSet->GetPointData()->GetScalars() == NULL) {
949        return false;
950    }
951    vtkIdType pt = _dataSet->FindPoint(x, y, z);
952    if (pt < 0)
953        return false;
954    *value = _dataSet->GetPointData()->GetScalars()->GetComponent(pt, 0);
955    return true;
956}
957
958/**
959 * \brief Get nearest vector data value given world coordinates x,y,z
960 *
961 * Note: no interpolation is performed on data
962 *
963 * \param[in] x World x coordinate to probe
964 * \param[in] y World y coordinate to probe
965 * \param[in] z World z coordinate to probe
966 * \param[out] vector On success, contains the data values
967 * \return boolean indicating success or failure
968 */
969bool DataSet::getVectorValue(double x, double y, double z, double vector[3]) const
970{
971    if (_dataSet == NULL)
972        return false;
973    if (_dataSet->GetPointData() == NULL ||
974        _dataSet->GetPointData()->GetVectors() == NULL) {
975        return false;
976    }
977    vtkIdType pt = _dataSet->FindPoint(x, y, z);
978    if (pt < 0)
979        return false;
980    assert(_dataSet->GetPointData()->GetVectors()->GetNumberOfComponents() == 3);
981    _dataSet->GetPointData()->GetVectors()->GetTuple(pt, vector);
982    return true;
983}
Note: See TracBrowser for help on using the repository browser.