Dune::Copasi
Loading...
Searching...
No Matches
model_single_compartment.hh
Go to the documentation of this file.
1#ifndef DUNE_COPASI_MODEL_DIFFUSION_REACTION_SINGLE_COMAPARTMENT_HH
2#define DUNE_COPASI_MODEL_DIFFUSION_REACTION_SINGLE_COMAPARTMENT_HH
3
4// file: diffusion reaction for single compartment models
5
11
12#include <dune/pdelab/basis/prebasis/composite.hh>
13#include <dune/pdelab/basis/prebasis/leaf.hh>
14
15#include <dune/grid/concepts/grid.hh>
16#include <dune/grid/concepts/gridview.hh>
17
19
20template<class Traits>
22 : public Model<typename Traits::Grid,
23 typename Traits::CompartmentEntitySet,
24 typename Traits::RangeQuatinty,
25 typename Traits::TimeQuantity>
26{
27 using Base = Model<typename Traits::Grid,
28 typename Traits::CompartmentEntitySet,
29 typename Traits::RangeQuatinty,
30 typename Traits::TimeQuantity>;
31
32public:
33 using State = typename Base::State;
34 using Grid = typename Traits::Grid;
35 using TimeQuantity = typename Traits::TimeQuantity;
36 using ScalarQuantity = typename Traits::RangeQuatinty;
37 using CompartmentEntitySet = typename Traits::CompartmentEntitySet;
38 using ScalarFiniteElementMap = typename Traits::ScalarFiniteElementMap;
39 using ScalarMergingStrategy = typename Traits::ScalarMergingStrategy;
41 PDELab::PreBasis<ScalarMergingStrategy, ScalarFiniteElementMap, Constraints<CompartmentEntitySet>>;
42 using CompartmentMergingStrategy = typename Traits::CompartmentMergingStrategy;
43 using CompartmentPreBasis = PDELab::PreBasisVector<CompartmentMergingStrategy, ScalarPreBasis>;
45
47
51 : _functor_factory{ std::move(functor_factory) }
52 , _cell_data{ std::move(cell_data) }
53 {
54 assert(_functor_factory);
55 }
56
57 std::unique_ptr<State> make_state(const std::shared_ptr<const Grid>&,
58 const ParameterTree&) const override;
59
60 void interpolate(State&, const std::unordered_map<std::string, GridFunction>&) const override;
61
62 std::unordered_map<std::string, GridFunction> make_initial(const Grid&,
63 const ParameterTree&) const override;
64
65 GridFunction make_compartment_function(const std::shared_ptr<const State>&,
66 std::string_view) const override;
67
68 std::unique_ptr<PDELab::OneStep<State>> make_step_operator(const State&,
69 const ParameterTree&) const override;
70
71 void write_vtk(const State&, const std::filesystem::path&, bool = true) const override;
72
73 std::map<std::string, double> reduce(const State&, const ParameterTree&) const override;
74
76 std::string_view,
77 const std::vector<std::string>&,
78 const ParameterTree& = {},
79 std::shared_ptr<const FunctorFactory<Grid::dimensionworld>> = nullptr);
80
81private:
82 static ScalarPreBasis make_scalar_field_pre_basis(std::shared_ptr<BoundaryEntityMapper<CompartmentEntitySet>>, const CompartmentEntitySet&, std::string_view, const ParameterTree&, std::shared_ptr<const FunctorFactory<Grid::dimensionworld>>);
83 static void setup_basis(State&, const Grid&, const ParameterTree&, std::shared_ptr<const FunctorFactory<Grid::dimensionworld>>);
84 static void setup_coefficient_vector(State&);
85 static CompartmentEntitySet get_entity_set(const Grid&, std::size_t);
86
87 mutable std::unordered_map<std::string, std::vector<double>> _writer_timesteps;
88 std::shared_ptr<const FunctorFactory<Grid::dimensionworld>> _functor_factory;
89 std::shared_ptr<const CellData<typename Grid::LeafGridView, ScalarQuantity>> _cell_data;
90};
91
92} // namespace Dune::Copasi::DiffusionReaction
93
94#ifndef DUNE_COPASI_PRECOMPILED_MODE
96#endif
97
98#endif // DUNE_COPASI_MODEL_DIFFUSION_REACTION_SINGLE_COMAPARTMENT_HH
Mapper of boundary entities.
Definition boundary_entity_mapper.hh:26
Container for cell data of a grid view.
Definition cell_data.hh:25
Definition model_single_compartment.hh:26
std::unique_ptr< State > make_state(const std::shared_ptr< const Grid > &, const ParameterTree &) const override
Definition model_single_compartment.impl.hh:154
PDELab::PreBasis< ScalarMergingStrategy, ScalarFiniteElementMap, Constraints< CompartmentEntitySet > > ScalarPreBasis
Definition model_single_compartment.hh:40
ScalarQuantity ResidualQuantity
Definition model_single_compartment.hh:44
typename Traits::CompartmentEntitySet CompartmentEntitySet
Definition model_single_compartment.hh:37
typename Traits::CompartmentMergingStrategy CompartmentMergingStrategy
Definition model_single_compartment.hh:42
typename Traits::ScalarFiniteElementMap ScalarFiniteElementMap
Definition model_single_compartment.hh:38
static CompartmentPreBasis make_compartment_pre_basis(const CompartmentEntitySet &, std::string_view, const std::vector< std::string > &, const ParameterTree &={}, std::shared_ptr< const FunctorFactory< Grid::dimensionworld > >=nullptr)
Definition model_single_compartment.impl.hh:82
typename Traits::Grid Grid
Definition model_single_compartment.hh:34
typename Base::State State
Definition model_single_compartment.hh:33
GridFunction make_compartment_function(const std::shared_ptr< const State > &, std::string_view) const override
Definition model_single_compartment.impl.hh:168
std::map< std::string, double > reduce(const State &, const ParameterTree &) const override
Definition model_single_compartment.impl.hh:236
PDELab::PreBasisVector< CompartmentMergingStrategy, ScalarPreBasis > CompartmentPreBasis
Definition model_single_compartment.hh:43
typename Base::GridFunction GridFunction
Definition model_single_compartment.hh:46
typename Traits::ScalarMergingStrategy ScalarMergingStrategy
Definition model_single_compartment.hh:39
void interpolate(State &, const std::unordered_map< std::string, GridFunction > &) const override
Definition model_single_compartment.impl.hh:50
void write_vtk(const State &, const std::filesystem::path &, bool=true) const override
Definition model_single_compartment.impl.hh:252
typename Traits::TimeQuantity TimeQuantity
Definition model_single_compartment.hh:35
std::unordered_map< std::string, GridFunction > make_initial(const Grid &, const ParameterTree &) const override
Definition model_single_compartment.impl.hh:189
typename Traits::RangeQuatinty ScalarQuantity
Definition model_single_compartment.hh:36
std::unique_ptr< PDELab::OneStep< State > > make_step_operator(const State &, const ParameterTree &) const override
Definition model_single_compartment.impl.hh:197
ModelSingleCompartment(std::shared_ptr< const FunctorFactory< Grid::dimensionworld > > functor_factory, std::shared_ptr< const CellData< typename Grid::LeafGridView, ScalarQuantity > > cell_data=nullptr)
Definition model_single_compartment.hh:48
Definition functor_factory.hh:24
Definition factory.hh:28
constexpr bool is_bitflags_v
Alias for Bitflag indicator.
Definition bit_flags.hh:24
Definition model.hh:37
Definition model.hh:29
Dune::Functions::GridViewFunction< RangeQuatinty(typename GridView::template Codim< 0 >::Geometry::GlobalCoordinate), GridView > GridFunction
Definition model.hh:44