12#ifndef DUMUX_IO_GRID_WRITER_HH 
   13#define DUMUX_IO_GRID_WRITER_HH 
   17#if DUMUX_HAVE_GRIDFORMAT 
   24#include <gridformat/common/type_traits.hpp> 
   25#include <gridformat/gridformat.hpp> 
   26#include <gridformat/traits/dune.hpp> 
   28#include <dune/common/exceptions.hh> 
   29#include <dune/common/timer.hh> 
   41template<
typename Grid, 
typename Format, 
typename Comm, 
typename... Args>
 
   42auto makeParallelWriter(
const Grid& grid, 
const Format& fmt, 
const Dune::Communication<Comm>& comm, Args&&... args)
 
   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)...};
 
   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)...}; }
 
   53template<
typename Grid, 
typename Format, 
typename... Args>
 
   54auto makeWriter(
const Grid& grid, 
const Format& fmt, Args&&... args)
 
   56    const auto& comm = GridFormat::Dune::Traits::GridView<Grid>::get(grid).comm();
 
   57    return makeParallelWriter(grid, fmt, comm, std::forward<Args>(args)...);
 
   61concept Container = 
requires(
const T& t) {
 
   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; }
 
   75    using GridFormat::float32;
 
   76    using GridFormat::float64;
 
   78    using GridFormat::uint64;
 
   79    using GridFormat::uint32;
 
   80    using GridFormat::uint16;
 
   81    using GridFormat::uint8;
 
   83    using GridFormat::int64;
 
   84    using GridFormat::int32;
 
   85    using GridFormat::int16;
 
   86    using GridFormat::int8;
 
   95struct Order { 
static_assert(order > 0, 
"order must be > 0"); };
 
   98inline constexpr auto order = Order<o>{};
 
  109template<Gr
idFormat::Concepts::Gr
id Gr
idView, 
int order = 1>
 
  112    using Grid = std::conditional_t<
 
  114        GridFormat::Dune::LagrangePolynomialGrid<GridView>,
 
  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>;
 
  128    template<
typename Format>
 
  130                        const GridView& gridView,
 
  131                        const Order<order>& = {})
 
  132    : gridView_{gridView}
 
  133    , grid_{makeGrid_(gridView)}
 
  134    , writer_{Detail::makeWriter(grid_, fmt)}
 
  136        if (gridView.comm().size() > 0 && 
getParam<bool>(
"IO.GridWriter.AddProcessRank", 
true))
 
  137            setCellField(
"rank", [&](
const Cell&) { 
return gridView.comm().rank(); }, Precision::uint64);
 
  144    template<
typename Format>
 
  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)}
 
  153        if (gridView.comm().size() > 0 && 
getParam<bool>(
"IO.GridWriter.AddProcessRank", 
true))
 
  154            setCellField(
"rank", [&](
const Cell&) { 
return gridView.comm().rank(); }, Precision::uint64);
 
  161    std::string write(
const std::string& name)
 const 
  162    { 
return writer_.write(name); }
 
  168    template<std::
floating_po
int T>
 
  169    std::string write(T time)
 const 
  170    { 
return writer_.write(time); }
 
  173    template<Gr
idFormat::Concepts::CellFunction<Gr
idView> F,
 
  174             Gr
idFormat::Concepts::Scalar T = Gr
idFormat::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); }
 
  179    template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F>
 
  180    void setCellField(
const std::string& name, F&& f)
 
  181    { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name); }
 
  184    template<
typename F, Gr
idFormat::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); }
 
  189    template<Gr
idFormat::Concepts::Po
intFunction<Gr
idView> F,
 
  190             Gr
idFormat::Concepts::Scalar T = Gr
idFormat::FieldScalar<std::invoke_result_t<F, Vertex>>>
 
  191    void setPointField(
const std::string& name, F&& f, 
const GridFormat::Precision<T>& prec = {})
 
  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);
 
  198    template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F>
 
  199    void setPointField(
const std::string& name, F&& f)
 
  200    { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name); }
 
  203    template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F, Gr
idFormat::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); }
 
  214        if constexpr (order > 1)
 
  215            grid_.update(gridView_);
 
  219    Grid makeGrid_(
const GridView& gv)
 const 
  221        if constexpr (order > 1)
 
  222            return Grid{gv, order};
 
  236template<
typename Gr
idVariables, 
typename SolutionVector>
 
  238    using ParentType = GridWriter<typename GridVariables::GridGeometry::GridView, 1>;
 
  239    using GridView = 
typename GridVariables::GridGeometry::GridView;
 
  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;
 
  248    class VolVarFieldStorage;
 
  250    static constexpr int defaultVerbosity_ = 1;
 
  253    using VolumeVariables = VolVar;
 
  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
 
  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}
 
  273        setVerbosity(defaultVerbosity_);
 
  280    template<
