source: vtkvis/trunk/DataSet.cpp @ 6538

Last change on this file since 6538 was 6538, checked in by ldelgass, 8 years ago

Support reading/writing VTK XML format in vtkvis

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