[cig-commits] r12161 - cs/benchmark/cigma/trunk/src

luis at geodynamics.org luis at geodynamics.org
Wed Jun 11 15:37:29 PDT 2008


Author: luis
Date: 2008-06-11 15:37:29 -0700 (Wed, 11 Jun 2008)
New Revision: 12161

Added:
   cs/benchmark/cigma/trunk/src/VtkReader2.cpp
   cs/benchmark/cigma/trunk/src/VtkReader2.h
Modified:
   cs/benchmark/cigma/trunk/src/Makefile.am
   cs/benchmark/cigma/trunk/src/Reader.cpp
   cs/benchmark/cigma/trunk/src/VtkList.cpp
Log:
Alternate implementation of VtkReader
 * Supports both legacy and xml vtk formats (vtk, vtu, vts, vtr)
 * Support for xml parallel vtk formats (pvtu, pvts, pvtr)


Modified: cs/benchmark/cigma/trunk/src/Makefile.am
===================================================================
--- cs/benchmark/cigma/trunk/src/Makefile.am	2008-06-11 22:37:27 UTC (rev 12160)
+++ cs/benchmark/cigma/trunk/src/Makefile.am	2008-06-11 22:37:29 UTC (rev 12161)
@@ -130,8 +130,8 @@
 	HdfReader.cpp \
 	TextReader.h \
 	TextReader.cpp \
-	VtkReader.h \
-	VtkReader.cpp \
+	VtkReader2.h \
+	VtkReader2.cpp \
 	VtkList.h \
 	VtkList.cpp \
 	\

Modified: cs/benchmark/cigma/trunk/src/Reader.cpp
===================================================================
--- cs/benchmark/cigma/trunk/src/Reader.cpp	2008-06-11 22:37:27 UTC (rev 12160)
+++ cs/benchmark/cigma/trunk/src/Reader.cpp	2008-06-11 22:37:29 UTC (rev 12161)
@@ -6,7 +6,7 @@
 #include "NullReader.h"
 #include "HdfReader.h"
 #include "TextReader.h"
-#include "VtkReader.h"
+#include "VtkReader2.h"
 
 #include "PathUtils.h"
 

Modified: cs/benchmark/cigma/trunk/src/VtkList.cpp
===================================================================
--- cs/benchmark/cigma/trunk/src/VtkList.cpp	2008-06-11 22:37:27 UTC (rev 12160)
+++ cs/benchmark/cigma/trunk/src/VtkList.cpp	2008-06-11 22:37:29 UTC (rev 12161)
@@ -2,7 +2,7 @@
 #include <cstdlib>
 
 #include "VtkList.h"
-#include "VtkReader.h"
+#include "VtkReader2.h"
 
 #include "vtkDataSet.h"
 #include "vtkPointData.h"