typename Format>
 
  282                          const GridVariables& gridVariables,
 
  283                          const SolutionVector& sol)
 
  284    : ParentType{fmt, gridVariables.gridGeometry().gridView(), order<1>}
 
  285    , gridVariables_{gridVariables}
 
  286    , solutionVector_{sol}
 
  288        setVerbosity(defaultVerbosity_);
 
  295    template<
typename Format>
 
  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}
 
  304        setVerbosity(defaultVerbosity_);
 
  310    template<std::invocable<const VolumeVariables&> VolVarFunction>
 
  311    void addVolumeVariable(VolVarFunction&& f, 
const std::string& name)
 
  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));
 
  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));
 
  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));
 
  329                Dune::AlwaysFalse<VolVarFunction>::value,
 
  330                "Could not identify the given volume variable as scalar, vector or tensor." 
  339    template<Detail::Container C>
 
  340    void addField(
const C& values, 
const std::string& name)
 
  341    { addField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
 
  347    template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
 
  348    void addField(
const C& values, 
const std::string& name, 
const GridFormat::Precision<T>& prec)
 
  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.");
 
  358            addCellField(values, name, prec);
 
  360            addPointField(values, name, prec);
 
  366    template<Gr
idFormat::Concepts::Po
intFunction<Gr
idView> DofFunction>
 
  367    void addPointField(DofFunction&& f, 
const std::string& name)
 
  368    { this->setPointField(name, std::forward<DofFunction>(f)); }
 
  373    template<Detail::Container C>
 
  374    void addPointField(
const C& values, 
const std::string& name)
 
  375    { addPointField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
 
  380    template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
 
  381    void addPointField(
const C& values, 
const std::string& name, 
const GridFormat::Precision<T>& prec)
 
  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");
 
  386        addPointField_(values, name, prec);
 
  392    template<Gr
idFormat::Concepts::CellFunction<Gr
idView> DofFunction>
 
  393    void addCellField(DofFunction&& f, 
const std::string& name)
 
  394    { this->setCellField(name, std::forward<DofFunction>(f)); }
 
  399    template<Detail::Container C>
 
  400    void addCellField(
const C& values, 
const std::string& name)
 
  401    { addCellField(values, name, GridFormat::Precision<GridFormat::MDRangeScalar<C>>{}); }
 
  406    template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
 
  407    void addCellField(
const C& values, 
const std::string& name, 
const GridFormat::Precision<T>& prec)
 
  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");
 
  412        addCellField_(values, name, prec);
 
  418    std::string write(
const std::string& name)
 
  422        volVarFields_.updateFieldData(gridVariables_, solutionVector_);
 
  423        auto filename = ParentType::write(name);
 
  424        volVarFields_.clearFieldData();
 
  428            std::cout << Fmt::format(
 
  429                "Writing output to \"{}\". Took {:.2g} seconds.\n",
 
  430                filename, timer.elapsed()
 
  439    template<std::
floating_po
int T>
 
  440    std::string write(T time)
 
  444        volVarFields_.updateFieldData(gridVariables_, solutionVector_);
 
  445        auto filename = ParentType::write(time);
 
  446        volVarFields_.clearFieldData();
 
  450            std::cout << Fmt::format(
 
  451                "Writing output to \"{}\". Took {:.2g} seconds.\n",
 
  452                filename, timer.elapsed()
 
  461        volVarFields_.clear();
 
  464    void setVerbosity(
int verbosity)
 
  466        if (gridVariables_.gridGeometry().gridView().comm().rank() == 0)
 
  467            verbosity_ = verbosity;
 
  472    template<Detail::Container C, Gr
idFormat::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)
 
  476        this->setPointField(name, [&] (
const auto& vertex) -> T {
 
  477            return values[gridVariables_.gridGeometry().vertexMapper().index(vertex)][0];
 
  481    template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
 
  482    void addPointField_(
const C& values, 
const std::string& name, 
const GridFormat::Precision<T>& prec)
 
  484        this->setPointField(name, [&] (
const auto& vertex) -> std::ranges::range_value_t<C> {
 
  485            return values[gridVariables_.gridGeometry().vertexMapper().index(vertex)];
 
  490    template<Detail::Container C, Gr
idFormat::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)
 
  494        this->setCellField(name, [&] (
const auto& element) -> T {
 
  495            return values[gridVariables_.gridGeometry().elementMapper().index(element)][0];
 
  499    template<Detail::Container C, Gr
idFormat::Concepts::Scalar T>
 
  500    void addCellField_(
const C& values, 
const std::string& name, 
const GridFormat::Precision<T>& prec)
 
  502        this->setCellField(name, [&] (
const auto& element) -> std::ranges::range_value_t<C> {
 
  503            return values[gridVariables_.gridGeometry().elementMapper().index(element)];
 
  507    template<
typename ResultType, 
typename Id>
 
  508    void setVolVarField_(
const std::string& name, Id&& volVarFieldId)
 
  510        auto dofEntityField = [&, _id=std::move(volVarFieldId)] (
const auto& entity) {
 
  511            return volVarFields_.getValue(_id, gridVariables_.gridGeometry().dofMapper().index(entity));
 
  514            this->setPointField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
 
  516            this->setCellField(name, std::move(dofEntityField), GridFormat::Precision<ResultType>{});
 
  519    const GridVariables& gridVariables_;
 
  520    const SolutionVector& solutionVector_;
 
  521    VolVarFieldStorage volVarFields_;
 
  526template<
typename Gr
idVariables, 
typename SolutionVector>
 
  527class OutputModule<GridVariables, SolutionVector>::VolVarFieldStorage
 
  530    { scalar, vector, tensor };
 
  543        std::function<T(
const VolVar&)> getter;
 
  547    template<FieldType ft>
 
  548    struct FieldId { std::size_t index; };
 
  550    template<std::ranges::range R>
 
  551    static constexpr auto toStorageVector(R&& in)
 
  554        std::ranges::copy(in, result.begin());
 
  558    template<Gr
idFormat::Concepts::MDRange<2> R>
 
  559    static constexpr auto toStorageTensor(R&& in)
 
  562        std::ranges::for_each(in, [&, i=0] (
const auto& row) 
mutable {
 
  563            std::ranges::copy(row, result[i++].begin());
 
  568    template<FieldType ft>
 
  569    const auto& getValue(
const FieldId<ft>& 
id, std::size_t idx)
 const 
  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);
 
  576            return tensorFieldStorage_.at(
id.index).data.at(idx);
 
  579    auto registerScalarField(std::string name, std::function<Scalar(
const VolVar&)> f)
 
  580    { 
return register_<FieldType::scalar>(std::move(name), scalarFieldStorage_, std::move(f)); }
 
  582    auto registerVectorField(std::string name, std::function<Vector(
const VolVar&)> f)
 
  583    { 
return register_<FieldType::vector>(std::move(name), vectorFieldStorage_, std::move(f)); }
 
  585    auto registerTensorField(std::string name, std::function<Tensor(
const VolVar&)> f)
 
  586    { 
return register_<FieldType::tensor>(std::move(name), tensorFieldStorage_, std::move(f)); }
 
  588    void updateFieldData(
const GridVariables& gridVars, 
const SolutionVector& x)
 
  590        resizeFieldData_(gridVars.gridGeometry().numDofs());
 
  591        const auto range = GridFormat::cells(gridVars.gridGeometry().gridView());
 
  593#
if __cpp_lib_parallel_algorithm >= 201603L
 
  594            std::execution::par_unseq,
 
  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))
 
  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); }
 
  611    void clearFieldData()
 
  613        for (
auto& s : scalarFieldStorage_) { s.data.clear(); }
 
  614        for (
auto& s : vectorFieldStorage_) { s.data.clear(); }
 
  615        for (
auto& s : tensorFieldStorage_) { s.data.clear(); }
 
  621        scalarFieldStorage_.clear();
 
  622        vectorFieldStorage_.clear();
 
  623        tensorFieldStorage_.clear();
 
  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)
 
  632        if (exists_<ft>(name))
 
  633            DUNE_THROW(Dune::InvalidStateException, 
"Volume variables field '" << name << 
"' is already defined.");
 
  635        FieldId<ft> 
id{storage.size()};
 
  636        fields_.emplace_back(FieldInfo{std::move(name), ft, 
id.index});
 
  637        storage.push_back({{}, std::move(f)});
 
  641    template<FieldType ft>
 
  642    bool exists_(
const std::string& name)
 const 
  644        return std::ranges::any_of(fields_, [&] (
const FieldInfo& info) {
 
  645            return info.type == ft && info.name == name;
 
  649    void resizeFieldData_(std::size_t size)
 
  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); });
 
  656    std::vector<FieldInfo> fields_;
 
  657    std::vector<FieldStorage<Scalar>> scalarFieldStorage_;
 
  658    std::vector<FieldStorage<Vector>> vectorFieldStorage_;
 
  659    std::vector<FieldStorage<Tensor>> tensorFieldStorage_;
 
  668template<
class... Args>
 
  672    template<
class... _Args>
 
  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)." 
 
 
  684template<
class... Args>
 
  688    template<
class... _Args>
 
  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)." 
 
 
 
Definition gridwriter.hh:670
GridWriter(_Args &&...)
Definition gridwriter.hh:673
Definition gridwriter.hh:686
OutputModule(_Args &&...)
Definition gridwriter.hh:689
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition localview.hh:26
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.