12#ifndef DUMUX_IO_GRID_WRITER_HH
13#define DUMUX_IO_GRID_WRITER_HH
17#if DUMUX_HAVE_GRIDFORMAT
23#include <gridformat/gridformat.hpp>
24#include <gridformat/traits/dune.hpp>
26#include <dune/common/exceptions.hh>
35template<
typename Grid,
typename Format,
typename Comm,
typename... Args>
36auto makeParallelWriter(
const Grid& grid,
const Format& fmt,
const Dune::Communication<Comm>& comm, Args&&... args)
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)...};
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)...}; }
47template<
typename Grid,
typename Format,
typename... Args>
48auto makeWriter(
const Grid& grid,
const Format& fmt, Args&&... args)
50 const auto& comm = GridFormat::Dune::Traits::GridView<Grid>::get(grid).comm();
51 return makeParallelWriter(grid, fmt, comm, std::forward<Args>(args)...);
55concept Container =
requires(
const T& t) {
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; }
69 using GridFormat::float32;
70 using GridFormat::float64;
72 using GridFormat::uint64;
73 using GridFormat::uint32;
74 using GridFormat::uint16;
75 using GridFormat::uint8;
77 using GridFormat::int64;
78 using GridFormat::int32;
79 using GridFormat::int16;
80 using GridFormat::int8;
89struct Order {
static_assert(order > 0,
"order must be > 0"); };
92inline constexpr auto order = Order<o>{};
103template<Gr
idFormat::Concepts::Gr
id Gr
idView,
int order = 1>
106 using Grid = std::conditional_t<
108 GridFormat::Dune::LagrangePolynomialGrid<GridView>,
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>;
122 template<
typename Format>
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()); }
135 template<
typename Format>
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()); }
149 std::string write(
const std::string& name)
const
150 {
return writer_.write(name); }
156 template<std::
floating_po
int T>
157 std::string write(T time)
const
158 {
return writer_.write(time); }
161 template<Gr
idFormat::Concepts::CellFunction<Gr
idView> F,
162 Gr
idFormat::Concepts::Scalar T = Gr
idFormat::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); }
167 template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F>
168 void setCellField(
const std::string& name, F&& f)
169 { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name); }
172 template<
typename F, Gr
idFormat::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); }
177 template<Gr
idFormat::Concepts::Po
intFunction<Gr
idView> F,
178 Gr
idFormat::Concepts::Scalar T = Gr
idFormat::FieldScalar<std::invoke_result_t<F, Vertex>>>
179 void setPointField(
const std::string& name, F&& f,
const GridFormat::Precision<T>& prec = {})
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);
186 template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F>
187 void setPointField(
const std::string& name, F&& f)
188 { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name); }
191 template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F, Gr
idFormat::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); }
202 if constexpr (order > 1)
203 grid_.update(gridView_);
207 Grid makeGrid_(
const GridView& gv)
const
209 if constexpr (order > 1)
210 return Grid{gv, order};
224template<
typename Gr
idVariables,
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;
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;
236 class VolVarFieldStorage;
239 using VolumeVariables = VolVar;
241 static constexpr auto defaultFileFormat = IO::Format::pvd_with(
242 IO::Format::vtu.with({
243 .encoder = IO::Encoding::ascii,
245 .data_format = VTK::DataFormat::inlined
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}
264 template<
typename Format>
266 const GridVariables& gridVariables,
267 const SolutionVector& sol)
268 : ParentType{fmt, gridVariables.gridGeometry().gridView(), order<1>}
269 , gridVariables_{gridVariables}
270 , solutionVector_{sol}
277 template<
typename Format>
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}
290 template<std::invocable<const VolumeVariables&> VolVarFunction>
291 void addVolumeVariable(VolVarFunction&& f,
const std::string& name)
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));
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));
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));
309 Dune::AlwaysFalse<VolVarFunction>::value,
310 "Could not identify the given volume variable as scalar, vector or tensor."
319 template<Detail::Container C>
320 void addField(
const C& values,
const std::string& name)
321 { addField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
327 template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
328 void addField(
const C& values,
const std::string& name,
const GridFormat::Precision<T>& prec)
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.");
338 addCellField(values, name, prec);
340 addPointField(values, name, prec);
346 template<Gr
idFormat::Concepts::Po
intFunction<Gr
idView> DofFunction>
347 void addPointField(DofFunction&& f,
const std::string& name)
348 { this->setPointField(name, std::forward<DofFunction>(f)); }
353 template<Detail::Container C>
354 void addPointField(
const C& values,
const std::string& name)
355 { addPointField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
360 template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
361 void addPointField(
const C& values,
const std::string& name,
const GridFormat::Precision<T>& prec)
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");
366 this->setPointField(name, [&] (
const auto& vertex) {
367 return values[gridVariables_.gridGeometry().vertexMapper().index(vertex)];
374 template<Gr
idFormat::Concepts::CellFunction<Gr
idView> DofFunction>
375 void addCellField(DofFunction&& f,
const std::string& name)
376 { this->setCellField(name, std::forward<DofFunction>(f)); }
381 template<Detail::Container C>
382 void addCellField(
const C& values,
const std::string& name)
383 { addCellField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
388 template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
389 void addCellField(
const C& values,
const std::string& name,
const GridFormat::Precision<T>& prec)
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");
394 this->setCellField(name, [&] (
const auto& element) {
395 return values[gridVariables_.gridGeometry().elementMapper().index(element)];
402 std::string write(
const std::string& name)
404 volVarFields_.updateFieldData(gridVariables_, solutionVector_);
405 auto filename = ParentType::write(name);
406 volVarFields_.clearFieldData();
413 template<std::
floating_po
int T>
414 std::string write(T time)
416 volVarFields_.updateFieldData(gridVariables_, solutionVector_);
417 auto filename = ParentType::write(time);
418 volVarFields_.clearFieldData();
425 volVarFields_.clear();
429 template<
typename ResultType,
typename Id>
430 void setVolVarField_(
const std::string& name, Id&& volVarFieldId)
432 auto dofEntityField = [&, _id=std::move(volVarFieldId)] (
const auto& entity) {
433 return volVarFields_.getValue(_id, gridVariables_.gridGeometry().dofMapper().index(entity));
436 this->setPointField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
438 this->setCellField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
441 const GridVariables& gridVariables_;
442 const SolutionVector& solutionVector_;
443 VolVarFieldStorage volVarFields_;
447template<
typename Gr
idVariables,
typename SolutionVector>
448class OutputModule<GridVariables, SolutionVector>::VolVarFieldStorage
451 { scalar, vector, tensor };
464 std::function<T(
const VolVar&)> getter;
468 template<FieldType ft>
469 struct FieldId { std::size_t index; };
471 template<std::ranges::range R>
472 static constexpr auto toStorageVector(R&& in)
475 std::ranges::copy(in, result.begin());
479 template<Gr
idFormat::Concepts::MDRange<2> R>
480 static constexpr auto toStorageTensor(R&& in)
483 std::ranges::for_each(in, [&, i=0] (
const auto& row)
mutable {
484 std::ranges::copy(row, result[i++].begin());
489 template<FieldType ft>
490 const auto& getValue(
const FieldId<ft>&
id, std::size_t idx)
const
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);
497 return tensorFieldStorage_.at(
id.index).data.at(idx);
500 auto registerScalarField(std::string name, std::function<Scalar(
const VolVar&)> f)
501 {
return register_<FieldType::scalar>(std::move(name), scalarFieldStorage_, std::move(f)); }
503 auto registerVectorField(std::string name, std::function<Vector(
const VolVar&)> f)
504 {
return register_<FieldType::vector>(std::move(name), vectorFieldStorage_, std::move(f)); }
506 auto registerTensorField(std::string name, std::function<Tensor(
const VolVar&)> f)
507 {
return register_<FieldType::tensor>(std::move(name), tensorFieldStorage_, std::move(f)); }
509 void updateFieldData(
const GridVariables& gridVars,
const SolutionVector& x)
511 resizeFieldData_(gridVars.gridGeometry().numDofs());
512 const auto range = GridFormat::cells(gridVars.gridGeometry().gridView());
514#
if __cpp_lib_parallel_algorithm >= 201603L
515 std::execution::par_unseq,
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))
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); }
532 void clearFieldData()
534 for (
auto& s : scalarFieldStorage_) { s.data.clear(); }
535 for (
auto& s : vectorFieldStorage_) { s.data.clear(); }
536 for (
auto& s : tensorFieldStorage_) { s.data.clear(); }
542 scalarFieldStorage_.clear();
543 vectorFieldStorage_.clear();
544 tensorFieldStorage_.clear();
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)
553 if (exists_<ft>(name))
554 DUNE_THROW(Dune::InvalidStateException,
"Volume variables field '" << name <<
"' is already defined.");
556 FieldId<ft>
id{storage.size()};
557 fields_.emplace_back(FieldInfo{std::move(name), ft,
id.index});
558 storage.push_back({{}, std::move(f)});
562 template<FieldType ft>
563 bool exists_(
const std::string& name)
const
565 return std::ranges::any_of(fields_, [&] (
const FieldInfo& info) {
566 return info.type == ft && info.name == name;
570 void resizeFieldData_(std::size_t size)
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); });
577 std::vector<FieldInfo> fields_;
578 std::vector<FieldStorage<Scalar>> scalarFieldStorage_;
579 std::vector<FieldStorage<Vector>> vectorFieldStorage_;
580 std::vector<FieldStorage<Tensor>> tensorFieldStorage_;
589template<
class... Args>
593 template<
class... _Args>
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)."
605template<
class... Args>
609 template<
class... _Args>
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)."
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
The available discretization methods in Dumux.
constexpr None none
Definition: method.hh:153
constexpr bool isCVFE
Definition: method.hh:67
Definition: gridwriter.hh:587