21#include <dune/grid/common/rangegenerators.hh>
22#include <dune/grid/io/file/vtk/common.hh>
24#include <opm/common/ErrorMacros.hpp>
34template <
class Gr
idView,
unsigned int partitions>
38 , dunePartition_(dunePartition)
40 dimw_ = GridView::dimension;
45template <
class Gr
idView,
unsigned int partitions>
48 const auto& elems = elements(gridView_, dunePartition_);
49 return std::any_of(elems.begin(), elems.end(),
51 { return Dune::VTK::geometryType(cit.geometry().type()) == Dune::VTK::polyhedron; });
54template <
class Gr
idView,
unsigned int partitions>
58 const auto& vert_partition = vertices(gridView_, Dune::Partitions::all);
59 nvertices_ = std::distance(vert_partition.begin(), vert_partition.end());
61 const auto& cell_partition = elements(gridView_, dunePartition_);
64 for (
const auto& cit : cell_partition) {
65 auto corner_geom = cit.geometry();
66 ncorners_ += corner_geom.corners();
71template <
class Gr
idView,
unsigned int partitions>
76 if (max_size < nvertices_) {
77 OPM_THROW(std::runtime_error,
78 "Opm::GridDataOutput::writeGridPoints( T& x_inout, T& "
79 "y_inout, T& z_inout ) "
81 +
") is not sufficient to fit the nvertices_ values (" +
std::to_string(nvertices_) +
")");
86 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
87 auto xyz_local = vit.geometry().corner(0);
88 x_inout[i] =
static_cast<T
>(xyz_local[0]);
89 y_inout[i] =
static_cast<T
>(xyz_local[1]);
90 z_inout[i] =
static_cast<T
>(xyz_local[2]);
93 }
else if (dimw_ == 2) {
94 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
95 auto xyz_local = vit.geometry().corner(0);
96 x_inout[i] =
static_cast<T
>(xyz_local[0]);
97 y_inout[i] =
static_cast<T
>(xyz_local[1]);
98 z_inout[i] =
static_cast<T
>(0.0);
102 assert(i == nvertices_);
107template <
class Gr
idView,
unsigned int partitions>
108template <
typename VectType>
110writeGridPoints(VectType& x_inout, VectType& y_inout, VectType& z_inout)
const
112 const std::size_t check_size_x = x_inout.size();
113 const std::size_t check_size_y = y_inout.size();
114 const std::size_t check_size_z = z_inout.size();
116 using VT =
decltype(x_inout.data()[0]);
118 if ((check_size_x <
static_cast<std::size_t
>(nvertices_)) ||
119 (check_size_y < static_cast<std::size_t>(nvertices_)) ||
120 (check_size_z <
static_cast<std::size_t
>(nvertices_))) {
122 OPM_THROW(std::runtime_error,
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)) {
133 auto xyz_local = vit.geometry().corner(0);
134 x_inout.data()[i] =
static_cast<VT
>(xyz_local[0]);
135 y_inout.data()[i] =
static_cast<VT
>(xyz_local[1]);
136 z_inout.data()[i] =
static_cast<VT
>(xyz_local[2]);
139 }
else if (dimw_ == 2) {
141 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
142 auto xyz_local = vit.geometry().corner(0);
143 x_inout.data()[i] =
static_cast<VT
>(xyz_local[0]);
144 y_inout.data()[i] =
static_cast<VT
>(xyz_local[1]);
145 z_inout.data()[i] =
static_cast<VT
>(td);
149 assert(i == nvertices_);
154template <
class Gr
idView,
unsigned int partitions>
159 if (max_size < nvertices_ * 3) {
160 assert(max_size >= nvertices_ * 3);
161 OPM_THROW(std::runtime_error,
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 ("
169 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
170 auto xyz_local = vit.geometry().corner(0);
171 xyz_inout[i++] =
static_cast<T
>(xyz_local[0]);
172 xyz_inout[i++] =
static_cast<T
>(xyz_local[1]);
173 xyz_inout[i++] =
static_cast<T
>(xyz_local[2]);
175 }
else if (dimw_ == 2) {
176 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
177 auto xyz_local = vit.geometry().corner(0);
178 xyz_inout[i++] =
static_cast<T
>(xyz_local[0]);
179 xyz_inout[i++] =
static_cast<T
>(xyz_local[1]);
180 xyz_inout[i++] =
static_cast<T
>(0.0);
186template <
class Gr
idView,
unsigned int partitions>
187template <
typename VectType>
191 const std::size_t check_size = xyz_inout.size();
193 using VT =
decltype(xyz_inout.data()[0]);
195 if (check_size <
static_cast<std::size_t
>(nvertices_ * 3)) {
196 assert(check_size >= nvertices_ * 3);
197 OPM_THROW(std::runtime_error,
198 "Opm::GridDataOutput::writeGridPoints_AOS( VectType& xyz_inout ) "
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)) {
207 auto xyz_local = vit.geometry().corner(0);
208 xyz_inout.data()[i++] =
static_cast<VT
>(xyz_local[0]);
209 xyz_inout[i++] =
static_cast<VT
>(xyz_local[1]);
210 xyz_inout[i++] =
static_cast<VT
>(xyz_local[2]);
212 }
else if (dimw_ == 2) {
214 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
215 auto xyz_local = vit.geometry().corner(0);
216 xyz_inout[i++] =
static_cast<VT
>(xyz_local[0]);
217 xyz_inout[i++] =
static_cast<VT
>(xyz_local[1]);
218 xyz_inout[i++] =
static_cast<VT
>(td);
224template <
class Gr
idView,
unsigned int partitions>
229 if (max_size < nvertices_ * 3) {
231 OPM_THROW(std::runtime_error,
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 ("
239 T* xyz_inout_y = xyz_inout + nvertices_;
240 T* xyz_inout_z = xyz_inout + (2 * nvertices_);
243 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
244 auto xyz_local = vit.geometry().corner(0);
245 xyz_inout[i] =
static_cast<T
>(xyz_local[0]);
246 xyz_inout_y[i] =
static_cast<T
>(xyz_local[1]);
247 xyz_inout_z[i] =
static_cast<T
>(xyz_local[2]);
250 }
else if (dimw_ == 2) {
251 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
252 auto xyz_local = vit.geometry().corner(0);
253 xyz_inout[i] =
static_cast<T
>(xyz_local[0]);
254 xyz_inout_y[i] =
static_cast<T
>(xyz_local[1]);
255 xyz_inout_z[i] =
static_cast<T
>(0.0);
262template <
class Gr
idView,
unsigned int partitions>
263template <
typename VectType>
267 const std::size_t check_size = xyz_inout.size();
269 if (check_size <
static_cast<std::size_t
>(nvertices_ * 3)) {
271 OPM_THROW(std::runtime_error,
272 "Opm::GridDataOutput::writeGridPoints_SOA( VectType& xyz_inout ) "
274 +
") is not sufficient to fit the nvertices_ * 3 values (" +
std::to_string(nvertices_ * 3)
278 using VT =
decltype(xyz_inout.data()[0]);
282 VT* xyz_inout_y = xyz_inout.data() + nvertices_;
283 VT* xyz_inout_z = xyz_inout.data() + (2 * nvertices_);
286 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
287 auto xyz_local = vit.geometry().corner(0);
288 xyz_inout.data()[i] =
static_cast<VT
>(xyz_local[0]);
289 xyz_inout_y[i] =
static_cast<VT
>(xyz_local[1]);
290 xyz_inout_z[i] =
static_cast<VT
>(xyz_local[2]);
293 }
else if (dimw_ == 2) {
295 for (
const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
296 auto xyz_local = vit.geometry().corner(0);
297 xyz_inout.data()[i] =
static_cast<VT
>(xyz_local[0]);
298 xyz_inout_y[i] =
static_cast<VT
>(xyz_local[1]);
299 xyz_inout_z[i] =
static_cast<VT
>(td);
306template <
class Gr
idView,
unsigned int partitions>
307template <
typename Integer>
312 if (max_size < ncorners_) {
314 OPM_THROW(std::runtime_error,
315 "Opm::GridDataOutput::writeConnectivity( T* connectivity_inout,... ) "
317 +
") is not sufficient to fit the ncorners_ values (" +
std::to_string(ncorners_) +
")");
321 if (whichOrder ==
DUNE) {
323 for (
const auto& cit : elements(gridView_, dunePartition_)) {
324 auto cell_corners = cit.geometry().corners();
325 for (
auto vx = 0; vx < cell_corners; ++vx) {
326 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
327 connectivity_inout[i + vx] = vxIdx;
333 for (
const auto& cit : elements(gridView_, dunePartition_)) {
334 auto cell_corners = cit.geometry().corners();
335 for (
auto vx = 0; vx < cell_corners; ++vx) {
336 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
338 vtkOrder = Dune::VTK::renumber(cit.type(), vx);
339 connectivity_inout[i + vtkOrder] = vxIdx;
347template <
class Gr
idView,
unsigned int partitions>
348template <
typename VectType>
353 const std::size_t check_size = connectivity_inout.size();
355 if (check_size <
static_cast<std::size_t
>(ncorners_)) {
357 OPM_THROW(std::runtime_error,
358 "Opm::GridDataOutput::writeConnectivity( VectType& "
359 "connectivity_inout ) "
361 +
") is not sufficient to fit the ncorners_ values (" +
std::to_string(ncorners_) +
")");
365 if (whichOrder ==
DUNE) {
367 for (
const auto& cit : elements(gridView_, dunePartition_)) {
368 auto cell_corners = cit.geometry().corners();
369 for (
auto vx = 0; vx < cell_corners; ++vx) {
370 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
371 connectivity_inout.data()[i + vx] = vxIdx;
377 for (
const auto& cit : elements(gridView_, dunePartition_)) {
378 auto cell_corners = cit.geometry().corners();
379 for (
auto vx = 0; vx < cell_corners; ++vx) {
380 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
382 vtkOrder = Dune::VTK::renumber(cit.type(), vx);
383 connectivity_inout.data()[i + vtkOrder] = vxIdx;
391template <
class Gr
idView,
unsigned int partitions>
392template <
typename Integer>
396 if (max_size < ncells_) {
398 OPM_THROW(std::runtime_error,
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 ("
404 offsets_inout[0] = 0;
405 for (
const auto& cit : elements(gridView_, dunePartition_)) {
406 auto cell_corners = cit.geometry().corners();
407 offsets_inout[i] = offsets_inout[i - 1] + cell_corners;
413template <
class Gr
idView,
unsigned int partitions>
414template <
typename VectType>
418 const std::size_t check_size = offsets_inout.size();
419 if (check_size <
static_cast<std::size_t
>(ncells_)) {
421 OPM_THROW(std::runtime_error,
422 "Opm::GridDataOutput::writeOffsetsCells( VectType& "
425 +
") is not sufficient to fit the ncells_ values (" +
std::to_string(ncells_) +
")");
431 offsets_inout.data()[0] = 0;
432 for (
const auto& cit : elements(gridView_, dunePartition_)) {
433 auto cell_corners = cit.geometry().corners();
434 offsets_inout.data()[i] = offsets_inout.data()[i - 1] + cell_corners;
440template <
class Gr
idView,
unsigned int partitions>
441template <
typename Integer>
445 if (max_size < ncells_) {
447 OPM_THROW(std::runtime_error,
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 ("
453 for (
const auto& cit : elements(gridView_, dunePartition_)) {
454 Integer vtktype =
static_cast<Integer
>(Dune::VTK::geometryType(cit.type()));
455 types_inout[i++] = vtktype;
460template <
class Gr
idView,
unsigned int partitions>
461template <
typename VectType>
465 const std::size_t check_size = types_inout.size();
467 if (check_size <
static_cast<std::size_t
>(ncells_)) {
468 OPM_THROW(std::runtime_error,
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 ("
475 for (
const auto& cit : elements(gridView_, dunePartition_)) {
476 int vtktype =
static_cast<int>(Dune::VTK::geometryType(cit.type()));
477 types_inout.data()[i++] = vtktype;
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>
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.
long writeGridPoints_AOS(T *xyz_inout, long max_size=0) const
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()
Definition: GridDataOutput_impl.hpp:55
long writeGridPoints(T *x_inout, T *y_inout, T *z_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:74
long writeOffsetsCells(Integer *offsets_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:394
long writeConnectivity(Integer *connectivity_inout, ConnectivityVertexOrder whichOrder, long max_size=0) const
Definition: GridDataOutput_impl.hpp:309
unsigned int partition_value_
Definition: GridDataOutput.hpp:359
void printGridDetails(std::ostream &outstr) const
Definition: GridDataOutput_impl.hpp:516
bool polyhedralCellPresent() const
Definition: GridDataOutput_impl.hpp:46
long writeCellTypes(Integer *types_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:443
int dimw_
Definition: GridDataOutput.hpp:374
long writeGridPoints_SOA(T *xyz_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:227
std::string getPartitionTypeString() const
Definition: GridDataOutput_impl.hpp:484
Definition: GridDataOutput.hpp:89
ConnectivityVertexOrder
Definition: GridDataOutput.hpp:93
@ DUNE
Definition: GridDataOutput.hpp:93
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)