12#ifndef DUMUX_IO_GRID_WRITER_HH
13#define DUMUX_IO_GRID_WRITER_HH
17#if DUMUX_HAVE_GRIDFORMAT
21#include <gridformat/gridformat.hpp>
22#include <gridformat/traits/dune.hpp>
30template<
typename Grid,
typename Format,
typename Comm,
typename... Args>
31auto makeParallelWriter(
const Grid& grid,
const Format& fmt,
const Dune::Communication<Comm>& comm, Args&&... args)
33 return comm.size() > 1
34 ? GridFormat::Writer<Grid>{fmt, grid,
static_cast<Comm
>(comm), std::forward<Args>(args)...}
35 : GridFormat::Writer<Grid>{fmt, grid, std::forward<Args>(args)...};
38template<
typename Grid,
typename Format,
typename... Args>
39auto makeParallelWriter(
const Grid& grid,
const Format& fmt,
const Dune::Communication<Dune::No_Comm>&, Args&&... args)
40{
return GridFormat::Writer<Grid>{fmt, grid, std::forward<Args>(args)...}; }
42template<
typename Grid,
typename Format,
typename... Args>
43auto makeWriter(
const Grid& grid,
const Format& fmt, Args&&... args)
45 const auto& comm = GridFormat::Dune::Traits::GridView<Grid>::get(grid).comm();
46 return makeParallelWriter(grid, fmt, comm, std::forward<Args>(args)...);
53namespace VTK {
using namespace GridFormat::VTK; }
54namespace Format {
using namespace GridFormat::Formats; }
55namespace Encoding {
using namespace GridFormat::Encoding; }
56namespace Compression {
using namespace GridFormat::Compression;
using GridFormat::none; }
58 using GridFormat::float32;
59 using GridFormat::float64;
61 using GridFormat::uint64;
62 using GridFormat::uint32;
63 using GridFormat::uint16;
64 using GridFormat::uint8;
66 using GridFormat::int64;
67 using GridFormat::int32;
68 using GridFormat::int16;
69 using GridFormat::int8;
78struct Order {
static_assert(order > 0,
"order must be > 0"); };
81inline constexpr auto order = Order<o>{};
92template<Gr
idFormat::Concepts::Gr
id Gr
idView,
int order = 1>
95 using Grid = std::conditional_t<
97 GridFormat::Dune::LagrangePolynomialGrid<GridView>,
100 using Cell = GridFormat::Cell<Grid>;
101 using Vertex =
typename GridView::template Codim<GridView::dimension>::Entity;
102 using Element =
typename GridView::template Codim<0>::Entity;
103 using Coordinate =
typename Element::Geometry::GlobalCoordinate;
104 using Writer = GridFormat::Writer<Grid>;
111 template<
typename Format>
113 const GridView& gridView,
114 const Order<order>& = {})
115 : gridView_{gridView}
116 , grid_{makeGrid_(gridView)}
117 , writer_{Detail::makeWriter(grid_, fmt)}
118 { writer_.set_meta_data(
"rank", gridView_.comm().rank()); }
124 template<
typename Format>
126 const GridView& gridView,
127 const std::string& filename,
128 const Order<order>& = {})
129 : gridView_{gridView}
130 , grid_{makeGrid_(gridView)}
131 , writer_{Detail::makeWriter(grid_, fmt, filename)}
132 { writer_.set_meta_data(
"rank", gridView_.comm().rank()); }
139 std::string write(
const std::string& name)
const
140 {
return writer_.write(name); }
147 template<std::
floating_po
int T>
148 std::string write(T time)
const
149 {
return writer_.write(time); }
152 template<Gr
idFormat::Concepts::CellFunction<Gr
idView> F,
153 Gr
idFormat::Concepts::Scalar T = Gr
idFormat::FieldScalar<std::invoke_result_t<F, Element>>>
154 void setCellField(
const std::string& name, F&& f,
const GridFormat::Precision<T>& prec = {})
155 { writer_.set_cell_field(name, std::move(f), prec); }
158 template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F>
159 void setCellField(
const std::string& name, F&& f)
160 { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name); }
163 template<
typename F, Gr
idFormat::Concepts::Scalar T>
164 void setCellField(
const std::string& name, F&& f,
const GridFormat::Precision<T>& prec)
165 { GridFormat::Dune::set_cell_function(std::forward<F>(f), writer_, name, prec); }
168 template<Gr
idFormat::Concepts::Po
intFunction<Gr
idView> F,
169 Gr
idFormat::Concepts::Scalar T = Gr
idFormat::FieldScalar<std::invoke_result_t<F, Vertex>>>
170 void setPointField(
const std::string& name, F&& f,
const GridFormat::Precision<T>& prec = {})
172 static_assert(order == 1,
"Point lambdas can only be used for order == 1. Use Dune::Functions instead.");
173 writer_.set_point_field(name, std::move(f), prec);
177 template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F>
178 void setPointField(
const std::string& name, F&& f)
179 { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name); }
182 template<Gr
idFormat::Dune::Concepts::Function<Gr
idView> F, Gr
idFormat::Concepts::Scalar T>
183 void setPointField(
const std::string& name, F&& f,
const GridFormat::Precision<T>& prec)
184 { GridFormat::Dune::set_point_function(std::forward<F>(f), writer_, name, prec); }
193 if constexpr (order > 1)
194 grid_.update(gridView_);
198 Grid makeGrid_(
const GridView& gv)
const
200 if constexpr (order > 1)
201 return Grid{gv, order};
217template<
class... Args>
221 template<
class... _Args>
226 "GridWriter only available when the GridFormat library is available. "
227 "Use `git submodule init && git submodule update` to pull it."
Definition: gridwriter.hh:219
GridWriter(_Args &&...)
Definition: gridwriter.hh:222
constexpr None none
Definition: method.hh:153
Definition: gridwriter.hh:215