Added: cs/benchmark/cigma/trunk/src/VtkReader2.cpp
===================================================================
--- cs/benchmark/cigma/trunk/src/VtkReader2.cpp	                        (rev 0)
+++ cs/benchmark/cigma/trunk/src/VtkReader2.cpp	2008-06-11 22:37:29 UTC (rev 12161)
@@ -0,0 +1,346 @@
+#include <iostream>
+#include <string>
+#include <cstdlib>
+#include <cassert>
+
+#include "StringUtils.h"
+#include "PathUtils.h"
+
+#include "VtkReader2.h"
+
+#include "vtkDataSet.h"
+#include "vtkPointData.h"
+
+#include "vtkDataArray.h"
+#include "vtkDataArraySelection.h"
+
+#include "vtkIdList.h"
+#include "vtkCell.h"
+#include "vtkGenericCell.h"
+
+
+
+using namespace std;
+using namespace cigma;
+
+
+// ----------------------------------------------------------------------------
+
+VtkReader::VtkReader()
+{
+    dataset = 0;
+    legacy_reader = 0;
+    xml_reader = 0;
+    gridReaderType = NULL_GRID_READER;
+}
+
+VtkReader::~VtkReader()
+{
+    close();
+}
+
+// ----------------------------------------------------------------------------
+
+int VtkReader::open(const char *filename)
+{
+    //cout << ">> VtkReader::open(): Reading " << filename << endl; // XXX: need journal logging
+
+    string filepath = filename;
+    string fileroot, ext;
+    path_splitext(filepath, fileroot, ext);
+
+    if (!VtkExtension(ext.c_str()))
+    {
+        return -1;
+    }
+
+    int status;
+
+    if (ext == ".vtk")
+    {
+        legacy_reader = vtkDataSetReader::New();
+
+        legacy_reader->SetFileName(filename);
+
+        status = legacy_reader->OpenVTKFile();
+
+        if (status == 0)
+        {
+            //cerr << "Could not read VTK file " << filename << endl;
+            return -1;
+        }
+
+        gridReaderType = VTK;
+
+        legacy_reader->Update();
+
+        dataset = legacy_reader->GetOutput();
+        dataset->Update();
+
+    }
+    else
+    {
+        if (ext == ".vtu")
+        {
+            xml_reader = vtkXMLUnstructuredGridReader::New();
+            gridReaderType = VTU;
+        }
+        else if (ext == ".vts")
+        {
+            xml_reader = vtkXMLStructuredGridReader::New();
+            gridReaderType = VTS;
+        }
+        else if (ext == ".vtr")
+        {
+            xml_reader = vtkXMLRectilinearGridReader::New();
+            gridReaderType = VTR;
+        }
+        else if (ext == ".pvtu")
+        {
+            xml_reader = vtkXMLPUnstructuredGridReader::New();
+            gridReaderType = PVTU;
+        }
+        else if (ext == ".pvts")
+        {
+            xml_reader = vtkXMLPStructuredGridReader::New();
+            gridReaderType = PVTS;
+        }
+        else if (ext == ".pvtr")
+        {
+            xml_reader = vtkXMLRectilinearGridReader::New();
+            gridReaderType = PVTR;
+        }
+        else
+        {
+            cerr << "Unsupported VTK file (" << ext << ")" << endl;
+            return -1;
+        }
+
+        status = xml_reader->CanReadFile(filename);
+
+        if (status == 0)
+        {
+            cerr << "Could not read XML VTK file " << filename << endl;
+            gridReaderType = NULL_GRID_READER;
+            return -1;
+        }
+
+        xml_reader->SetFileName(filename);
+        xml_reader->Update();
+        xml_reader->GetPointDataArraySelection()->EnableAllArrays();
+
+        dataset = xml_reader->GetOutputAsDataSet();
+        dataset->Update();
+    }
+    return 0;
+}
+
+
+int VtkReader::close()
+{
+    gridReaderType = NULL_GRID_READER;
+
+    if (dataset != 0)
+    {
+        dataset = 0;
+    }
+
+    if (legacy_reader)
+    {
+        legacy_reader->Delete();
+        legacy_reader = 0;
+    }
+
+    if (xml_reader)
+    {
+        xml_reader->Delete();
+        xml_reader = 0;
+    }
+
+    return 0;
+}
+
+
+
+// ----------------------------------------------------------------------------
+
+static int GetDimensionFromCellType(int cellType)
+{
+    if ((cellType == VTK_TRIANGLE) ||
+        (cellType == VTK_QUADRATIC_TRIANGLE) ||
+        (cellType == VTK_PARAMETRIC_TRI_SURFACE) ||
+        (cellType == VTK_HIGHER_ORDER_TRIANGLE))
+    {
+        return 2;
+    }
+
+    if ((cellType == VTK_QUAD) ||
+        (cellType == VTK_QUADRATIC_QUAD) ||
+        (cellType == VTK_PARAMETRIC_QUAD_SURFACE) ||
+        (cellType == VTK_HIGHER_ORDER_QUAD))
+    {
+        return 2;
+    }
+
+    if ((cellType == VTK_TETRA) ||
+        (cellType == VTK_QUADRATIC_TETRA) ||
+        (cellType == VTK_PARAMETRIC_TETRA_REGION) ||
+        (cellType == VTK_HIGHER_ORDER_TETRAHEDRON))
+    {
+        return 3;
+    }
+
+    if ((cellType == VTK_HEXAHEDRON) ||
+        (cellType == VTK_QUADRATIC_HEXAHEDRON) ||
+        (cellType == VTK_PARAMETRIC_HEX_REGION) ||
+        (cellType == VTK_HIGHER_ORDER_HEXAHEDRON))
+    {
+        return 3;
+    }
+
+    return -1;
+}
+
+int VtkReader::get_dataset(const char *loc, double **data, int *num, int *dim)
+{
+    //cout << ">> VtkReader::get_dataset(): requesting location " << loc << endl;
+
+    int i,j;
+    int n_points;
+    int rank;
+    double *array = 0;
+
+    if (dataset == 0)
+    {
+        return -1;
+    }
+
+    vtkPointData *pointData = dataset->GetPointData();
+
+    assert(pointData != 0);
+
+    if (pointData->HasArray(loc) == 1)
+    {
+        vtkDataArray *dataArray = pointData->GetArray(loc);
+        assert(dataArray != 0);
+
+        /* XXX: verify that we are loading a dataset of the correct type
+        int dataType = dataArray->GetDataType();
+        if ((dataType != VTK_DOUBLE) || (dataType != VTK_FLOAT))
+        {
+            assert(false); // XXX
+            return -1;
+        }*/
+
+        n_points = dataArray->GetNumberOfTuples();
+        rank = dataArray->GetNumberOfComponents();
+        array = new double[n_points * rank];
+
+        for (i = 0; i < n_points; i++)
+        {
+            double dofval[rank];
+            dataArray->GetTuple(i, dofval);
+            for (j = 0; j < rank; j++)
+            {
+                array[i*rank + j] = dofval[j];
+            }
+        }
+    }
+    else
+    {
+        assert(false); // XXX
+        return -1;
+    }
+
+    *data = array;
+    *num = n_points;
+    *dim = rank;
+
+    return 0;
+}
+
+int VtkReader::get_coordinates(const char *loc, double **coordinates, int *nno, int *nsd)
+{
+    //cout << ">> VtkReader::get_coordinates(): requesting location " << loc << endl;
+
+    int i,j;
+    int n_points;
+    int n_dims;
+    double *array = 0;
+
+    if (dataset == 0)
+    {
+        return -1;
+    }
+
+    int cellType = dataset->GetCellType(0); // XXX: assumes uniform cell types
+
+    n_points = dataset->GetNumberOfPoints();
+    n_dims = GetDimensionFromCellType(cellType);
+
+    assert(n_dims <= 3);
+
+    if (n_dims < 0)
+    {
+        assert(false);
+        return -1;
+    }
+
+    array = new double[n_points * n_dims];
+
+    for (i = 0; i < n_points; i++)
+    {
+        double point[3];
+        dataset->GetPoint(i, point);
+        for (j = 0; j < n_dims; j++)
+        {
+            array[i*n_dims + j] = point[j];
+        }
+    }
+
+    *coordinates = array;
+    *nno = n_points;
+    *nsd = n_dims;
+
+    return 0;
+}
+
+int VtkReader::get_connectivity(const char *loc, int **connectivity, int *nel, int *ndofs)
+{
+    //cout << ">> VtkReader::get_connectivity(): requesting location " << loc << endl;
+
+    int i,j;
+    int n_cells;
+    int n_dofs;
+    int *array = 0;
+
+    if (dataset == 0)
+    {
+        return -1;
+    }
+
+    vtkCell *firstCell = dataset->GetCell(0);
+
+    n_cells = dataset->GetNumberOfCells();
+    n_dofs = firstCell->GetNumberOfPoints();    // XXX: assumes uniform cell types
+    array = new int[n_cells * n_dofs];
+
+    vtkGenericCell *cell = vtkGenericCell::New();
+    for (i = 0; i < n_cells; i++)
+    {
+        dataset->GetCell(i, cell);
+        for (j = 0; j < n_dofs; j++)
+        {
+            array[i*n_dofs + j] = cell->GetPointIds()->GetId(j);
+        }
+    }
+
+    *connectivity = array;
+    *nel = n_cells;
+    *ndofs = n_dofs;
+
+    return 0;
+}
+
+
+// ----------------------------------------------------------------------------

