Skip to content
Snippets Groups Projects
Commit 1dfcc606 authored by Stéphane Del Pino's avatar Stéphane Del Pino
Browse files

Improve VTKWriter

Add encapsulate friend function to build an Array from a Kokkos::View. This is
meant to be used exceptionnaly
parent 782991e2
No related branches found
No related tags found
1 merge request!11Feature/mpi
......@@ -2,6 +2,7 @@
#define CONNECTIVITY_MATRIX_HPP
#include <PastisUtils.hpp>
#include <Array.hpp>
#include <Kokkos_StaticCrsGraph.hpp>
class ConnectivityMatrix
......@@ -20,9 +21,11 @@ class ConnectivityMatrix
return m_is_built;
}
const auto& entries() const
auto entries() const
{
return m_host_matrix.entries;
return encapsulate(m_host_matrix.entries);
// using DataType = typename decltype(m_host_matrix.entries)::value_type;
// return Array<DataType>(m_host_matrix.entries);
}
const auto& rowsMap() const
......
......@@ -67,7 +67,7 @@ class ItemToItemMatrix
return m_connectivity_matrix.numEntries();
}
const auto& entries() const
auto entries() const
{
return m_connectivity_matrix.entries();
}
......
......@@ -279,14 +279,14 @@ class VTKWriter
using Rd = TinyVector<MeshType::dimension>;
const NodeValue<const Rd>& xr = mesh.xr();
Array<TinyVector<3>> positions(mesh.numberOfNodes());
for (NodeId r=0; r<mesh.numberOfNodes(); ++r) {
parallel_for(mesh.numberOfNodes(), PASTIS_LAMBDA(NodeId r) {
for (unsigned short i=0; i<MeshType::dimension; ++i) {
positions[r][i] = xr[r][i];
}
for (unsigned short i=MeshType::dimension; i<3; ++i) {
positions[r][i] = 0;
}
}
});
_write_array(fout, "Positions", positions);
}
fout << "</Points>\n";
......@@ -295,12 +295,8 @@ class VTKWriter
{
const auto& cell_to_node_matrix
= mesh.connectivity().cellToNodeMatrix();
Array<int> connectivity(cell_to_node_matrix.numberOfEntries());
for (size_t i=0; i<cell_to_node_matrix.numberOfEntries(); ++i) {
connectivity[i] = cell_to_node_matrix.entries()[i];
}
_write_array(fout, "connectivity", connectivity);
_write_array(fout, "connectivity", cell_to_node_matrix.entries());
}
{
......@@ -317,38 +313,45 @@ class VTKWriter
}
{
const auto& cell_to_node_matrix
= mesh.connectivity().cellToNodeMatrix();
Array<int8_t> types(mesh.numberOfCells());
for (CellId j=0; j<mesh.numberOfCells(); ++j) {
const auto& cell_nodes = cell_to_node_matrix[j];
switch (cell_nodes.size()) {
case 2: {
const auto& cell_type
= mesh.connectivity().cellType();
parallel_for(mesh.numberOfCells(), PASTIS_LAMBDA(const CellId& j) {
switch (cell_type[j]) {
case CellType::Line: {
types[j]=3;
break;
}
case 3: {
case CellType::Triangle: {
types[j]=5;
break;
}
case 4: {
if (mesh.meshDimension() == 3) {
types[j]=10;
} else {
case CellType::Quadrangle: {
types[j]=9;
break;
}
case CellType::Tetrahedron: {
types[j]=10;
break;
}
case 8: {
types[j]=12;
case CellType::Pyramid: {
types[j]=14;
break;
}
default: {
types[j]=7;
case CellType::Prism: {
types[j]=13;
break;
}
case CellType::Hexahedron: {
types[j]=12;
break;
}
default: {
std::cerr << __FILE__ << ':' << __LINE__ << ": unknown cell type\n";
std::exit(1);
}
}
});
_write_array(fout, "types", types);
}
......
......@@ -38,6 +38,10 @@ class Array
return image;
}
template <typename DataType2, typename ... RT>
friend PASTIS_INLINE
Array<DataType2> encapsulate(const Kokkos::View<DataType2*, RT...>& values);
PASTIS_INLINE
DataType& operator[](const index_type& i) const noexcept(NO_ASSERT)
{
......@@ -106,6 +110,15 @@ class Array
~Array() = default;
};
template <typename DataType, typename ... RT>
PASTIS_INLINE
Array<DataType> encapsulate(const Kokkos::View<DataType*, RT...>& values)
{
Array<DataType> array;
array.m_values = values;
return array;
}
// map, multimap, unordered_map and stack cannot be copied this way
template <typename Container>
PASTIS_INLINE
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment