Dune::Copasi
Loading...
Searching...
No Matches
model_multi_compartment.impl.hh
Go to the documentation of this file.
1#ifndef DUNE_COPASI_MODEL_DIFFUSION_REACTION_MULTI_COMPARTMENT_IMPL_HH
2#define DUNE_COPASI_MODEL_DIFFUSION_REACTION_MULTI_COMPARTMENT_IMPL_HH
3
11
12#include <dune/pdelab/basis/backend/istl.hh>
13#include <dune/pdelab/basis/basis.hh>
14#include <dune/pdelab/basis/discrete_global_function.hh>
15#include <dune/pdelab/common/trace.hh>
16#include <dune/pdelab/common/partition/identity.hh>
17#include <dune/pdelab/concepts/multiindex.hh>
18
19#include <dune/grid/io/file/vtk.hh>
20
21#include <spdlog/spdlog.h>
22
23#include <fmt/core.h>
24
25#ifdef DUNE_COPASI_PRECOMPILED_MODE
26#warning "Including this file in pre-compiled mode may defeat the purpose of pre-compilation"
27#endif
28
30
31template<class Traits>
32void
34 State& state,
35 const std::unordered_map<std::string, GridFunction>& initial) const
36{
38 PDELab::Basis<PDELab::EntitySetPartitioner::Identity<MultiCompartmentEntitySet>, MultiCompartmentPreBasis, TypeTree::HybridTreePath<>>;
39
40 using CoefficientsBackend = PDELab::ISTLUniformBackend<ScalarQuantity>;
41 using Coefficients = typename MultiCompartmentBasis::template Container<CoefficientsBackend>;
43 PDELab::Basis<PDELab::EntitySetPartitioner::Identity<MultiCompartmentEntitySet>, MultiCompartmentPreBasis, TypeTree::HybridTreePath<>>;
44
45 using CoefficientsBackend = PDELab::ISTLUniformBackend<ScalarQuantity>;
46 using Coefficients = typename MultiCompartmentBasis::template Container<CoefficientsBackend>;
47 const auto& basis = any_cast<const MultiCompartmentBasis&>(state.basis);
48 auto& coefficients = any_cast<Coefficients&>(state.coefficients);
49
50 Dune::Copasi::interpolate(basis, coefficients, initial);
51}
52
53template<class Traits>
54auto
56 const Grid& grid,
57 const ParameterTree& config,
58 std::shared_ptr<const FunctorFactory<Grid::dimensionworld>> functor_factory) -> MultiCompartmentPreBasis
59{
60 TRACE_EVENT("dune", "Basis::SetUp");
61 spdlog::info("Setup basis functions");
62
63 const auto& compartments_config = config.sub("compartments", true);
64
65 std::vector<CompartmentPreBasis> compartment_pre_basis_vec;
66
67 std::map<std::string, std::vector<std::string>> compartment2componets;
68 const auto& scalar_fields_config = config.sub("scalar_field");
69 for (const auto& component : scalar_fields_config.getSubKeys())
70 compartment2componets[config[fmt::format("scalar_field.{}.compartment", component)]].push_back(component);
71
72 for (const auto& compartment : compartments_config.getSubKeys()) {
73 using SubDomainIndex = typename Grid::SubDomainIndex;
76 compartments_config.sub(compartment, true).template get<SubDomainIndex>("id");
77
78 CompartmentEntitySet sub_grid_view = grid.subDomain(domain_id).leafGridView();
79
80 CompartmentPreBasis compartment_pre_basis =
84 }
85
86 spdlog::info("Setup of multi-compartment grid function space");
87
88 MultiCompartmentPreBasis multi_compartment_pre_basis =
89 composite(MultiCompartmentMergingStrategy{}, compartment_pre_basis_vec);
90
91 multi_compartment_pre_basis.name("compartments");
93}
94
95template<class Traits>
96void
97ModelMultiCompartment<Traits>::setup_coefficient_vector(State& state)
98{
100 PDELab::Basis<PDELab::EntitySetPartitioner::Identity<MultiCompartmentEntitySet>, MultiCompartmentPreBasis, TypeTree::HybridTreePath<>>;
101 using CoefficientsBackend = PDELab::ISTLUniformBackend<ScalarQuantity>;
102 using Coefficients = typename MultiCompartmentBasis::template Container<CoefficientsBackend>;
103
104 spdlog::info("Setup coefficient vector");
105 const auto& basis = any_cast<const MultiCompartmentBasis&>(state.basis);
106 state.coefficients = Coefficients{ basis.makeContainer(CoefficientsBackend{}) };
107}
108
109template<class Traits>
110auto
111ModelMultiCompartment<Traits>::make_state(const std::shared_ptr<const Grid>& grid,
112 const ParameterTree& config) const
113 -> std::unique_ptr<State>
114{
115 auto state_ptr = std::make_unique<State>();
116 state_ptr->basis = makeBasis(grid->leafGridView(),
117 make_multi_compartment_pre_basis(*grid, config, _functor_factory));
118 setup_coefficient_vector(*state_ptr);
119 state_ptr->grid = grid;
120 state_ptr->time = TimeQuantity{ 0. };
121 return state_ptr;
122}
123
124template<class Traits>
125auto
127 const std::shared_ptr<const State>& state,
128 std::string_view name) const -> GridFunction
129{
131 PDELab::Basis<PDELab::EntitySetPartitioner::Identity<MultiCompartmentEntitySet>, MultiCompartmentPreBasis, TypeTree::HybridTreePath<>>;
132 using CoefficientsBackend = PDELab::ISTLUniformBackend<ScalarQuantity>;
133 using Coefficients = typename MultiCompartmentBasis::template Container<CoefficientsBackend>;
134
135 const auto& basis = any_cast<const MultiCompartmentBasis&>(state->basis);
136 const auto& coeff = any_cast<const Coefficients&>(state->coefficients);
137 std::shared_ptr<const Coefficients> coeff_ptr(state, &coeff);
138
139 MultiCompartmentBasis multi_compartment_basis = basis.subSpace(TypeTree::treePath());
140 for (std::size_t compartment = 0; compartment != multi_compartment_basis.degree();
141 ++compartment) {
142 auto compartment_basis =
143 multi_compartment_basis.subSpace(TypeTree::treePath(compartment));
144 for (std::size_t species = 0; species != compartment_basis.degree(); ++species) {
145 auto species_basis = compartment_basis.subSpace(TypeTree::treePath(species));
146 if (species_basis.name() == name)
148 }
149 }
150 throw format_exception(RangeError{}, "State doesn't contain any function with name: {}", name);
151}
152
153template<class Traits>
154auto
156 const ParameterTree& config) const
157 -> std::unordered_map<std::string, GridFunction>
158{
159 return Dune::Copasi::make_initial<GridFunction>(grid, config, *_functor_factory);
160}
161
162template<class Traits>
163auto
165 -> std::map<std::string, double>
166{
168 PDELab::Basis<PDELab::EntitySetPartitioner::Identity<MultiCompartmentEntitySet>, MultiCompartmentPreBasis, TypeTree::HybridTreePath<>>;
169 using CoefficientsBackend = PDELab::ISTLUniformBackend<ScalarQuantity>;
170 using Coefficients = typename MultiCompartmentBasis::template Container<CoefficientsBackend>;
171
172 // TODO(sospinar): set/get basis with concurrent entity set partition
173 const auto& basis = any_cast<const MultiCompartmentBasis&>(state.basis);
174 const auto& coeff = any_cast<const Coefficients&>(state.coefficients);
175
176 return Dune::Copasi::DiffusionReaction::reduce(PDELab::Execution::seq, basis, coeff, state.time, config, _functor_factory);
177}
178
179template<class Traits>
180auto
182 const State& state,
183 const ParameterTree& config) const -> std::unique_ptr<PDELab::OneStep<State>>
184{
186 PDELab::Basis<PDELab::EntitySetPartitioner::Identity<MultiCompartmentEntitySet>, MultiCompartmentPreBasis, TypeTree::HybridTreePath<>>;
187 using CoefficientsBackend = PDELab::ISTLUniformBackend<ScalarQuantity>;
188 using Coefficients = typename MultiCompartmentBasis::template Container<CoefficientsBackend>;
189 using ResidualBackend = PDELab::ISTLUniformBackend<ResidualQuantity>;
190 using Residual = typename MultiCompartmentBasis::template Container<ResidualBackend>;
191 using LocalBasisTraits = typename ScalarFiniteElementMap::Traits::FiniteElement::Traits::LocalBasisType::Traits;
192
193 const auto& basis = any_cast<const MultiCompartmentBasis&>(state.basis);
194 std::shared_ptr one_step = DiffusionReaction::make_step_operator<LocalBasisTraits, Coefficients, Residual, ResidualQuantity, TimeQuantity>(config, basis, 2, _functor_factory, _cell_data);
195
196 // type erase the original runge kutta operator
197 auto type_erased_one_step = std::make_unique<PDELab::OperatorAdapter<State, State>>(
198 [one_step](PDELab::Operator<State, State>& self, State& domain, State& range) mutable {
199 auto log_guard = ostream2spdlog();
200 // copy contents of this operator into runge-kutta operator
201 static_cast<PDELab::PropertyTree&>(*one_step) =
202 static_cast<const PDELab::PropertyTree&>(self);
205 return one_step->apply(domain_coeff, range_coeff);
206 });
207
208 // assign properties of the original one step to the type erased one
209 static_cast<PDELab::PropertyTree&>(*type_erased_one_step) =
210 static_cast<const PDELab::PropertyTree&>(*one_step);
211
212 auto residual_ptr = std::make_shared<Residual>(basis.makeContainer(ResidualBackend{}));
213 type_erased_one_step->get("initial_residual") = residual_ptr;
214 type_erased_one_step->get("time") = state.time;
216}
217
218template<class Traits>
219void
221 const std::filesystem::path& path,
222 bool append) const
223{
225 PDELab::Basis<PDELab::EntitySetPartitioner::Identity<MultiCompartmentEntitySet>, MultiCompartmentPreBasis, TypeTree::HybridTreePath<>>;
226 using CoefficientsBackend = PDELab::ISTLUniformBackend<ScalarQuantity>;
227 using Coefficients = typename MultiCompartmentBasis::template Container<CoefficientsBackend>;
228
229 const auto& basis = any_cast<const MultiCompartmentBasis&>(state.basis);
230 const auto& coeff = any_cast<const Coefficients&>(state.coefficients);
231 // warning: we use Dune::stackobject_to_shared_ptr to avoid copying the coefficients vector,
232 // be we need to check that no one else took ownership of this pointer when we leave this function!
233 std::shared_ptr<const Coefficients> coeff_ptr = Dune::stackobject_to_shared_ptr(coeff);
234
235 // create directory if necessary
236 auto path_entry = std::filesystem::directory_entry{ path };
237 if (not path_entry.exists()) {
238 spdlog::info("Creating output directory '{}'", path_entry.path().string());
239 std::error_code ec{ 0, std::generic_category() };
240 std::filesystem::create_directories(path_entry.path(), ec);
241 if (ec)
243 "Category: {}\n"
244 "Value: {}\n"
245 "Message: {}",
246 ec.category().name(),
247 ec.value(),
248 ec.message());
249 }
250
251 spdlog::info("Writing solution for {:.2f}s time stamp", state.time);
252 MultiCompartmentBasis multi_compartment_basis = basis.subSpace(TypeTree::treePath());
253
254 // Recover old timestesps in case something was written before
255 auto& timesteps = _writer_timesteps[path.string()];
256 std::vector<double> tmp_timesteps;
257
258 for (std::size_t compartment = 0; compartment != multi_compartment_basis.degree();
259 ++compartment) {
260 auto compartment_basis =
261 multi_compartment_basis.subSpace(TypeTree::treePath(compartment));
262 std::string name = fmt::format("{}-{}", path.filename().string(), compartment_basis.name());
263 if (not append) {
264 timesteps.clear();
265 spdlog::info("Creating a time sequence file: '{}.pvd'", name);
266 } else {
267 spdlog::info("Overriding time sequence file: '{}.pvd'", name);
268 }
269
270 // setup writer again with old timesteps if necessary
271 auto writer = std::make_shared<VTKWriter<CompartmentEntitySet>>(compartment_basis.entitySet(),
272 Dune::VTK::conforming);
273 auto sequential_writer =
274 VTKSequenceWriter<CompartmentEntitySet>{ writer, name, path.string(), path.string() };
275 sequential_writer.setTimeSteps(timesteps);
276
277 for (std::size_t species = 0; species != compartment_basis.degree(); ++species) {
278 auto species_basis =
279 compartment_basis.subSpace(TypeTree::treePath(std::size_t{ species }));
280 sequential_writer.vtkWriter()->addVertexData(
282 VTK::FieldInfo{ species_basis.name(), VTK::FieldInfo::Type::scalar, 1 });
283 }
284 spdlog::info("Writing vtu file: '{0}/{0}-{1:0>5}.vtu'", name, timesteps.size());
285 TRACE_EVENT("dune", "WriteVTK");
286 sequential_writer.write(state.time, Dune::VTK::base64);
287 sequential_writer.vtkWriter()->clear();
288
289 tmp_timesteps = sequential_writer.getTimeSteps();
290 }
291
293
294 if (coeff_ptr.use_count() != 1)
295 throw format_exception(
297 "Fake shared pointer from coefficient vector may have been leaked outsie of this function!");
298}
299
300} // namespace Dune::Copasi::DiffusionReaction
301
302#endif // DUNE_COPASI_MODEL_DIFFUSION_REACTION_MULTI_COMPARTMENT_IMPL_HH
Definition model_multi_compartment.hh:28
void write_vtk(const State &, const std::filesystem::path &, bool) const override
Definition model_multi_compartment.impl.hh:220
GridFunction make_compartment_function(const std::shared_ptr< const State > &, std::string_view) const override
Definition model_multi_compartment.impl.hh:126
typename Base::State State
Definition model_multi_compartment.hh:35
void interpolate(State &, const std::unordered_map< std::string, GridFunction > &) const override
Definition model_multi_compartment.impl.hh:33
PDELab::PreBasisVector< MultiCompartmentMergingStrategy, CompartmentPreBasis > MultiCompartmentPreBasis
Definition model_multi_compartment.hh:51
std::unordered_map< std::string, GridFunction > make_initial(const Grid &, const ParameterTree &) const override
Definition model_multi_compartment.impl.hh:155
typename Traits::Grid Grid
Definition model_multi_compartment.hh:36
typename Base::GridFunction GridFunction
Definition model_multi_compartment.hh:56
typename Traits::TimeQuantity TimeQuantity
Definition model_multi_compartment.hh:37
std::unique_ptr< State > make_state(const std::shared_ptr< const Grid > &, const ParameterTree &) const override
Definition model_multi_compartment.impl.hh:111
std::unique_ptr< PDELab::OneStep< State > > make_step_operator(const State &, const ParameterTree &) const override
Definition model_multi_compartment.impl.hh:181
std::map< std::string, double > reduce(const State &, const ParameterTree &) const override
Definition model_multi_compartment.impl.hh:164
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
Definition functor_factory.hh:24
Definition factory.hh:28
auto format_exception(Exception &&e, fmt::format_string< Args... > format, Args &&... args)
Definition exceptions.hh:23
constexpr bool is_bitflags_v
Alias for Bitflag indicator.
Definition bit_flags.hh:24
auto ostream2spdlog()
Standard output redirection to spdlog.
Definition ostream_to_spdlog.hh:20