swh:1:snp:f521c49ab17ef7db6ec70b2430e1ed203f50383f
Raw File
Tip revision: e8c101c3f8dfd53e6697096b77ad5ed99aa6fa9a authored by Lars Bilke on 20 July 2021, 07:25:50 UTC
Merge branch 'fix-qt' into 'master'
Tip revision: e8c101c
VtkMeshConverter.cpp
/**
 * \file
 * \author Karsten Rink
 * \date   2011-08-23
 * \brief  Implementation of the VtkMeshConverter class.
 *
 * \copyright
 * Copyright (c) 2012-2021, OpenGeoSys Community (http://www.opengeosys.org)
 *            Distributed under a Modified BSD License.
 *              See accompanying file LICENSE.txt or
 *              http://www.opengeosys.org/project/license
 *
 */

#include "VtkMeshConverter.h"

#include "MeshLib/Elements/Elements.h"
#include "MeshLib/Mesh.h"
#include "MeshLib/Node.h"

// Conversion from Image to QuadMesh
#include <vtkBitArray.h>
#include <vtkCharArray.h>
#include <vtkDoubleArray.h>
#include <vtkFloatArray.h>
#include <vtkImageData.h>
#include <vtkIntArray.h>
#include <vtkLongArray.h>
#include <vtkLongLongArray.h>
#include <vtkPointData.h>
#include <vtkSmartPointer.h>
#include <vtkUnsignedCharArray.h>
#include <vtkUnsignedLongArray.h>
#include <vtkUnsignedLongLongArray.h>

// Conversion from vtkUnstructuredGrid
#include <vtkCell.h>
#include <vtkCellData.h>
#include <vtkUnsignedIntArray.h>
#include <vtkUnstructuredGrid.h>

