version 3.10.0
Loading...
Searching...
No Matches
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-FileCopyrightText: 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 <concepts>
22#include <type_traits>
23
24#include <gridformat/common/type_traits.hpp>
25#include <gridformat/gridformat.hpp>
26#include <gridformat/traits/dune.hpp>
27
28#include <dune/common/exceptions.hh>
29#include <dune/common/timer.hh>
30
31#include <dumux/io/format.hh>
34
35namespace Dumux::IO {
36
37
38#ifndef DOXYGEN
39namespace Detail {
40
41template<typename Grid, typename Format, typename Comm, typename... Args>
42auto makeParallelWriter(const Grid& grid, const Format& fmt, const Dune::Communication<Comm>& comm, Args&&... args)
43{
44 return comm.size() > 1
45 ? GridFormat::Writer<Grid>{fmt, grid, static_cast<Comm>(comm), std::forward<Args>(args)...}
46 : GridFormat::Writer<Grid>{fmt, grid, std::forward<Args>(args)...};
47}
48
49template<typename Grid, typename Format, typename... Args>
50auto makeParallelWriter(const Grid& grid, const Format& fmt, const Dune::Communication<Dune::No_Comm>&, Args&&... args)
51{ return GridFormat::Writer<Grid>{fmt, grid, std::forward<Args>(args)...}; }
52
53template<typename Grid, typename Format, typename... Args>
54auto makeWriter(const Grid& grid, const Format& fmt, Args&&... args)
55{
56 const auto& comm = GridFormat::Dune::Traits::GridView<Grid>::get(grid).comm();
57 return makeParallelWriter(grid, fmt, comm, std::forward<Args>(args)...);
58}
59
60template<typename T>
61concept Container = requires(const T& t) {
62 { t.size() };
63 { t[std::size_t{}] };
64};
65
66} // namespace Detail
67#endif // DOXYGEN
68
69
70namespace VTK { using namespace GridFormat::VTK; }
71namespace Format { using namespace GridFormat::Formats; }
72namespace Encoding { using namespace GridFormat::Encoding; }
73namespace Compression { using namespace GridFormat::Compression; using GridFormat::none; }
74namespace Precision {
75 using GridFormat::float32;
76 using GridFormat::float64;
77
78 using GridFormat::uint64;
79 using GridFormat::uint32;
80 using GridFormat::uint16;
81 using GridFormat::uint8;
82
83 using GridFormat::int64;
84 using GridFormat::int32;
85 using GridFormat::int16;
86 using GridFormat::int8;
87} // namespace Precision
88
89
94template<int order>
95struct Order { static_assert(order > 0, "order must be > 0"); };
96
97template<int o>
98inline constexpr auto order = Order<o>{};
99
109template<GridFormat::Concepts::Grid GridView, int order = 1>
110class GridWriter
111{
112 using Grid = std::conditional_t<
113 (order > 1),
114 GridFormat::Dune::LagrangePolynomialGrid<GridView>,
115 GridView
116 >;
117 using Cell = GridFormat::Cell<Grid>;
118 using Vertex = typename GridView::template Codim<GridView::dimension>::Entity;
119 using Element = typename GridView::template Codim<0>::Entity;
120 using Coordinate = typename Element::Geometry::GlobalCoordinate;
121 using Writer = GridFormat::Writer<Grid>;
122
123 public:
128 template<typename Format>
129 explicit GridWriter(const Format& fmt,
130 const GridView& gridView,
131 const Order<order>& = {})
132 : gridView_{gridView}
133 , grid_{makeGrid_(gridView)}
134 , writer_{Detail::makeWriter(grid_, fmt)}
135 {
136 if (gridView.comm().size() > 0 && getParam<bool>("IO.GridWriter.AddProcessRank", true))
137 setCellField("rank", [&](const Cell&) { return gridView.comm().rank(); }, Precision::uint64);
138 }
139
144 template<typename Format>
145 explicit GridWriter(const Format& fmt,
146 const GridView& gridView,
147 const std::string& filename,
148 const Order<order>& = {})
149 : gridView_{gridView}
150 , grid_{makeGrid_(gridView)}
151 , writer_{Detail::makeWriter(grid_, fmt, filename)}
152 {
153 if (gridView.comm().size() > 0 && getParam<bool>("IO.GridWriter.AddProcessRank", true))
154 setCellField("rank", [&](const Cell&) { return gridView.comm().rank(); }, Precision::uint64);
155 }
156
161 std::string write(const std::string& name) const
162 { return writer_.write(name); }
163
168 template<std::floating_point T>
169 std::string write(T time) const
170 { return writer_.write(time); }
171
173 template<GridFormat::Concepts::CellFunction<GridView> F,
174 GridFormat::Concepts::Scalar T = GridFormat::FieldScalar<std::invoke_result_t<F, Element>>>
175 void setCellField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec = {})
176 { writer_.set_cell_field(name, std::move(f), prec); }
177
179 template<GridFormat::Dune::Concepts::Function<GridView> F>
180 void setCellField(const std::string& name, F&& f)
181 { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name); }
182
184 template<typename F, GridFormat::Concepts::Scalar T>
185 void setCellField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec)
186 { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name, prec); }
187
189 template<GridFormat::Concepts::PointFunction<GridView> F,
190 GridFormat::Concepts::Scalar T = GridFormat::FieldScalar<std::invoke_result_t<F, Vertex>>>
191 void setPointField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec = {})
192 {
193 static_assert(order == 1, "Point lambdas can only be used for order == 1. Use Dune::Functions instead.");
194 writer_.set_point_field(name, std::move(f), prec);
195 }
196
198 template<GridFormat::Dune::Concepts::Function<GridView> F>
199 void setPointField(const std::string& name, F&& f)
200 { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name); }
201
203 template<GridFormat::Dune::Concepts::Function<GridView> F, GridFormat::Concepts::Scalar T>
204 void setPointField(const std::string& name, F&& f, const GridFormat::Precision<T>& prec)
205 { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name, prec); }
206
208 void clear()
209 { writer_.clear(); }
210
212 void update()
213 {
214 if constexpr (order > 1)
215 grid_.update(gridView_);
216 }
217
218 private:
219 Grid makeGrid_(const GridView& gv) const
220 {
221 if constexpr (order > 1)
222 return Grid{gv, order};
223 else
224 return gv;
225 }
226
227 GridView gridView_;
228 Grid grid_;
229 Writer writer_;
230};
231
236template<typename GridVariables, typename SolutionVector>
237class OutputModule : private GridWriter<typename GridVariables::GridGeometry::GridView, 1> {
238 using ParentType = GridWriter<typename GridVariables::GridGeometry::GridView, 1>;
239 using GridView = typename GridVariables::GridGeometry::GridView;
240
242 static constexpr int dimWorld = GridVariables::GridGeometry::GridView::dimensionworld;
243 using Scalar = typename GridVariables::Scalar;
244 using Vector = Dune::FieldVector<Scalar, dimWorld>;
245 using Tensor = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
246 using VolVar = typename GridVariables::VolumeVariables;
247
248 class VolVarFieldStorage;
249
250 static constexpr int defaultVerbosity_ = 1;
251
252 public:
253 using VolumeVariables = VolVar;
254
255 static constexpr auto defaultFileFormat = IO::Format::pvd_with(
256 IO::Format::vtu.with({
257 .encoder = IO::Encoding::ascii,
258 .compressor = IO::Compression::none,
259 .data_format = VTK::DataFormat::inlined
260 })
261 );
262
266 explicit OutputModule(const GridVariables& gridVariables,
267 const SolutionVector& sol,
268 const std::string& filename)
269 : ParentType{defaultFileFormat, gridVariables.gridGeometry().gridView(), filename, order<1>}
270 , gridVariables_{gridVariables}
271 , solutionVector_{sol}
272 {
273 setVerbosity(defaultVerbosity_);
274 }
275
280 template<typename Format>
281 explicit OutputModule(const Format& fmt,
282 const GridVariables& gridVariables,
283 const SolutionVector& sol)
284 : ParentType{fmt, gridVariables.gridGeometry().gridView(), order<1>}
285 , gridVariables_{gridVariables}
286 , solutionVector_{sol}
287 {
288 setVerbosity(defaultVerbosity_);
289 }
290
295 template<typename Format>
296 explicit OutputModule(const Format& fmt,
297 const GridVariables& gridVariables,
298 const SolutionVector& sol,
299 const std::string& filename)
300 : ParentType{fmt, gridVariables.gridGeometry().gridView(), filename, order<1>}
301 , gridVariables_{gridVariables}
302 , solutionVector_{sol}
303 {
304 setVerbosity(defaultVerbosity_);
305 }
306
310 template<std::invocable<const VolumeVariables&> VolVarFunction>
311 void addVolumeVariable(VolVarFunction&& f, const std::string& name)
312 {
313 using ResultType = std::invoke_result_t<std::remove_cvref_t<VolVarFunction>, const VolumeVariables&>;
314 if constexpr (GridFormat::Concepts::Scalar<ResultType>)
315 setVolVarField_<ResultType>(name, volVarFields_.registerScalarField(name, [_f=std::move(f)] (const auto& vv) {
316 return static_cast<Scalar>(_f(vv));
317 }));
318 else if constexpr (GridFormat::mdrange_dimension<ResultType> == 1)
319 setVolVarField_<GridFormat::MDRangeScalar<ResultType>>(name, volVarFields_.registerVectorField(name, [_f=std::move(f)] (const auto& vv) {
320 return VolVarFieldStorage::toStorageVector(_f(vv));
321 }));
322 else if constexpr (GridFormat::mdrange_dimension<ResultType> == 2)
323 setVolVarField_<GridFormat::MDRangeScalar<ResultType>>(name, volVarFields_.registerTensorField(name, [_f=std::move(f)] (const auto& vv) {
324 return VolVarFieldStorage::toStorageTensor(_f(vv));
325 }));
326 else
327 {
328 static_assert(
329 Dune::AlwaysFalse<VolVarFunction>::value,
330 "Could not identify the given volume variable as scalar, vector or tensor."
331 );
332 }
333 }
334
339 template<Detail::Container C>
340 void addField(const C& values, const std::string& name)
341 { addField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
342
347 template<Detail::Container C, GridFormat::Concepts::Scalar T>
348 void addField(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
349 {
350 const bool hasCellSize = values.size() == gridVariables_.gridGeometry().elementMapper().size();
351 const bool hasPointSize = values.size() == gridVariables_.gridGeometry().vertexMapper().size();
352 if (hasCellSize && hasPointSize)
353 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of field type failed. Please use addCellField or addPointField instead.");
354 if (!hasCellSize && !hasPointSize)
355 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of field type failed. Given container size does not match neither the number of points nor cells.");
356
357 if (hasCellSize)
358 addCellField(values, name, prec);
359 else
360 addPointField(values, name, prec);
361 }
362
366 template<GridFormat::Concepts::PointFunction<GridView> DofFunction>
367 void addPointField(DofFunction&& f, const std::string& name)
368 { this->setPointField(name, std::forward<DofFunction>(f)); }
369
373 template<Detail::Container C>
374 void addPointField(const C& values, const std::string& name)
375 { addPointField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
376
380 template<Detail::Container C, GridFormat::Concepts::Scalar T>
381 void addPointField(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
382 {
383 if (values.size() != gridVariables_.gridGeometry().vertexMapper().size())
384 DUNE_THROW(Dune::InvalidStateException, "Given container does not match the number of points in the grid");
385
386 addPointField_(values, name, prec);
387 }
388
392 template<GridFormat::Concepts::CellFunction<GridView> DofFunction>
393 void addCellField(DofFunction&& f, const std::string& name)
394 { this->setCellField(name, std::forward<DofFunction>(f)); }
395
399 template<Detail::Container C>
400 void addCellField(const C& values, const std::string& name)
401 { addCellField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
402
406 template<Detail::Container C, GridFormat::Concepts::Scalar T>
407 void addCellField(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
408 {
409 if (values.size() != gridVariables_.gridGeometry().elementMapper().size())
410 DUNE_THROW(Dune::InvalidStateException, "Given container does not match the number of cells in the grid");
411
412 addCellField_(values, name, prec);
413 }
414
418 std::string write(const std::string& name)
419 {
420 Dune::Timer timer;
421
422 volVarFields_.updateFieldData(gridVariables_, solutionVector_);
423 auto filename = ParentType::write(name);
424 volVarFields_.clearFieldData();
425
426 timer.stop();
427 if (verbosity_ > 0)
428 std::cout << Fmt::format(
429 "Writing output to \"{}\". Took {:.2g} seconds.\n",
430 filename, timer.elapsed()
431 );
432
433 return filename;
434 }
435
439 template<std::floating_point T>
440 std::string write(T time)
441 {
442 Dune::Timer timer;
443
444 volVarFields_.updateFieldData(gridVariables_, solutionVector_);
445 auto filename = ParentType::write(time);
446 volVarFields_.clearFieldData();
447
448 timer.stop();
449 if (verbosity_ > 0)
450 std::cout << Fmt::format(
451 "Writing output to \"{}\". Took {:.2g} seconds.\n",
452 filename, timer.elapsed()
453 );
454
455 return filename;
456 }
457
459 void clear() {
460 ParentType::clear();
461 volVarFields_.clear();
462 }
463
464 void setVerbosity(int verbosity)
465 {
466 if (gridVariables_.gridGeometry().gridView().comm().rank() == 0)
467 verbosity_ = verbosity;
468 }
469
470private:
471 // reproduce behaviour of VtkOutputModule that flattens vectors of size 1
472 template<Detail::Container C, GridFormat::Concepts::Scalar T>
473 requires(GridFormat::has_sub_range<C> && std::ranges::range_value_t<C>::size() == 1)
474 void addPointField_(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
475 {
476 this->setPointField(name, [&] (const auto& vertex) -> T {
477 return values[gridVariables_.gridGeometry().vertexMapper().index(vertex)][0];
478 }, prec);
479 }
480
481 template<Detail::Container C, GridFormat::Concepts::Scalar T>
482 void addPointField_(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
483 {
484 this->setPointField(name, [&] (const auto& vertex) -> std::ranges::range_value_t<C> {
485 return values[gridVariables_.gridGeometry().vertexMapper().index(vertex)];
486 }, prec);
487 }
488
489 // reproduce behaviour of VtkOutputModule that flattens vectors of size 1
490 template<Detail::Container C, GridFormat::Concepts::Scalar T>
491 requires(GridFormat::has_sub_range<C> && std::ranges::range_value_t<C>::size() == 1)
492 void addCellField_(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
493 {
494 this->setCellField(name, [&] (const auto& element) -> T {
495 return values[gridVariables_.gridGeometry().elementMapper().index(element)][0];
496 }, prec);
497 }
498
499 template<Detail::Container C, GridFormat::Concepts::Scalar T>
500 void addCellField_(const C& values, const std::string& name, const GridFormat::Precision<T>& prec)
501 {
502 this->setCellField(name, [&] (const auto& element) -> std::ranges::range_value_t<C> {
503 return values[gridVariables_.gridGeometry().elementMapper().index(element)];
504 }, prec);
505 }
506
507 template<typename ResultType, typename Id>
508 void setVolVarField_(const std::string& name, Id&& volVarFieldId)
509 {
510 auto dofEntityField = [&, _id=std::move(volVarFieldId)] (const auto& entity) {
511 return volVarFields_.getValue(_id, gridVariables_.gridGeometry().dofMapper().index(entity));
512 };
513 if constexpr (isCVFE)
514 this->setPointField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
515 else
516 this->setCellField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
517 }
518
519 const GridVariables& gridVariables_;
520 const SolutionVector& solutionVector_;
521 VolVarFieldStorage volVarFields_;
522 int verbosity_ = 0;
523};
524
525// Class to store vol var fields; implementation detail of the OutputModule class
526template<typename GridVariables, typename SolutionVector>
527class OutputModule<GridVariables, SolutionVector>::VolVarFieldStorage
528{
529 enum class FieldType
530 { scalar, vector, tensor };
531
532 struct FieldInfo
533 {
534 std::string name;
535 FieldType type;
536 std::size_t index;
537 };
538
539 template<typename T>
540 struct FieldStorage
541 {
542 std::vector<T> data;
543 std::function<T(const VolVar&)> getter;
544 };
545
546public:
547 template<FieldType ft>
548 struct FieldId { std::size_t index; };
549
550 template<std::ranges::range R>
551 static constexpr auto toStorageVector(R&& in)
552 {
553 Vector result;
554 std::ranges::copy(in, result.begin());
555 return result;
556 }
557
558 template<GridFormat::Concepts::MDRange<2> R>
559 static constexpr auto toStorageTensor(R&& in)
560 {
561 Tensor result;
562 std::ranges::for_each(in, [&, i=0] (const auto& row) mutable {
563 std::ranges::copy(row, result[i++].begin());
564 });
565 return result;
566 }
567
568 template<FieldType ft>
569 const auto& getValue(const FieldId<ft>& id, std::size_t idx) const
570 {
571 if constexpr (ft == FieldType::scalar)
572 return scalarFieldStorage_.at(id.index).data.at(idx);
573 else if constexpr (ft == FieldType::vector)
574 return vectorFieldStorage_.at(id.index).data.at(idx);
575 else
576 return tensorFieldStorage_.at(id.index).data.at(idx);
577 }
578
579 auto registerScalarField(std::string name, std::function<Scalar(const VolVar&)> f)
580 { return register_<FieldType::scalar>(std::move(name), scalarFieldStorage_, std::move(f)); }
581
582 auto registerVectorField(std::string name, std::function<Vector(const VolVar&)> f)
583 { return register_<FieldType::vector>(std::move(name), vectorFieldStorage_, std::move(f)); }
584
585 auto registerTensorField(std::string name, std::function<Tensor(const VolVar&)> f)
586 { return register_<FieldType::tensor>(std::move(name), tensorFieldStorage_, std::move(f)); }
587
588 void updateFieldData(const GridVariables& gridVars, const SolutionVector& x)
589 {
590 resizeFieldData_(gridVars.gridGeometry().numDofs());
591 const auto range = GridFormat::cells(gridVars.gridGeometry().gridView());
592 std::for_each(
593#if __cpp_lib_parallel_algorithm >= 201603L
594 std::execution::par_unseq,
595#endif
596 std::ranges::begin(range),
597 std::ranges::end(range),
598 [&] (const auto& element) {
599 auto fvGeometry = localView(gridVars.gridGeometry()).bindElement(element);
600 auto elemVolVars = localView(gridVars.curGridVolVars()).bindElement(element, fvGeometry, x);
601 for (const auto& scv : scvs(fvGeometry))
602 {
603 const auto& volVars = elemVolVars[scv];
604 for (auto& s : scalarFieldStorage_) { s.data.at(scv.dofIndex()) = s.getter(volVars); }
605 for (auto& s : vectorFieldStorage_) { s.data.at(scv.dofIndex()) = s.getter(volVars); }
606 for (auto& s : tensorFieldStorage_) { s.data.at(scv.dofIndex()) = s.getter(volVars); }
607 }
608 });
609 }
610
611 void clearFieldData()
612 {
613 for (auto& s : scalarFieldStorage_) { s.data.clear(); }
614 for (auto& s : vectorFieldStorage_) { s.data.clear(); }
615 for (auto& s : tensorFieldStorage_) { s.data.clear(); }
616 }
617
618 void clear()
619 {
620 fields_.clear();
621 scalarFieldStorage_.clear();
622 vectorFieldStorage_.clear();
623 tensorFieldStorage_.clear();
624 }
625
626private:
627 template<FieldType ft, typename T>
628 auto register_(std::string&& name,
629 std::vector<FieldStorage<T>>& storage,
630 std::function<T(const VolVar&)>&& f)
631 {
632 if (exists_<ft>(name))
633 DUNE_THROW(Dune::InvalidStateException, "Volume variables field '" << name << "' is already defined.");
634
635 FieldId<ft> id{storage.size()};
636 fields_.emplace_back(FieldInfo{std::move(name), ft, id.index});
637 storage.push_back({{}, std::move(f)});
638 return id;
639 }
640
641 template<FieldType ft>
642 bool exists_(const std::string& name) const
643 {
644 return std::ranges::any_of(fields_, [&] (const FieldInfo& info) {
645 return info.type == ft && info.name == name;
646 });
647 }
648
649 void resizeFieldData_(std::size_t size)
650 {
651 std::ranges::for_each(scalarFieldStorage_, [&] (auto& s) { s.data.resize(size); });
652 std::ranges::for_each(vectorFieldStorage_, [&] (auto& s) { s.data.resize(size); });
653 std::ranges::for_each(tensorFieldStorage_, [&] (auto& s) { s.data.resize(size); });
654 }
655
656 std::vector<FieldInfo> fields_;
657 std::vector<FieldStorage<Scalar>> scalarFieldStorage_;
658 std::vector<FieldStorage<Vector>> vectorFieldStorage_;
659 std::vector<FieldStorage<Tensor>> tensorFieldStorage_;
660};
661
662} // namespace Dumux::IO
663
664#else // DUMUX_HAVE_GRIDFORMAT
665
666namespace Dumux::IO {
667
668template<class... Args>
670{
671public:
672 template<class... _Args>
673 GridWriter(_Args&&...)
674 {
675 static_assert(
676 false,
677 "GridWriter only available when the GridFormat library is available. "
678 "Use `git submodule update --init` to pull it and reconfigure the project "
679 "(note: C++20 is required)."
680 );
681 }
682};
683
684template<class... Args>
686{
687public:
688 template<class... _Args>
689 OutputModule(_Args&&...)
690 {
691 static_assert(
692 false,
693 "OutputModule only available when the GridFormat library is available. "
694 "Use `git submodule update --init` to pull it and reconfigure the project "
695 "(note: C++20 is required)."
696 );
697 }
698};
699
700} // namespace Dumux::IO
701
702#endif // DUMUX_HAVE_GRIDFORMAT
703#endif // DUMUX_IO_GRID_HH
Definition gridwriter.hh:670
GridWriter(_Args &&...)
Definition gridwriter.hh:673
Definition gridwriter.hh:686
OutputModule(_Args &&...)
Definition gridwriter.hh:689
Formatting based on the fmt-library which implements std::format of C++20.
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
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:139
The available discretization methods in Dumux.
constexpr bool isCVFE
Definition method.hh:67
Definition gridwriter.hh:666
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.