21#include <dune/grid/common/rangegenerators.hh>
22#include <dune/grid/io/file/vtk/common.hh>
24#include <opm/common/ErrorMacros.hpp>
32namespace Opm::GridDataOutput {
34template <
class Gr
idView,
unsigned int partitions>
40 dimw_ = GridView::dimension;
45template <
class Gr
idView,
unsigned int partitions>
51 { return Dune::VTK::geometryType(cit.geometry().type()) == Dune::VTK::polyhedron; });
54template <
class Gr
idView,
unsigned int partitions>
71template <
class Gr
idView,
unsigned int partitions>
76 if (max_size < nvertices_) {
78 "Opm::GridDataOutput::writeGridPoints( T& x_inout, T& "
79 "y_inout, T& z_inout ) "
80 +
" Input objects max_size (" + std::to_string(max_size)
81 +
") is not sufficient to fit the nvertices_ values (" + std::to_string(nvertices_) +
")");
86 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
93 }
else if (dimw_ == 2) {
94 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
98 z_inout[i] =
static_cast<T
>(0.0);
107template <
class Gr
idView,
unsigned int partitions>
108template <
typename VectType>
118 if ((
check_size_x <
static_cast<std::size_t
>(nvertices_)) ||
120 (
check_size_z <
static_cast<std::size_t
>(nvertices_))) {
123 "Opm::GridDataOutput::writeGridPoints( VectType& x_inout, VectType& "
124 "y_inout, VectType& z_inout ) At least one of the inputs"
125 +
" object x size " + std::to_string(
check_size_x) +
" object y size "
127 +
" is not sufficient to fit the nvertices_ values( " + std::to_string(nvertices_) +
" )");
132 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
139 }
else if (dimw_ == 2) {
141 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
154template <
class Gr
idView,
unsigned int partitions>
159 if (max_size < nvertices_ * 3) {
160 assert(max_size >= nvertices_ * 3);
162 "Opm::GridDataOutput::writeGridPoints_AOS( T* xyz_inout ) " +
" Input objects max_size ("
163 + std::to_string(max_size) +
") is not sufficient to fit the nvertices_ * 3 values ("
164 + std::to_string(nvertices_ * 3) +
")");
169 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
175 }
else if (dimw_ == 2) {
176 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
186template <
class Gr
idView,
unsigned int partitions>
187template <
typename VectType>
195 if (
check_size <
static_cast<std::size_t
>(nvertices_ * 3)) {
198 "Opm::GridDataOutput::writeGridPoints_AOS( VectType& xyz_inout ) "
199 +
" Input objects check_size (" + std::to_string(
check_size)
200 +
") is not sufficient to fit the nvertices_ * 3 values (" + std::to_string(nvertices_ * 3)
206 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
212 }
else if (dimw_ == 2) {
214 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
224template <
class Gr
idView,
unsigned int partitions>
229 if (max_size < nvertices_ * 3) {
232 "Opm::GridDataOutput::writeGridPoints_SOA( T& xyz_inout ) " +
" Input objects max_size ("
233 + std::to_string(max_size) +
") is not sufficient to fit the nvertices_ * 3 values ("
234 + std::to_string(nvertices_ * 3) +
")");
243 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
250 }
else if (dimw_ == 2) {
251 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
262template <
class Gr
idView,
unsigned int partitions>
263template <
typename VectType>
269 if (
check_size <
static_cast<std::size_t
>(nvertices_ * 3)) {
272 "Opm::GridDataOutput::writeGridPoints_SOA( VectType& xyz_inout ) "
273 +
" Input objects check_size (" + std::to_string(
check_size)
274 +
") is not sufficient to fit the nvertices_ * 3 values (" + std::to_string(nvertices_ * 3)
286 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
293 }
else if (dimw_ == 2) {
295 for (
const auto&
vit :
vertices(gridView_, Dune::Partitions::all)) {
306template <
class Gr
idView,
unsigned int partitions>
307template <
typename Integer>
312 if (max_size < ncorners_) {
315 "Opm::GridDataOutput::writeConnectivity( T* connectivity_inout,... ) "
316 +
" Input max_size value (" + std::to_string(max_size)
317 +
") is not sufficient to fit the ncorners_ values (" + std::to_string(ncorners_) +
")");
323 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
326 const int vxIdx = gridView_.indexSet().subIndex(
cit,
vx, 3);
333 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
336 const int vxIdx = gridView_.indexSet().subIndex(
cit,
vx, 3);
347template <
class Gr
idView,
unsigned int partitions>
348template <
typename VectType>
355 if (
check_size <
static_cast<std::size_t
>(ncorners_)) {
358 "Opm::GridDataOutput::writeConnectivity( VectType& "
359 "connectivity_inout ) "
360 +
" Input objects size (" + std::to_string(
check_size)
361 +
") is not sufficient to fit the ncorners_ values (" + std::to_string(ncorners_) +
")");
367 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
370 const int vxIdx = gridView_.indexSet().subIndex(
cit,
vx, 3);
377 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
380 const int vxIdx = gridView_.indexSet().subIndex(
cit,
vx, 3);
391template <
class Gr
idView,
unsigned int partitions>
392template <
typename Integer>
396 if (max_size < ncells_) {
399 "Opm::GridDataOutput::writeOffsetsCells( T* offsets_inout ) " +
" Input objects max_size ("
400 + std::to_string(max_size) +
") is not sufficient to fit the ncells_ values ("
401 + std::to_string(ncells_) +
")");
405 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
413template <
class Gr
idView,
unsigned int partitions>
414template <
typename VectType>
419 if (
check_size <
static_cast<std::size_t
>(ncells_)) {
422 "Opm::GridDataOutput::writeOffsetsCells( VectType& "
424 +
" Input objects size (" + std::to_string(
offsets_inout.size())
425 +
") is not sufficient to fit the ncells_ values (" + std::to_string(ncells_) +
")");
432 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
440template <
class Gr
idView,
unsigned int partitions>
441template <
typename Integer>
445 if (max_size < ncells_) {
448 "Opm::GridDataOutput::writeCellTypes( T* types_inout ) " +
" Input objects max_size ("
449 + std::to_string(max_size) +
") is not sufficient to fit the ncells_ values ("
450 + std::to_string(ncells_) +
")");
453 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
460template <
class Gr
idView,
unsigned int partitions>
461template <
typename VectType>
467 if (
check_size <
static_cast<std::size_t
>(ncells_)) {
469 "Opm::GridDataOutput::writeCellTypes( VectType& types_inout ) " +
" Input objects check_size ("
470 + std::to_string(
check_size) +
") is not sufficient to fit the ncells_ values ("
471 + std::to_string(ncells_) +
")");
475 for (
const auto&
cit :
elements(gridView_, dunePartition_)) {
476 int vtktype =
static_cast<int>(Dune::VTK::geometryType(
cit.type()));
482template <
class Gr
idView,
unsigned int partitions>
486 if (this->dunePartition_ == Dune::Partitions::all) {
487 return "Dune::Partitions::all";
489 if (this->dunePartition_ == Dune::Partitions::interior) {
490 return "Dune::Partitions::interior";
492 if (this->dunePartition_ == Dune::Partitions::interiorBorder) {
493 return "Dune::Partitions::interiorBorder";
495 if (this->dunePartition_ == Dune::Partitions::interiorBorderOverlap) {
496 return "Dune::Partitions::interiorBorderOverlap";
498 if (this->dunePartition_ == Dune::Partitions::front) {
499 return "Dune::Partitions::front";
501 if (this->dunePartition_ == Dune::Partitions::interiorBorderOverlapFront) {
502 return "Dune::Partitions::InteriorBorderOverlapFront";
504 if (this->dunePartition_ == Dune::Partitions::border) {
505 return "Dune::Partitions::border";
507 if (this->dunePartition_ == Dune::Partitions::ghost) {
508 return "Dune::Partitions::ghost";
511 return "Unknown Dune::PartitionSet<>";
514template <
class Gr
idView,
unsigned int partitions>
515void SimMeshDataAccessor<GridView,partitions>::
516printGridDetails(std::ostream&
outstr)
const
518 outstr <<
"Dune Partition = " << partition_value_ <<
", " << getPartitionTypeString() << std::endl;
519 outstr <<
"ncells_: " << getNCells() << std::endl;
520 outstr <<
"nvertices_: " << getNVertices() << std::endl;
521 outstr <<
"ncorners_: " << getNCorners() << std::endl;
Allows model geometry data to be passed to external code - via a copy direct to input pointers.
ConnectivityVertexOrder
Allows selection of order of vertices in writeConnectivity()
Definition GridDataOutput.hpp:93
Definition GridDataOutput.hpp:97
long writeGridPoints_AOS(T *xyz_inout, long max_size=0) const
Write the positions of vertices - directly to the pointers given in parameters as Array of Structures...
Definition GridDataOutput_impl.hpp:157
SimMeshDataAccessor(const GridView &gridView, Dune::PartitionSet< partitions > dunePartition)
Construct a SimMeshDataAccessor working on a specific GridView and specialize to a Dune::PartitionSet...
Definition GridDataOutput_impl.hpp:36
void countEntities()
Count the vertices, cells and corners.
Definition GridDataOutput_impl.hpp:55
long writeGridPoints(T *x_inout, T *y_inout, T *z_inout, long max_size=0) const
Write the positions of vertices - directly to the pointers given in parameters.
Definition GridDataOutput_impl.hpp:74
long writeOffsetsCells(Integer *offsets_inout, long max_size=0) const
Write the offsets values - directly to the pointer given in parameter 1.
Definition GridDataOutput_impl.hpp:394
long writeConnectivity(Integer *connectivity_inout, ConnectivityVertexOrder whichOrder, long max_size=0) const
Write the connectivity array - directly to the pointer given in parameter 1 Reorders the indices as s...
Definition GridDataOutput_impl.hpp:309
bool polyhedralCellPresent() const
Checks for cells that have polyhedral type within the current partition of cells.
Definition GridDataOutput_impl.hpp:46
long writeCellTypes(Integer *types_inout, long max_size=0) const
Write the cell types values - directly to the pointer given in parameter 1.
Definition GridDataOutput_impl.hpp:443
long writeGridPoints_SOA(T *xyz_inout, long max_size=0) const
Write the positions of vertices - directly to the pointers given in parameters as Structure of Arrays...
Definition GridDataOutput_impl.hpp:227
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240