namespace MeshLib
{
namespace detail
{
template <class T_ELEMENT>
MeshLib::Element* createElementWithSameNodeOrder(
    const std::vector<MeshLib::Node*>& nodes, vtkIdList* const node_ids,
    std::size_t const element_id)
{
    auto** ele_nodes = new MeshLib::Node*[T_ELEMENT::n_all_nodes];
    for (unsigned k(0); k < T_ELEMENT::n_all_nodes; k++)
    {
        ele_nodes[k] = nodes[node_ids->GetId(k)];
    }
    return new T_ELEMENT(ele_nodes, element_id);
}
}  // namespace detail

MeshLib::Mesh* VtkMeshConverter::convertUnstructuredGrid(
    vtkUnstructuredGrid* grid, std::string const& mesh_name)
{
    if (!grid)
    {
        return nullptr;
    }

    // set mesh nodes
    const std::size_t nNodes = grid->GetPoints()->GetNumberOfPoints();
    std::vector<MeshLib::Node*> nodes(nNodes);
    double* coords = nullptr;
    for (std::size_t i = 0; i < nNodes; i++)
    {
        coords = grid->GetPoints()->GetPoint(i);
        nodes[i] = new MeshLib::Node(coords[0], coords[1], coords[2], i);
    }

    // set mesh elements
    const std::size_t nElems = grid->GetNumberOfCells();
    std::vector<MeshLib::Element*> elements(nElems);
    auto node_ids = vtkSmartPointer<vtkIdList>::New();
    for (std::size_t i = 0; i < nElems; i++)
    {
        MeshLib::Element* elem;
        grid->GetCellPoints(i, node_ids);

        int cell_type = grid->GetCellType(i);
        switch (cell_type)
        {
            case VTK_VERTEX:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Point>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_LINE:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Line>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_TRIANGLE:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Tri>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_QUAD:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Quad>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_PIXEL:
            {
                auto** quad_nodes = new MeshLib::Node*[4];
                quad_nodes[0] = nodes[node_ids->GetId(0)];
                quad_nodes[1] = nodes[node_ids->GetId(1)];
                quad_nodes[2] = nodes[node_ids->GetId(3)];
                quad_nodes[3] = nodes[node_ids->GetId(2)];
                elem = new MeshLib::Quad(quad_nodes, i);
                break;
            }
            case VTK_TETRA:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Tet>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_HEXAHEDRON:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Hex>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_VOXEL:
            {
                auto** voxel_nodes = new MeshLib::Node*[8];
                voxel_nodes[0] = nodes[node_ids->GetId(0)];
                voxel_nodes[1] = nodes[node_ids->GetId(1)];
                voxel_nodes[2] = nodes[node_ids->GetId(3)];
                voxel_nodes[3] = nodes[node_ids->GetId(2)];
                voxel_nodes[4] = nodes[node_ids->GetId(4)];
                voxel_nodes[5] = nodes[node_ids->GetId(5)];
                voxel_nodes[6] = nodes[node_ids->GetId(7)];
                voxel_nodes[7] = nodes[node_ids->GetId(6)];
                elem = new MeshLib::Hex(voxel_nodes, i);
                break;
            }
            case VTK_PYRAMID:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Pyramid>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_WEDGE:
            {
                auto** prism_nodes = new MeshLib::Node*[6];
                for (unsigned j = 0; j < 3; ++j)
                {
                    prism_nodes[j] = nodes[node_ids->GetId(j + 3)];
                    prism_nodes[j + 3] = nodes[node_ids->GetId(j)];
                }
                elem = new MeshLib::Prism(prism_nodes, i);
                break;
            }
            case VTK_QUADRATIC_EDGE:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Line3>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_QUADRATIC_TRIANGLE:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Tri6>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_QUADRATIC_QUAD:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Quad8>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_BIQUADRATIC_QUAD:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Quad9>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_QUADRATIC_TETRA:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Tet10>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_QUADRATIC_HEXAHEDRON:
            {
                elem = detail::createElementWithSameNodeOrder<MeshLib::Hex20>(
                    nodes, node_ids, i);
                break;
            }
            case VTK_QUADRATIC_PYRAMID:
            {
                elem =
                    detail::createElementWithSameNodeOrder<MeshLib::Pyramid13>(
                        nodes, node_ids, i);
                break;
            }
            case VTK_QUADRATIC_WEDGE:
            {
                auto** prism_nodes = new MeshLib::Node*[15];
                for (unsigned j = 0; j < 3; ++j)
                {
                    prism_nodes[j] = nodes[node_ids->GetId(j + 3)];
                    prism_nodes[j + 3] = nodes[node_ids->GetId(j)];
                }
                for (unsigned j = 0; j < 3; ++j)
                {
                    prism_nodes[6 + j] = nodes[node_ids->GetId(8 - j)];
                }
                prism_nodes[9] = nodes[node_ids->GetId(12)];
                prism_nodes[10] = nodes[node_ids->GetId(14)];
                prism_nodes[11] = nodes[node_ids->GetId(13)];
                for (unsigned j = 0; j < 3; ++j)
                {
                    prism_nodes[12 + j] = nodes[node_ids->GetId(11 - j)];
                }
                elem = new MeshLib::Prism15(prism_nodes, i);
                break;
            }
            default:
                ERR("VtkMeshConverter::convertUnstructuredGrid(): Unknown mesh "
                    "element type '{:d}'.",
                    cell_type);
                return nullptr;
        }

        elements[i] = elem;
    }

    MeshLib::Mesh* mesh = new MeshLib::Mesh(mesh_name, nodes, elements);
    convertScalarArrays(*grid, *mesh);

    return mesh;
}

void VtkMeshConverter::convertScalarArrays(vtkUnstructuredGrid& grid,
                                           MeshLib::Mesh& mesh)
{
    vtkPointData* point_data = grid.GetPointData();
    auto const n_point_arrays =
        static_cast<unsigned>(point_data->GetNumberOfArrays());
    for (unsigned i = 0; i < n_point_arrays; ++i)
    {
        convertArray(*point_data->GetArray(i),
                     mesh.getProperties(),
                     MeshLib::MeshItemType::Node);
    }

    vtkCellData* cell_data = grid.GetCellData();
    auto const n_cell_arrays =
        static_cast<unsigned>(cell_data->GetNumberOfArrays());
    for (unsigned i = 0; i < n_cell_arrays; ++i)
    {
        convertArray(*cell_data->GetArray(i),
                     mesh.getProperties(),
                     MeshLib::MeshItemType::Cell);
    }

    vtkFieldData* field_data = grid.GetFieldData();
    auto const n_field_arrays =
        static_cast<unsigned>(field_data->GetNumberOfArrays());
    for (unsigned i = 0; i < n_field_arrays; ++i)
    {
        convertArray(
            *vtkDataArray::SafeDownCast(field_data->GetAbstractArray(i)),
            mesh.getProperties(),
            MeshLib::MeshItemType::IntegrationPoint);
    }
}

void VtkMeshConverter::convertArray(vtkDataArray& array,
                                    MeshLib::Properties& properties,
                                    MeshLib::MeshItemType type)
{
    if (vtkDoubleArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<double>(array, properties, type);
        return;
    }

    if (vtkFloatArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<float>(array, properties, type);
        return;
    }

    if (vtkIntArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<int>(array, properties, type);
        return;
    }

    if (vtkBitArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<bool>(array, properties, type);
        return;
    }

    if (vtkCharArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<char>(array, properties, type);
        return;
    }

    if (vtkLongArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<long>(array, properties, type);
        return;
    }

    if (vtkLongLongArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<long long>(array, properties, type);
        return;
    }

    if (vtkUnsignedLongArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<unsigned long>(array, properties,
                                                           type);
        return;
    }

    if (vtkUnsignedLongLongArray::SafeDownCast(&array))
    {
        VtkMeshConverter::convertTypedArray<unsigned long long>(
            array, properties, type);
        return;
    }

    if (vtkUnsignedIntArray::SafeDownCast(&array))
    {
        // MaterialIDs are assumed to be integers
        if (std::strncmp(array.GetName(), "MaterialIDs", 11) == 0)
        {
            VtkMeshConverter::convertTypedArray<int>(array, properties, type);
        }
        else
        {
            VtkMeshConverter::convertTypedArray<unsigned>(array, properties,
                                                          type);
        }

        return;
    }

    WARN(
        "Array '{:s}' in VTU file uses unsupported data type '{:s}'. The data "
        "array will not be available.",
        array.GetName(), array.GetDataTypeAsString());
}

}  // end namespace MeshLib
back to top