GridDataOutput_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2023 Inria, Bretagne–Atlantique Research Center
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
20
21#include <dune/grid/common/rangegenerators.hh>
22#include <dune/grid/io/file/vtk/common.hh>
23
24#include <opm/common/ErrorMacros.hpp>
25
26#include <cassert>
27#include <cstddef>
28#include <iterator>
29#include <ostream>
30
31namespace Opm::GridDataOutput {
32
33template <class GridView, unsigned int partitions>
35 SimMeshDataAccessor(const GridView& gridView, Dune::PartitionSet<partitions> dunePartition)
36 : gridView_(gridView)
37 , dunePartition_(dunePartition)
38{
39 dimw_ = GridView::dimension; // this is an enum
40 partition_value_ = dunePartition.value;
42}
43
44template <class GridView, unsigned int partitions>
46{
47 for (const auto& cit : elements(gridView_, dunePartition_)) {
48 auto corner_geom = cit.geometry();
49 if (Dune::VTK::geometryType(corner_geom.type()) == Dune::VTK::polyhedron) {
50 return true;
51 }
52 }
53 return false;
54}
55
56template <class GridView, unsigned int partitions>
58{
59 // We include all the vertices for this ranks partition
60 const auto& vert_partition = vertices(gridView_, Dune::Partitions::all);
61 nvertices_ = std::distance(vert_partition.begin(), vert_partition.end());
62
63 const auto& cell_partition = elements(gridView_, dunePartition_);
64 ncells_ = 0;
65 ncorners_ = 0;
66 for (const auto& cit : cell_partition) {
67 auto corner_geom = cit.geometry();
68 ncorners_ += corner_geom.corners();
69 ++ncells_;
70 }
71}
72
73template <class GridView, unsigned int partitions>
74template <typename T>
76writeGridPoints(T* x_inout, T* y_inout, T* z_inout, long max_size) const
77{
78 if (max_size < nvertices_) {
79 OPM_THROW(std::runtime_error,
80 "Opm::GridDataOutput::writeGridPoints( T& x_inout, T& "
81 "y_inout, T& z_inout ) "
82 + " Input objects max_size (" + std::to_string(max_size)
83 + ") is not sufficient to fit the nvertices_ values (" + std::to_string(nvertices_) + ")");
84 }
85
86 long i = 0;
87 if (dimw_ == 3) {
88 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
89 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
90 x_inout[i] = static_cast<T>(xyz_local[0]);
91 y_inout[i] = static_cast<T>(xyz_local[1]);
92 z_inout[i] = static_cast<T>(xyz_local[2]);
93 i++;
94 }
95 } else if (dimw_ == 2) {
96 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
97 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
98 x_inout[i] = static_cast<T>(xyz_local[0]);
99 y_inout[i] = static_cast<T>(xyz_local[1]);
100 z_inout[i] = static_cast<T>(0.0);
101 i++;
102 }
103 }
104 assert(i == nvertices_); // As we are templated on the
105 // Dune::PartitionSet<partitions>, this cannot change
106 return i;
107}
108
109template <class GridView, unsigned int partitions>
110template <typename VectType>
112writeGridPoints(VectType& x_inout, VectType& y_inout, VectType& z_inout) const
113{
114 const std::size_t check_size_x = x_inout.size();
115 const std::size_t check_size_y = y_inout.size();
116 const std::size_t check_size_z = z_inout.size();
117
118 using VT = decltype(x_inout.data()[0]);
119
120 if ((check_size_x < static_cast<std::size_t>(nvertices_)) ||
121 (check_size_y < static_cast<std::size_t>(nvertices_)) ||
122 (check_size_z < static_cast<std::size_t>(nvertices_))) {
123 // assert(check_size >= nvertices_);
124 OPM_THROW(std::runtime_error,
125 "Opm::GridDataOutput::writeGridPoints( VectType& x_inout, VectType& "
126 "y_inout, VectType& z_inout ) At least one of the inputs"
127 + " object x size " + std::to_string(check_size_x) + " object y size "
128 + std::to_string(check_size_y) + " object z size " + std::to_string(check_size_z)
129 + " is not sufficient to fit the nvertices_ values( " + std::to_string(nvertices_) + " )");
130 }
131
132 long i = 0;
133 if (dimw_ == 3) {
134 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
135 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
136 x_inout.data()[i] = static_cast<VT>(xyz_local[0]);
137 y_inout.data()[i] = static_cast<VT>(xyz_local[1]);
138 z_inout.data()[i] = static_cast<VT>(xyz_local[2]);
139 i++;
140 }
141 } else if (dimw_ == 2) {
142 double td = 0.0;
143 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
144 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
145 x_inout.data()[i] = static_cast<VT>(xyz_local[0]);
146 y_inout.data()[i] = static_cast<VT>(xyz_local[1]);
147 z_inout.data()[i] = static_cast<VT>(td);
148 i++;
149 }
150 }
151 assert(i == nvertices_); // As we are templated on the
152 // Dune::PartitionSet<partitions>, this cannot change
153 return i;
154}
155
156template <class GridView, unsigned int partitions>
157template <typename T>
159writeGridPoints_AOS(T* xyz_inout, long max_size) const
160{
161 if (max_size < nvertices_ * 3) {
162 assert(max_size >= nvertices_ * 3);
163 OPM_THROW(std::runtime_error,
164 "Opm::GridDataOutput::writeGridPoints_AOS( T* xyz_inout ) " + " Input objects max_size ("
165 + std::to_string(max_size) + ") is not sufficient to fit the nvertices_ * 3 values ("
166 + std::to_string(nvertices_ * 3) + ")");
167 }
168
169 long i = 0;
170 if (dimw_ == 3) {
171 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
172 auto xyz_local = vit.geometry().corner(0);
173 xyz_inout[i++] = static_cast<T>(xyz_local[0]);
174 xyz_inout[i++] = static_cast<T>(xyz_local[1]);
175 xyz_inout[i++] = static_cast<T>(xyz_local[2]);
176 }
177 } else if (dimw_ == 2) {
178 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
179 auto xyz_local = vit.geometry().corner(0);
180 xyz_inout[i++] = static_cast<T>(xyz_local[0]);
181 xyz_inout[i++] = static_cast<T>(xyz_local[1]);
182 xyz_inout[i++] = static_cast<T>(0.0);
183 }
184 }
185 return i / 3;
186}
187
188template <class GridView, unsigned int partitions>
189template <typename VectType>
191writeGridPoints_AOS(VectType& xyz_inout) const
192{
193 const std::size_t check_size = xyz_inout.size();
194
195 using VT = decltype(xyz_inout.data()[0]);
196
197 if (check_size < static_cast<std::size_t>(nvertices_ * 3)) {
198 assert(check_size >= nvertices_ * 3);
199 OPM_THROW(std::runtime_error,
200 "Opm::GridDataOutput::writeGridPoints_AOS( VectType& xyz_inout ) "
201 + " Input objects check_size (" + std::to_string(check_size)
202 + ") is not sufficient to fit the nvertices_ * 3 values (" + std::to_string(nvertices_ * 3)
203 + ")");
204 }
205
206 long i = 0;
207 if (dimw_ == 3) {
208 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
209 auto xyz_local = vit.geometry().corner(0);
210 xyz_inout.data()[i++] = static_cast<VT>(xyz_local[0]);
211 xyz_inout[i++] = static_cast<VT>(xyz_local[1]);
212 xyz_inout[i++] = static_cast<VT>(xyz_local[2]);
213 }
214 } else if (dimw_ == 2) {
215 double td = 0.0;
216 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
217 auto xyz_local = vit.geometry().corner(0);
218 xyz_inout[i++] = static_cast<VT>(xyz_local[0]);
219 xyz_inout[i++] = static_cast<VT>(xyz_local[1]);
220 xyz_inout[i++] = static_cast<VT>(td);
221 }
222 }
223 return i / 3;
224}
225
226template <class GridView, unsigned int partitions>
227template <typename T>
229writeGridPoints_SOA(T* xyz_inout, long max_size) const
230{
231 if (max_size < nvertices_ * 3) {
232 // assert(max_size >= nvertices_ * 3);
233 OPM_THROW(std::runtime_error,
234 "Opm::GridDataOutput::writeGridPoints_SOA( T& xyz_inout ) " + " Input objects max_size ("
235 + std::to_string(max_size) + ") is not sufficient to fit the nvertices_ * 3 values ("
236 + std::to_string(nvertices_ * 3) + ")");
237 }
238
239 long i = 0;
240 // Get offsets into structure
241 T* xyz_inout_y = xyz_inout + nvertices_;
242 T* xyz_inout_z = xyz_inout + (2 * nvertices_);
243
244 if (dimw_ == 3) {
245 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
246 auto xyz_local = vit.geometry().corner(0);
247 xyz_inout[i] = static_cast<T>(xyz_local[0]);
248 xyz_inout_y[i] = static_cast<T>(xyz_local[1]);
249 xyz_inout_z[i] = static_cast<T>(xyz_local[2]);
250 i++;
251 }
252 } else if (dimw_ == 2) {
253 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
254 auto xyz_local = vit.geometry().corner(0);
255 xyz_inout[i] = static_cast<T>(xyz_local[0]);
256 xyz_inout_y[i] = static_cast<T>(xyz_local[1]);
257 xyz_inout_z[i] = static_cast<T>(0.0);
258 i++;
259 }
260 }
261 return i;
262}
263
264template <class GridView, unsigned int partitions>
265template <typename VectType>
267writeGridPoints_SOA(VectType& xyz_inout) const
268{
269 const std::size_t check_size = xyz_inout.size();
270
271 if (check_size < static_cast<std::size_t>(nvertices_ * 3)) {
272 // assert(check_size >= nvertices_ * 3);
273 OPM_THROW(std::runtime_error,
274 "Opm::GridDataOutput::writeGridPoints_SOA( VectType& xyz_inout ) "
275 + " Input objects check_size (" + std::to_string(check_size)
276 + ") is not sufficient to fit the nvertices_ * 3 values (" + std::to_string(nvertices_ * 3)
277 + ")");
278 }
279
280 using VT = decltype(xyz_inout.data()[0]);
281
282 long i = 0;
283 // Get offsets into structure
284 VT* xyz_inout_y = xyz_inout.data() + nvertices_;
285 VT* xyz_inout_z = xyz_inout.data() + (2 * nvertices_);
286
287 if (dimw_ == 3) {
288 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
289 auto xyz_local = vit.geometry().corner(0);
290 xyz_inout.data()[i] = static_cast<VT>(xyz_local[0]);
291 xyz_inout_y[i] = static_cast<VT>(xyz_local[1]);
292 xyz_inout_z[i] = static_cast<VT>(xyz_local[2]);
293 i++;
294 }
295 } else if (dimw_ == 2) {
296 double td = 0.0;
297 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
298 auto xyz_local = vit.geometry().corner(0);
299 xyz_inout.data()[i] = static_cast<VT>(xyz_local[0]);
300 xyz_inout_y[i] = static_cast<VT>(xyz_local[1]);
301 xyz_inout_z[i] = static_cast<VT>(td);
302 i++;
303 }
304 }
305 return i;
306}
307
308template <class GridView, unsigned int partitions>
309template <typename Integer>
311writeConnectivity(Integer* connectivity_inout,
312 ConnectivityVertexOrder whichOrder, long max_size) const
313{
314 if (max_size < ncorners_) {
315
316 OPM_THROW(std::runtime_error,
317 "Opm::GridDataOutput::writeConnectivity( T* connectivity_inout,... ) "
318 + " Input max_size value (" + std::to_string(max_size)
319 + ") is not sufficient to fit the ncorners_ values (" + std::to_string(ncorners_) + ")");
320 }
321
322 long i = 0;
323 if (whichOrder == DUNE) {
324 // DUNE order
325 for (const auto& cit : elements(gridView_, dunePartition_)) {
326 auto cell_corners = cit.geometry().corners();
327 for (auto vx = 0; vx < cell_corners; ++vx) {
328 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
329 connectivity_inout[i + vx] = vxIdx;
330 }
331 i += cell_corners;
332 }
333 } else {
334 // VTK order
335 for (const auto& cit : elements(gridView_, dunePartition_)) {
336 auto cell_corners = cit.geometry().corners();
337 for (auto vx = 0; vx < cell_corners; ++vx) {
338 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
339 int vtkOrder;
340 vtkOrder = Dune::VTK::renumber(cit.type(), vx);
341 connectivity_inout[i + vtkOrder] = vxIdx;
342 }
343 i += cell_corners;
344 }
345 }
346 return i;
347}
348
349template <class GridView, unsigned int partitions>
350template <typename VectType>
352writeConnectivity(VectType& connectivity_inout,
353 ConnectivityVertexOrder whichOrder) const
354{
355 const std::size_t check_size = connectivity_inout.size();
356
357 if (check_size < static_cast<std::size_t>(ncorners_)) {
358 // assert(check_size >= ncorners_);
359 OPM_THROW(std::runtime_error,
360 "Opm::GridDataOutput::writeConnectivity( VectType& "
361 "connectivity_inout ) "
362 + " Input objects size (" + std::to_string(check_size)
363 + ") is not sufficient to fit the ncorners_ values (" + std::to_string(ncorners_) + ")");
364 }
365
366 long i = 0;
367 if (whichOrder == DUNE) {
368 // DUNE order
369 for (const auto& cit : elements(gridView_, dunePartition_)) {
370 auto cell_corners = cit.geometry().corners();
371 for (auto vx = 0; vx < cell_corners; ++vx) {
372 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
373 connectivity_inout.data()[i + vx] = vxIdx;
374 }
375 i += cell_corners;
376 }
377 } else {
378 // VTK order
379 for (const auto& cit : elements(gridView_, dunePartition_)) {
380 auto cell_corners = cit.geometry().corners();
381 for (auto vx = 0; vx < cell_corners; ++vx) {
382 const int vxIdx = gridView_.indexSet().subIndex(cit, vx, 3);
383 int vtkOrder;
384 vtkOrder = Dune::VTK::renumber(cit.type(), vx);
385 connectivity_inout.data()[i + vtkOrder] = vxIdx;
386 }
387 i += cell_corners;
388 }
389 }
390 return i;
391}
392
393template <class GridView, unsigned int partitions>
394template <typename Integer>
396writeOffsetsCells(Integer* offsets_inout, long max_size) const
397{
398 if (max_size < ncells_) {
399 // assert(max_size >= ncells_);
400 OPM_THROW(std::runtime_error,
401 "Opm::GridDataOutput::writeOffsetsCells( T* offsets_inout ) " + " Input objects max_size ("
402 + std::to_string(max_size) + ") is not sufficient to fit the ncells_ values ("
403 + std::to_string(ncells_) + ")");
404 }
405 long i = 1;
406 offsets_inout[0] = 0;
407 for (const auto& cit : elements(gridView_, dunePartition_)) {
408 auto cell_corners = cit.geometry().corners();
409 offsets_inout[i] = offsets_inout[i - 1] + cell_corners;
410 i++;
411 }
412 return i; // This should be 1 greater than ncells_
413}
414
415template <class GridView, unsigned int partitions>
416template <typename VectType>
418writeOffsetsCells(VectType& offsets_inout) const
419{
420 const std::size_t check_size = offsets_inout.size();
421 if (check_size < static_cast<std::size_t>(ncells_)) {
422 // assert(check_size >= ncells_);
423 OPM_THROW(std::runtime_error,
424 "Opm::GridDataOutput::writeOffsetsCells( VectType& "
425 "offsets_inout ) "
426 + " Input objects size (" + std::to_string(offsets_inout.size())
427 + ") is not sufficient to fit the ncells_ values (" + std::to_string(ncells_) + ")");
428 }
429
430 // using VT = decltype(offsets_inout.data()[0]);
431
432 long i = 1;
433 offsets_inout.data()[0] = 0;
434 for (const auto& cit : elements(gridView_, dunePartition_)) {
435 auto cell_corners = cit.geometry().corners();
436 offsets_inout.data()[i] = offsets_inout.data()[i - 1] + cell_corners;
437 i++;
438 }
439 return i; // This should be 1 greater than ncells_
440}
441
442template <class GridView, unsigned int partitions>
443template <typename Integer>
445writeCellTypes(Integer* types_inout, long max_size) const
446{
447 if (max_size < ncells_) {
448 // assert(max_size >= ncells_);
449 OPM_THROW(std::runtime_error,
450 "Opm::GridDataOutput::writeCellTypes( T* types_inout ) " + " Input objects max_size ("
451 + std::to_string(max_size) + ") is not sufficient to fit the ncells_ values ("
452 + std::to_string(ncells_) + ")");
453 }
454 int i = 0;
455 for (const auto& cit : elements(gridView_, dunePartition_)) {
456 Integer vtktype = static_cast<Integer>(Dune::VTK::geometryType(cit.type()));
457 types_inout[i++] = vtktype;
458 }
459 return i;
460}
461
462template <class GridView, unsigned int partitions>
463template <typename VectType>
465writeCellTypes(VectType& types_inout) const
466{
467 const std::size_t check_size = types_inout.size();
468
469 if (check_size < static_cast<std::size_t>(ncells_)) {
470 OPM_THROW(std::runtime_error,
471 "Opm::GridDataOutput::writeCellTypes( VectType& types_inout ) " + " Input objects check_size ("
472 + std::to_string(check_size) + ") is not sufficient to fit the ncells_ values ("
473 + std::to_string(ncells_) + ")");
474 }
475
476 int i = 0;
477 for (const auto& cit : elements(gridView_, dunePartition_)) {
478 int vtktype = static_cast<int>(Dune::VTK::geometryType(cit.type()));
479 types_inout.data()[i++] = vtktype;
480 }
481 return i;
482}
483
484template <class GridView, unsigned int partitions>
487{
488 if (this->dunePartition_ == Dune::Partitions::all) {
489 return "Dune::Partitions::all";
490 }
491 if (this->dunePartition_ == Dune::Partitions::interior) {
492 return "Dune::Partitions::interior";
493 }
494 if (this->dunePartition_ == Dune::Partitions::interiorBorder) {
495 return "Dune::Partitions::interiorBorder";
496 }
497 if (this->dunePartition_ == Dune::Partitions::interiorBorderOverlap) {
498 return "Dune::Partitions::interiorBorderOverlap";
499 }
500 if (this->dunePartition_ == Dune::Partitions::front) {
501 return "Dune::Partitions::front";
502 }
503 if (this->dunePartition_ == Dune::Partitions::interiorBorderOverlapFront) {
504 return "Dune::Partitions::InteriorBorderOverlapFront";
505 }
506 if (this->dunePartition_ == Dune::Partitions::border) {
507 return "Dune::Partitions::border";
508 }
509 if (this->dunePartition_ == Dune::Partitions::ghost) {
510 return "Dune::Partitions::ghost";
511 }
512
513 return "Unknown Dune::PartitionSet<>";
514}
515
516template <class GridView, unsigned int partitions>
518printGridDetails(std::ostream& outstr) const
519{
520 outstr << "Dune Partition = " << partition_value_ << ", " << getPartitionTypeString() << std::endl;
521 outstr << "ncells_: " << getNCells() << std::endl;
522 outstr << "nvertices_: " << getNVertices() << std::endl;
523 outstr << "ncorners_: " << getNCorners() << std::endl;
524}
525
526} // namespace Opm::GridDataOutput
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:159
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:35
void countEntities()
Definition: GridDataOutput_impl.hpp:57
long writeGridPoints(T *x_inout, T *y_inout, T *z_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:76
long writeOffsetsCells(Integer *offsets_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:396
long writeConnectivity(Integer *connectivity_inout, ConnectivityVertexOrder whichOrder, long max_size=0) const
Definition: GridDataOutput_impl.hpp:311
unsigned int partition_value_
Definition: GridDataOutput.hpp:359
void printGridDetails(std::ostream &outstr) const
Definition: GridDataOutput_impl.hpp:518
bool polyhedralCellPresent() const
Definition: GridDataOutput_impl.hpp:45
long writeCellTypes(Integer *types_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:445
int dimw_
Definition: GridDataOutput.hpp:374
long writeGridPoints_SOA(T *xyz_inout, long max_size=0) const
Definition: GridDataOutput_impl.hpp:229
std::string getPartitionTypeString() const
Definition: GridDataOutput_impl.hpp:486
Definition: GridDataOutput.hpp:89
ConnectivityVertexOrder
Definition: GridDataOutput.hpp:93
@ DUNE
Definition: GridDataOutput.hpp:93
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)