source: vtkvis/trunk/DataSet.cpp @ 6372

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

Require VTK >= 6.0.0

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