1#ifndef DUNE_COPASI_MODEL_MODEL_HH
2#define DUNE_COPASI_MODEL_MODEL_HH
6#include <dune/pdelab/operator/operator.hh>
8#include <dune/functions/gridfunctions/gridviewfunction.hh>
10#include <dune/grid/concepts/grid.hh>
11#include <dune/grid/concepts/gridview.hh>
13#include <dune/common/parametertree.hh>
20#include <unordered_map>
24template<Dune::Concept::Grid Grid_,
25 Dune::Concept::GridView GridView_ =
typename Grid_::LeafGridView,
26 class RangeQuatinty_ = double,
27 class TimeQuantity_ =
double>
38 std::shared_ptr<const Grid>
grid;
45 RangeQuatinty(
typename GridView::template Codim<0>::Geometry::GlobalCoordinate),
57 [[nodiscard]]
virtual std::unique_ptr<State>
make_state(
const std::shared_ptr<const Grid>&,
58 const ParameterTree&)
const = 0;
60 virtual void interpolate(
State&,
const std::unordered_map<std::string, GridFunction>&)
const = 0;
62 [[nodiscard]]
virtual std::unordered_map<std::string, GridFunction>
make_initial(
64 const ParameterTree&)
const = 0;
68 std::string_view name)
const
75 std::string_view)
const = 0;
83 std::string_view)
const
88 virtual void write_vtk(
const State&,
const std::filesystem::path&,
bool)
const
90 throw format_exception(NotImplemented{},
"Model write has not been implemented");
95 const ParameterTree&)
const = 0;
97 virtual std::map<std::string, double>
reduce(
const State&,
const ParameterTree&,
const ParameterTree& = {})
const
99 throw format_exception(NotImplemented{},
"Model reduce has not been implemented");
Definition: axis_names.hh:7
auto format_exception(Exception &&e, fmt::format_string< Args... > format, Args &&... args)
Definition: exceptions.hh:23
std::any coefficients
Definition: model.hh:40
TimeQuantity time
Definition: model.hh:39
std::any basis
Definition: model.hh:41
std::shared_ptr< const Grid > grid
Definition: model.hh:38
virtual std::unique_ptr< PDELab::OneStep< State > > make_step_operator(const State &, const ParameterTree &) const =0
GridView_ GridView
Definition: model.hh:32
Model(const Model &)=delete
virtual std::nullptr_t make_membrane_function(const State &, std::string_view) const
Definition: model.hh:77
Model & operator=(const Model &)=delete
virtual GridFunction make_compartment_function(const std::shared_ptr< const State > &, std::string_view) const =0
TimeQuantity_ TimeQuantity
Definition: model.hh:33
virtual std::nullptr_t make_membrane_function(const std::shared_ptr< const State > &, std::string_view) const
Definition: model.hh:82
Grid_ Grid
Definition: model.hh:31
Model & operator=(Model &&)=delete
virtual void interpolate(State &, const std::unordered_map< std::string, GridFunction > &) const =0
virtual std::unordered_map< std::string, GridFunction > make_initial(const Grid &, const ParameterTree &) const =0
virtual std::unique_ptr< State > make_state(const std::shared_ptr< const Grid > &, const ParameterTree &) const =0
virtual std::map< std::string, double > reduce(const State &, const ParameterTree &, const ParameterTree &={}) const
Definition: model.hh:97
GridFunction make_compartment_function(const State &state, std::string_view name) const
Definition: model.hh:67
virtual void write_vtk(const State &, const std::filesystem::path &, bool) const
Definition: model.hh:88
Dune::Functions::GridViewFunction< RangeQuatinty(typename GridView::template Codim< 0 >::Geometry::GlobalCoordinate), GridView > GridFunction
Definition: model.hh:46
RangeQuatinty_ RangeQuatinty
Definition: model.hh:34