Added: cs/benchmark/cigma/trunk/src/VtkReader2.h
===================================================================
--- cs/benchmark/cigma/trunk/src/VtkReader2.h	                        (rev 0)
+++ cs/benchmark/cigma/trunk/src/VtkReader2.h	2008-06-11 22:37:29 UTC (rev 12161)
@@ -0,0 +1,62 @@
+#ifndef __VTK_READER2_H__
+#define __VTK_READER2_H__
+
+#include "Reader.h"
+
+#include "vtkDataSetReader.h"
+#include "vtkXMLDataReader.h"
+
+#include "vtkUnstructuredGridReader.h"
+#include "vtkStructuredGridReader.h"
+#include "vtkRectilinearGridReader.h"
+
+#include "vtkXMLUnstructuredGridReader.h"
+#include "vtkXMLStructuredGridReader.h"
+#include "vtkXMLRectilinearGridReader.h"
+
+#include "vtkXMLPUnstructuredGridReader.h"
+#include "vtkXMLPStructuredGridReader.h"
+#include "vtkXMLPRectilinearGridReader.h"
+
+
+namespace cigma
+{
+    class VtkReader;
+}
+
+
+class cigma::VtkReader : public Reader
+{
+public:
+    typedef enum {
+        NULL_GRID_READER, VTK, VTU, VTS, VTR, PVTU, PVTS, PVTR
+    } GridReaderType;
+
+public:
+    VtkReader();
+    ~VtkReader();
+
+public:
+    ReaderType getType() { return VTK_READER; }
+
+public:
+    int open(const char *filename);
+    int close();
+
+public:
+    int get_dataset(const char *loc, double **data, int *num, int *dim);
+    int get_coordinates(const char *loc, double **coordinates, int *nno, int *nsd);
+    int get_connectivity(const char *loc, int **connectivity, int *nel, int *ndofs);
+
+public:
+    GridReaderType gridReaderType;
+    vtkDataSetReader *legacy_reader;
+    vtkXMLReader *xml_reader;
+
+public:
+    vtkDataSet *dataset;
+
+};
+
+
+#endif



More information about the cig-commits mailing list