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 <algorithm>
27#include <cassert>
28#include <cstddef>
29#include <iterator>
30#include <ostream>
31
32namespace Opm::GridDataOutput {
33
34template <class GridView, unsigned int partitions>
36 SimMeshDataAccessor(const GridView& gridView, Dune::PartitionSet<partitions> dunePartition)
37 : gridView_(gridView)
38 , dunePartition_(dunePartition)
39{
40 dimw_ = GridView::dimension; // this is an enum
41 partition_value_ = dunePartition.value;
43}
44
45template <class GridView, unsigned int partitions>
47{
48 const auto& elems = elements(gridView_, dunePartition_);
49 return std::any_of(elems.begin(), elems.end(),
50 [](const auto& cit)
51 { return Dune::VTK::geometryType(cit.geometry().type()) == Dune::VTK::polyhedron; });
52}
53
54template <class GridView, unsigned int partitions>
56{
57 // We include all the vertices for this ranks partition
58 const auto& vert_partition = vertices(gridView_, Dune::Partitions::all);
59 nvertices_ = std::distance(vert_partition.begin(), vert_partition.end());
60
61 const auto& cell_partition = elements(gridView_, dunePartition_);
62 ncells_ = 0;
63 ncorners_ = 0;
64 for (const auto& cit : cell_partition) {
65 auto corner_geom = cit.geometry();
66 ncorners_ += corner_geom.corners();
67 ++ncells_;
68 }
69}
70
71template <class GridView, unsigned int partitions>
72template <typename T>
74writeGridPoints(T* x_inout, T* y_inout, T* z_inout, long max_size) const
75{
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 ) "
80 + " Input objects max_size (" + std::to_string(max_size)
81 + ") is not sufficient to fit the nvertices_ values (" + std::to_string(nvertices_) + ")");
82 }
83
84 long i = 0;
85 if (dimw_ == 3) {
86 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
87 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
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]);
91 i++;
92 }
93 } else if (dimw_ == 2) {
94 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
95 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
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);
99 i++;
100 }
101 }
102 assert(i == nvertices_); // As we are templated on the
103 // Dune::PartitionSet<partitions>, this cannot change
104 return i;
105}
106
107template <class GridView, unsigned int partitions>
108template <typename VectType>
110writeGridPoints(VectType& x_inout, VectType& y_inout, VectType& z_inout) const
111{
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();
115
116 using VT = decltype(x_inout.data()[0]);
117
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_))) {
121 // assert(check_size >= 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 "
126 + std::to_string(check_size_y) + " object z size " + std::to_string(check_size_z)
127 + " is not sufficient to fit the nvertices_ values( " + std::to_string(nvertices_) + " )");
128 }
129
130 long i = 0;
131 if (dimw_ == 3) {
132 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
133 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
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]);
137 i++;
138 }
139 } else if (dimw_ == 2) {
140 double td = 0.0;
141 for (const auto& vit : vertices(gridView_, Dune::Partitions::all)) {
142 auto xyz_local = vit.geometry().corner(0); // vertices only have one corner
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);
146 i++;
147 }
148 }
149 assert(i == nvertices_); // As we are templated on the
150 // Dune::PartitionSet<partitions>, this cannot change
151 return i;
152}
153
154template <class GridView, unsigned int partitions>
155template <typename T>
157writeGridPoints_AOS(T* xyz_inout, long max_size) const
158{
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 ("
164 + std::to_string(nvertices_ * 3) + ")");
165 }
166
167 long i = 0;
168 if (dimw_ == 3) {
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]);
174 }
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);
181 }
182 }
183 return i / 3;
184}
185
186template <class GridView, unsigned int partitions>
187template <typename VectType>
189writeGridPoints_AOS(VectType& xyz_inout) const
190{
191 const std::size_t check_size = xyz_inout.size();
192
193 using VT = decltype(xyz_inout.data()[0]);
194
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 ) "
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)
201 + ")");
202 }
203
204 long i = 0;
205 if (dimw_ == 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]);
211 }
212 } else if (dimw_ == 2) {
213 double td = 0.0;
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);
219 }
220 }
221 return i / 3;
222}
223
224template <class GridView, unsigned int partitions>
225template <typename T>
227writeGridPoints_SOA(T* xyz_inout, long max_size) const
228{
229 if (max_size < nvertices_ * 3) {
230 // assert(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 ("
234 + std::to_string(nvertices_ * 3) + ")");
235 }
236
237 long i = 0;
238 // Get offsets into structure
239 T* xyz_inout_y = xyz_inout + nvertices_;
240 T* xyz_inout_z = xyz_inout + (2 * nvertices_);
241
242 if (dimw_ == 3) {
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]);
248 i++;
249 }
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);
256 i++;
257 }
258 }
259 return i;
260}
261
262template <class GridView, unsigned int partitions>
263template <typename VectType>
265writeGridPoints_SOA(VectType& xyz_inout) const
266{
267 const std::size_t check_size = xyz_inout.size();
268
269 if (check_size < static_cast<std::size_t>(nvertices_ * 3)) {
270 // assert(check_size >= nvertices_ * 3);
271 OPM_THROW(std::runtime_error,
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)
275 + ")");
276 }
277
278 using VT = decltype(xyz_inout.data()[0]);
279
280 long i = 0;
281 // Get offsets into structure
282 VT* xyz_inout_y = xyz_inout.data() + nvertices_;
283 VT* xyz_inout_z = xyz_inout.data() + (2 * nvertices_);
284
285 if (dimw_ == 3) {
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]);
291 i++;
292 }
293 } else if (dimw_ == 2) {
294 double td = 0.0;
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);
300 i++;
301 }
302 }
303 return i;
304}
305
306template <class GridView, unsigned int partitions>
307template <typename Integer>
309writeConnectivity(Integer* connectivity_inout,
310 ConnectivityVertexOrder whichOrder, long max_size) const
311{
312 if (max_size < ncorners_) {
313
314 OPM_THROW(std::runtime_error,
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_) + ")");
318 }
319
320 long i = 0;
321 if (whichOrder == DUNE) {
322 // DUNE order
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;
328 }
329 i += cell_corners;
330 }
331 } else {
332 // VTK order
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);
337 int vtkOrder;
338 vtkOrder = Dune::VTK::renumber(cit.type(), vx);
339 connectivity_inout[i + vtkOrder] = vxIdx;
340 }
341 i += cell_corners;
342 }
343 }
344 return i;
345}
346
347template <class GridView, unsigned int partitions>
348template <typename VectType>
350writeConnectivity(VectType& connectivity_inout,
351 ConnectivityVertexOrder whichOrder) const
352{
353 const std::size_t check_size = connectivity_inout.size();
354
355 if (check_size < static_cast<std::size_t>(ncorners_)) {
356 // assert(check_size >= ncorners_);
357 OPM_THROW(std::runtime_error,
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_) + ")");
362 }
363
364 long i = 0;
365 if (whichOrder == DUNE) {
366 // DUNE order
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;
372 }
373 i += cell_corners;
374 }
375 } else {
376 // VTK order
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);
381 int vtkOrder;
382 vtkOrder = Dune::VTK::renumber(cit.type(), vx);
383 connectivity_inout.data()[i + vtkOrder] = vxIdx;
384 }
385 i += cell_corners;
386 }
387 }
388 return i;
389}
390
391template <class GridView, unsigned int partitions>
392template <typename Integer>
394writeOffsetsCells(Integer* offsets_inout, long max_size) const
395{
396 if (max_size < ncells_) {
397 // assert(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 ("
401 + std::to_string(ncells_) + ")");
402 }
403 long i = 1;
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;
408 i++;
409 }
410 return i; // This should be 1 greater than ncells_
411}
412
413template <class GridView, unsigned int partitions>
414template <typename VectType>
416writeOffsetsCells(VectType& offsets_inout) const
417{
418 const std::size_t check_size = offsets_inout.size();
419 if (check_size < static_cast<std::size_t>(ncells_)) {
420 // assert(check_size >= ncells_);
421 OPM_THROW(std::runtime_error,
422 "Opm::GridDataOutput::writeOffsetsCells( VectType& "
423 "offsets_inout ) "
424 + " Input objects size (" + std::to_string(offsets_inout.size())
425 + ") is not sufficient to fit the ncells_ values (" + std::to_string(ncells_) + ")");
426 }
427
428 // using VT = decltype(offsets_inout.data()[0]);
429
430 long i = 1;
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;
435 i++;
436 }
437 return i; // This should be 1 greater than ncells_
438}
439
440template <class GridView, unsigned int partitions>
441template <typename Integer>
443writeCellTypes(Integer* types_inout, long max_size) const
444{
445 if (max_size < ncells_) {
446 // assert(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 ("
450 + std::to_string(ncells_) + ")");
451 }
452 int i = 0;
453 for (const auto& cit : elements(gridView_, dunePartition_)) {
454 Integer vtktype = static_cast<Integer>(Dune::VTK::geometryType(cit.type()));
455 types_inout[i++] = vtktype;
456 }
457 return i;
458}
459
460template <class GridView, unsigned int partitions>
461template <typename VectType>
463writeCellTypes(VectType& types_inout) const
464{
465 const std::size_t check_size = types_inout.size();
466
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 ("
471 + std::to_string(ncells_) + ")");
472 }
473
474 int i = 0;
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;
478 }
479 return i;
480}
481
482template <class GridView, unsigned int partitions>
485{
486 if (this->dunePartition_ == Dune::Partitions::all) {
487 return "Dune::Partitions::all";
488 }
489 if (this->dunePartition_ == Dune::Partitions::interior) {
490 return "Dune::Partitions::interior";
491 }
492 if (this->dunePartition_ == Dune::Partitions::interiorBorder) {
493 return "Dune::Partitions::interiorBorder";
494 }
495 if (this->dunePartition_ == Dune::Partitions::interiorBorderOverlap) {
496 return "Dune::Partitions::interiorBorderOverlap";
497 }
498 if (this->dunePartition_ == Dune::Partitions::front) {
499 return "Dune::Partitions::front";
500 }
501 if (this->dunePartition_ == Dune::Partitions::interiorBorderOverlapFront) {
502 return "Dune::Partitions::InteriorBorderOverlapFront";
503 }
504 if (this->dunePartition_ == Dune::Partitions::border) {
505 return "Dune::Partitions::border";
506 }
507 if (this->dunePartition_ == Dune::Partitions::ghost) {
508 return "Dune::Partitions::ghost";
509 }
510
511 return "Unknown Dune::PartitionSet<>";
512}
513
514template <class GridView, unsigned int partitions>
516printGridDetails(std::ostream& outstr) const
517{
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;
522}
523
524} // 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: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)