version 3.10-dev
gridwriter.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_GRID_WRITER_HH
13#define DUMUX_IO_GRID_WRITER_HH
14
15#include <config.h>
16
17#if DUMUX_HAVE_GRIDFORMAT
18
19#include <ranges>
20#include <execution>
21#include <type_traits>
22
23#include <gridformat/gridformat.hpp>
24#include <gridformat/traits/dune.hpp>
25
26#include <dune/common/exceptions.hh>
28
29namespace Dumux::IO {
30
31
32#ifndef DOXYGEN
33namespace Detail {
34
35template<typename Grid, typename Format, typename Comm, typename... Args>
36auto makeParallelWriter(const Grid& grid, const Format& fmt, const Dune::Communication<Comm>& comm, Args&&... args)
37{
38 return comm.size() > 1
39 ? GridFormat::Writer<Grid>{fmt, grid, static_cast<Comm>(comm), std::forward<Args>(args)...}
40 : GridFormat::Writer<Grid>{fmt, grid, std::forward<Args>(args)...};
41}
42
43template<typename Grid, typename Format, typename... Args>
44auto makeParallelWriter(const Grid& grid, const Format& fmt, const Dune::Communication<Dune::No_Comm>&, Args&&... args)
45{ return GridFormat::Writer<Grid>{fmt, grid, std::forward<Args>(args)...}; }
46
47template<typename Grid, typename Format, typename... Args>
48auto makeWriter(const Grid& grid, const Format& fmt, Args&&... args)
49{
50 const auto& comm = GridFormat::Dune::Traits::GridView<Grid>::get(grid).comm();
51 return makeParallelWriter(grid, fmt, comm, std::forward<Args>(args)...);
52}
53
54template<typename T>
55concept Container = requires(const T& t) {
56 { t.size() };
57 { t[std::size_t{}] };
58};
59
60} // namespace Detail
61#endif // DOXYGEN
62
63
64namespace VTK { using namespace GridFormat::VTK; }
65namespace Format { using namespace GridFormat::Formats; }
66namespace Encoding { using namespace GridFormat::Encoding; }
67namespace Compression { using namespace GridFormat::Compression; using GridFormat::none; }
68namespace Precision {
69 using GridFormat::float32;
70 using GridFormat::float64;
71
72 using GridFormat::uint64;
73 using GridFormat::uint32;
74 using GridFormat::uint16;
75 using GridFormat::uint8;
76
77 using GridFormat::int64;
78 using GridFormat::int32;
79 using GridFormat::int16;
80 using GridFormat::int8;
81} // namespace Precision
82
83
88template<int order>
89struct Order { static_assert(order > 0, "order must be > 0"); };
90
91template<int o>
92inline constexpr auto order = Order<o>{};
93
103template<GridFormat::Concepts::Grid GridView, int order = 1>
104class GridWriter
105{
106 using Grid = std::conditional_t<
107 (order > 1),
108 GridFormat::Dune::LagrangePolynomialGrid<GridView>,
109 GridView
110 >;
111 using Cell = GridFormat::Cell<Grid>;
112 using Vertex = typename GridView::template Codim<GridView::dimension>::Entity;
113 using Element = typename GridView::template Codim<0>::Entity;
114 using Coordinate = typename Element::Geometry::GlobalCoordinate;
115 using Writer = GridFormat::Writer<Grid>;
116
117 public:
122 template<typename Format>
123 explicit GridWriter(const Format& fmt,
124 const GridView& gridView,
125 const Order<order>& = {})
126 : gridView_{gridView}
127 , grid_{makeGrid_(gridView)}
128 , writer_{Detail::makeWriter(grid_, fmt)}
129 { writer_.set_meta_data("rank", gridView_.comm().rank()); }
130
135 template<typename Format>
136 explicit GridWriter(const Format& fmt,
137 const GridView& gridView,
138 const std::string& filename,
139 const Order<order>& = {})
140 : gridView_{gridView}
141 , grid_{makeGrid_(gridView)}
142 , writer_{Detail::makeWriter(grid_, fmt, filename)}
143 { writer_.set_meta_data("rank", gridView_.comm().rank()); }
144
149 std::string write(const std::string& name) const
150 { return writer_.write(name); }
151
156 template<std::floating_point T>
157 std::string write(T time) const
158 { return writer_.write(time); }
159
161 template<GridFormat::Concepts::CellFunction<GridView> F,
162 GridFormat::Concepts::Scalar T = GridFormat::FieldScalar<std::invoke_result_t<F, Element>>>
163 void setCellField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec = {})
164 { writer_.set_cell_field(name, std::move(f), prec); }
165
167 template<GridFormat::Dune::Concepts::Function<GridView> F>
168 void setCellField(const std::string& name, F&& f)
169 { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name); }
170
172 template<typename F, GridFormat::Concepts::Scalar T>
173 void setCellField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec)
174 { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name, prec); }
175
177 template<GridFormat::Concepts::PointFunction<GridView> F,
178 GridFormat::Concepts::Scalar T = GridFormat::FieldScalar<std::invoke_result_t<F, Vertex>>>
179 void setPointField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec = {})
180 {
181 static_assert(order == 1, "Point lambdas can only be used for order == 1. Use Dune::Functions instead.");
182 writer_.set_point_field(name, std::move(f), prec);
183 }
184
186 template<GridFormat::Dune::Concepts::Function<GridView> F>
187 void setPointField(const std::string& name, F&& f)
188 { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name); }
189
191 template<GridFormat::Dune::Concepts::Function<GridView> F, GridFormat::Concepts::Scalar T>
192 void setPointField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec)
193 { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name, prec); }
194
196 void clear()
197 { writer_.clear(); }
198
200 void update()
201 {
202 if constexpr (order > 1)
203 grid_.update(gridView_);
204 }
205
206 private:
207 Grid makeGrid_(const GridView& gv) const
208 {
209 if constexpr (order > 1)
210 return Grid{gv, order};
211 else
212 return gv;
213 }
214
215 GridView gridView_;
216 Grid grid_;
217 Writer writer_;
218};
219
224template<typename GridVariables, typename SolutionVector>
225class OutputModule : private GridWriter<typename GridVariables::GridGeometry::GridView, 1> {
226 using ParentType = GridWriter<typename GridVariables::GridGeometry::GridView, 1>;
227 using GridView = typename GridVariables::GridGeometry::GridView;
228
229 static constexpr bool isCVFE = DiscretizationMethods::isCVFE<typename GridVariables::GridGeometry::DiscretizationMethod>;
230 static constexpr int dimWorld = GridVariables::GridGeometry::GridView::dimensionworld;
231 using Scalar = typename GridVariables::Scalar;
232 using Vector = Dune::FieldVector<Scalar, dimWorld>;
233 using Tensor = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
234 using VolVar = typename GridVariables::VolumeVariables;
235
236 class VolVarFieldStorage;
237
238 public:
239 using VolumeVariables = VolVar;
240
241 static constexpr auto defaultFileFormat = IO::Format::pvd_with(
242 IO::Format::vtu.with({
243 .encoder = IO::Encoding::ascii,
244 .compressor = IO::Compression::none,
245 .data_format = VTK::DataFormat::inlined
246 })
247 );
248
252 explicit OutputModule(const GridVariables& gridVariables,
253 const SolutionVector& sol,
254 const std::string& filename)
255 : ParentType{defaultFileFormat, gridVariables.gridGeometry().gridView(), filename, order<1>}
256 , gridVariables_{gridVariables}
257 , solutionVector_{sol}
258 {}
259
264 template<typename Format>
265 explicit OutputModule(const Format& fmt,
266 const GridVariables& gridVariables,
267 const SolutionVector& sol)
268 : ParentType{fmt, gridVariables.gridGeometry().gridView(), order<1>}
269 , gridVariables_{gridVariables}
270 , solutionVector_{sol}
271 {}
272
277 template<typename Format>
278 explicit OutputModule(const Format& fmt,
279 const GridVariables& gridVariables,
280 const SolutionVector& sol,
281 const std::string& filename)
282 : ParentType{fmt, gridVariables.gridGeometry().gridView(), filename, order<1>}
283 , gridVariables_{gridVariables}
284 , solutionVector_{sol}
285 {}
286
290 template<std::invocable<const VolumeVariables&> VolVarFunction>
291 void addVolumeVariable(VolVarFunction&& f, const std::string& name)
292 {
293 using ResultType = std::invoke_result_t<std::remove_cvref_t<VolVarFunction>, const VolumeVariables&>;
294 if constexpr (GridFormat::Concepts::Scalar<ResultType>)
295 setVolVarField_<ResultType>(name, volVarFields_.registerScalarField(name, [_f=std::move(f)] (const auto& vv) {
296 return static_cast<Scalar>(_f(vv));
297 }));
298 else if constexpr (GridFormat::mdrange_dimension<ResultType> == 1)
299 setVolVarField_<GridFormat::MDRangeScalar<ResultType>>(name, volVarFields_.registerVectorField(name, [_f=std::move(f)] (const auto& vv) {
300 return VolVarFieldStorage::toStorageVector(_f(vv));
301 }));
302 else if constexpr (GridFormat::mdrange_dimension<ResultType> == 2)
303 setVolVarField_<GridFormat::MDRangeScalar<ResultType>>(name, volVarFields_.registerTensorField(name, [_f=std::move(f)] (const auto& vv) {
304 return VolVarFieldStorage::toStorageTensor(_f(vv));
305 }));
306 else
307 {
308 static_assert(
309 Dune::AlwaysFalse<VolVarFunction>::value,
310 "Could not identify the given volume variable as scalar, vector or tensor."
311 );
312 }
313 }
314
319 template<Detail::Container C>
320 void addField(const C& values, const std::string& name)
321 { addField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
322
327 template<Detail::Container C, GridFormat::Concepts::Scalar T>
328 void addField(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
329 {
330 const bool hasCellSize = values.size() == gridVariables_.gridGeometry().elementMapper().size();
331 const bool hasPointSize = values.size() == gridVariables_.gridGeometry().vertexMapper().size();
332 if (hasCellSize && hasPointSize)
333 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of field type failed. Please use addCellField or addPointField instead.");
334 if (!hasCellSize && !hasPointSize)
335 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of field type failed. Given container size does not match neither the number of points nor cells.");
336
337 if (hasCellSize)
338 addCellField(values, name, prec);
339 else
340 addPointField(values, name, prec);
341 }
342
346 template<GridFormat::Concepts::PointFunction<GridView> DofFunction>
347 void addPointField(DofFunction&& f, const std::string& name)
348 { this->setPointField(name, std::forward<DofFunction>(f)); }
349
353 template<Detail::Container C>
354 void addPointField(const C& values, const std::string& name)
355 { addPointField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
356
360 template<Detail::Container C, GridFormat::Concepts::Scalar T>
361 void addPointField(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
362 {
363 if (values.size() != gridVariables_.gridGeometry().vertexMapper().size())
364 DUNE_THROW(Dune::InvalidStateException, "Given container does not match the number of points in the grid");
365
366 this->setPointField(name, [&] (const auto& vertex) {
367 return values[gridVariables_.gridGeometry().vertexMapper().index(vertex)];
368 }, prec);
369 }
370
374 template<GridFormat::Concepts::CellFunction<GridView> DofFunction>
375 void addCellField(DofFunction&& f, const std::string& name)
376 { this->setCellField(name, std::forward<DofFunction>(f)); }
377
381 template<Detail::Container C>
382 void addCellField(const C& values, const std::string& name)
383 { addCellField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
384
388 template<Detail::Container C, GridFormat::Concepts::Scalar T>
389 void addCellField(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
390 {
391 if (values.size() != gridVariables_.gridGeometry().elementMapper().size())
392 DUNE_THROW(Dune::InvalidStateException, "Given container does not match the number of cells in the grid");
393
394 this->setCellField(name, [&] (const auto& element) {
395 return values[gridVariables_.gridGeometry().elementMapper().index(element)];
396 }, prec);
397 }
398
402 std::string write(const std::string& name)
403 {
404 volVarFields_.updateFieldData(gridVariables_, solutionVector_);
405 auto filename = ParentType::write(name);
406 volVarFields_.clearFieldData();
407 return filename;
408 }
409
413 template<std::floating_point T>
414 std::string write(T time)
415 {
416 volVarFields_.updateFieldData(gridVariables_, solutionVector_);
417 auto filename = ParentType::write(time);
418 volVarFields_.clearFieldData();
419 return filename;
420 }
421
423 void clear() {
424 ParentType::clear();
425 volVarFields_.clear();
426 }
427
428 private:
429 template<typename ResultType, typename Id>
430 void setVolVarField_(const std::string& name, Id&& volVarFieldId)
431 {
432 auto dofEntityField = [&, _id=std::move(volVarFieldId)] (const auto& entity) {
433 return volVarFields_.getValue(_id, gridVariables_.gridGeometry().dofMapper().index(entity));
434 };
435 if constexpr (isCVFE)
436 this->setPointField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
437 else
438 this->setCellField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
439 }
440
441 const GridVariables& gridVariables_;
442 const SolutionVector& solutionVector_;
443 VolVarFieldStorage volVarFields_;
444};
445
446// Class to store vol var fields; implementation detail of the OutputModule class
447template<typename GridVariables, typename SolutionVector>
448class OutputModule<GridVariables, SolutionVector>::VolVarFieldStorage
449{
450 enum class FieldType
451 { scalar, vector, tensor };
452
453 struct FieldInfo
454 {
455 std::string name;
456 FieldType type;
457 std::size_t index;
458 };
459
460 template<typename T>
461 struct FieldStorage
462 {
463 std::vector<T> data;
464 std::function<T(const VolVar&)> getter;
465 };
466
467public:
468 template<FieldType ft>
469 struct FieldId { std::size_t index; };
470
471 template<std::ranges::range R>
472 static constexpr auto toStorageVector(R&& in)
473 {
474 Vector result;
475 std::ranges::copy(in, result.begin());
476 return result;
477 }
478
479 template<GridFormat::Concepts::MDRange<2> R>
480 static constexpr auto toStorageTensor(R&& in)
481 {
482 Tensor result;
483 std::ranges::for_each(in, [&, i=0] (const auto& row) mutable {
484 std::ranges::copy(row, result[i++].begin());
485 });
486 return result;
487 }
488
489 template<FieldType ft>
490 const auto& getValue(const FieldId<ft>& id, std::size_t idx) const
491 {
492 if constexpr (ft == FieldType::scalar)
493 return scalarFieldStorage_.at(id.index).data.at(idx);
494 else if constexpr (ft == FieldType::vector)
495 return vectorFieldStorage_.at(id.index).data.at(idx);
496 else
497 return tensorFieldStorage_.at(id.index).data.at(idx);
498 }
499
500 auto registerScalarField(std::string name, std::function<Scalar(const VolVar&)> f)
501 { return register_<FieldType::scalar>(std::move(name), scalarFieldStorage_, std::move(f)); }
502
503 auto registerVectorField(std::string name, std::function<Vector(const VolVar&)> f)
504 { return register_<FieldType::vector>(std::move(name), vectorFieldStorage_, std::move(f)); }
505
506 auto registerTensorField(std::string name, std::function<Tensor(const VolVar&)> f)
507 { return register_<FieldType::tensor>(std::move(name), tensorFieldStorage_, std::move(f)); }
508
509 void updateFieldData(const GridVariables& gridVars, const SolutionVector& x)
510 {
511 resizeFieldData_(gridVars.gridGeometry().numDofs());
512 const auto range = GridFormat::cells(gridVars.gridGeometry().gridView());
513 std::for_each(
514#if __cpp_lib_parallel_algorithm >= 201603L
515 std::execution::par_unseq,
516#endif
517 std::ranges::begin(range),
518 std::ranges::end(range),
519 [&] (const auto& element) {
520 auto fvGeometry = localView(gridVars.gridGeometry()).bindElement(element);
521 auto elemVolVars = localView(gridVars.curGridVolVars()).bindElement(element, fvGeometry, x);
522 for (const auto& scv : scvs(fvGeometry))
523 {
524 const auto& volVars = elemVolVars[scv];
525 for (auto& s : scalarFieldStorage_) { s.data.at(scv.dofIndex()) = s.getter(volVars); }
526 for (auto& s : vectorFieldStorage_) { s.data.at(scv.dofIndex()) = s.getter(volVars); }
527 for (auto& s : tensorFieldStorage_) { s.data.at(scv.dofIndex()) = s.getter(volVars); }
528 }
529 });
530 }
531
532 void clearFieldData()
533 {
534 for (auto& s : scalarFieldStorage_) { s.data.clear(); }
535 for (auto& s : vectorFieldStorage_) { s.data.clear(); }
536 for (auto& s : tensorFieldStorage_) { s.data.clear(); }
537 }
538
539 void clear()
540 {
541 fields_.clear();
542 scalarFieldStorage_.clear();
543 vectorFieldStorage_.clear();
544 tensorFieldStorage_.clear();
545 }
546
547private:
548 template<FieldType ft, typename T>
549 auto register_(std::string&& name,
550 std::vector<FieldStorage<T>>& storage,
551 std::function<T(const VolVar&)>&& f)
552 {
553 if (exists_<ft>(name))
554 DUNE_THROW(Dune::InvalidStateException, "Volume variables field '" << name << "' is already defined.");
555
556 FieldId<ft> id{storage.size()};
557 fields_.emplace_back(FieldInfo{std::move(name), ft, id.index});
558 storage.push_back({{}, std::move(f)});
559 return id;
560 }
561
562 template<FieldType ft>
563 bool exists_(const std::string& name) const
564 {
565 return std::ranges::any_of(fields_, [&] (const FieldInfo& info) {
566 return info.type == ft && info.name == name;
567 });
568 }
569
570 void resizeFieldData_(std::size_t size)
571 {
572 std::ranges::for_each(scalarFieldStorage_, [&] (auto& s) { s.data.resize(size); });
573 std::ranges::for_each(vectorFieldStorage_, [&] (auto& s) { s.data.resize(size); });
574 std::ranges::for_each(tensorFieldStorage_, [&] (auto& s) { s.data.resize(size); });
575 }
576
577 std::vector<FieldInfo> fields_;
578 std::vector<FieldStorage<Scalar>> scalarFieldStorage_;
579 std::vector<FieldStorage<Vector>> vectorFieldStorage_;
580 std::vector<FieldStorage<Tensor>> tensorFieldStorage_;
581};
582
583} // namespace Dumux::IO
584
585#else // DUMUX_HAVE_GRIDFORMAT
586
587namespace Dumux::IO {
588
589template<class... Args>
591{
592public:
593 template<class... _Args>
594 GridWriter(_Args&&...)
595 {
596 static_assert(
597 false,
598 "GridWriter only available when the GridFormat library is available. "
599 "Use `git submodule update --init` to pull it and reconfigure the project "
600 "(note: C++20 is required)."
601 );
602 }
603};
604
605template<class... Args>
607{
608public:
609 template<class... _Args>
610 OutputModule(_Args&&...)
611 {
612 static_assert(
613 false,
614 "OutputModule only available when the GridFormat library is available. "
615 "Use `git submodule update --init` to pull it and reconfigure the project "
616 "(note: C++20 is required)."
617 );
618 }
619};
620
621} // namespace Dumux::IO
622
623#endif // DUMUX_HAVE_GRIDFORMAT
624#endif // DUMUX_IO_GRID_HH
Definition: gridwriter.hh:591
GridWriter(_Args &&...)
Definition: gridwriter.hh:594
Definition: gridwriter.hh:607
OutputModule(_Args &&...)
Definition: gridwriter.hh:610
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
FieldType
Identifier for vtk field types.
Definition: fieldtype.hh:22
The available discretization methods in Dumux.
constexpr None none
Definition: method.hh:153
constexpr bool isCVFE
Definition: method.hh:67
Definition: gridwriter.hh:587