source: vtkvis/trunk/DataSet.cpp

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

add include

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