Dune::Copasi 2.1.0-git79717215+dune.gitlab.629933
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages Concepts
reduce.hh
Go to the documentation of this file.
1#ifndef DUNE_COPASI_MODEL_DIFFUSION_REACTION_REDUCE_HH
2#define DUNE_COPASI_MODEL_DIFFUSION_REACTION_REDUCE_HH
3
11
12#include <dune/pdelab/concepts/basis.hh>
13#include <dune/pdelab/common/local_container.hh>
14#include <dune/pdelab/common/trace.hh>
15#include <dune/pdelab/common/for_each_entity.hh>
16
17#include <dune/grid/common/partitionset.hh>
18
19#include <dune/common/float_cmp.hh>
20
21#include <function2/function2.hpp>
22#include <spdlog/spdlog.h>
23
24#include <fmt/core.h>
25
26#if HAVE_TBB
27#include <tbb/enumerable_thread_specific.h>
28#endif // HAVE_TBB
29
30#include <string>
31#include <map>
32
34
35class ReductionError : public Dune::Exception {};
36
37// Evaluates a evaluation/reduction/transformation algorigthm over the grid.
38// Starting with val = 0, the this function evaluates `val = reduction(val, evaluation(), weight)` for each quadrature point of the grid. Finally, the result gets transformed by `value = transformation(value)`
39// For each sub-tree in the parameter tree, the evaluation and reduce expressions are extracted to form the evaluation/reduce operation for its key.
40template<class ExecutionPolicy, PDELab::Concept::Basis Basis, class Coefficients>
42 PDELab::Execution::is_execution_policy_v<ExecutionPolicy> )
43static std::map<std::string, double>
44reduce(const ExecutionPolicy exec,
45 const Basis& basis,
46 const Coefficients& coefficients,
47 auto time,
48 const ParameterTree& config,
49 const ParameterTree& domain_cfg = {},
50 std::shared_ptr<const FunctorFactory<Basis::EntitySet::GridView::dimension>>
51 functor_factory = nullptr)
52{
53 constexpr std::size_t dim = Basis::EntitySet::GridView::dimension;
54 TRACE_EVENT("dune", "Reduce");
55
56 if (basis.entitySet().size(0) == 0) {
57 return {};
58 }
59
60 std::shared_ptr<const ParserContext> parser_context;
61 auto functor_factory_parser = std::dynamic_pointer_cast<const FunctorFactoryParser<dim>>(functor_factory);
62 if (functor_factory_parser) {
63 parser_context = functor_factory_parser->parser_context();
64 }
65 if (not parser_context) {
66 parser_context = std::make_shared<const ParserContext>();
67 }
68
69 auto parser_type = string2parser.at(config.get("parser_type", default_parser_str));
70
71 auto first_compartment_finite_elment = [](const Concept::CompartmentLocalBasisNode auto& lnode) noexcept {
72 for (std::size_t i = 0; i != lnode.degree(); ++i)
73 if (lnode.child(i).size() != 0)
74 return lnode.child(i).finiteElement();
75 std::terminate();
76 };
77
78 auto first_finite_element = overload(first_compartment_finite_elment,
79 [first_compartment_finite_elment](const Concept::MultiCompartmentLocalBasisNode auto& lnode) noexcept {
80 for (std::size_t i = 0; i != lnode.degree(); ++i)
81 if (lnode.child(i).size() != 0)
82 return first_compartment_finite_elment(lnode.child(i));
83 std::terminate();
84 });
85
86 using FEM = std::decay_t<decltype(first_finite_element(basis.localView().tree()))>;
87 using Geometry = typename Basis::EntitySet::template Codim<0>::Entity::Geometry;
88
89 struct ThreadLocalData
90 {
91 typename Basis::LocalView lbasis;
92 PDELab::LocalContainerBuffer<Basis, const Coefficients> lcoeff;
93 std::unique_ptr<LocalEquations<dim>> leqs;
94 std::optional<typename Geometry::JacobianInverse> geojacinv_opt = std::nullopt;
95 std::vector<double> values = {};
96 std::vector<fu2::unique_function<double() const>> evaluations = {};
97 std::vector<fu2::unique_function<double(double, double) const>> reductions = {};
98 };
99
100 auto init_thread_data = [&]() -> ThreadLocalData {
101 auto lbasis = basis.localView();
102 auto leqs_ptr = LocalEquations<dim>::make(lbasis, {}, domain_cfg, functor_factory);
103 leqs_ptr->time = time;
104 const auto& leqs = *leqs_ptr;
105 ThreadLocalData data{ std::move(lbasis), { basis, &coefficients }, std::move(leqs_ptr) };
106
107 auto sz = config.getSubKeys().size();
108 data.values.reserve(sz);
109 data.evaluations.reserve(sz);
110 data.reductions.reserve(sz);
111 for (const auto& key : config.getSubKeys()) {
112 const auto& sub_config = config.sub(key, true);
113 data.values.emplace_back(sub_config.get("initial.value", 0.));
114 auto& evaluation = data.evaluations.emplace_back();
115 if (sub_config.hasKey("evaluation.expression"))
116 evaluation =
117 functor_factory->make_scalar(key + ".evaluation", sub_config.sub("evaluation"), leqs);
118 if (not evaluation)
119 evaluation = [] { return 0.; };
120 auto& reduction = data.reductions.emplace_back(std::plus<>{});
121 if (sub_config.hasKey("reduction.expression")) {
122 auto [args, expr] =
123 parser_context->parse_function_expression(sub_config.sub("reduction")["expression"]);
124 if (args.size() == 3)
125 spdlog::warn(
126 "Reduction expression for '{0}' uses 3 arguments whereas 2 are now required.\nThis "
127 "functionality has changed: the third argument became the contextual keyword "
128 "'integration_factor' that may be used to express '{0}.evaluation.expression'.",
129 key);
130 if (args.size() != 2)
131 throw format_exception(IOError{}, "Reduction arguments must be exactly 2");
132 auto sub_parser_type = string2parser.at(
133 sub_config.get("reduction.parser_type", std::string{ parser2string.at(parser_type) }));
134 auto str_args = std::array{ std::string{ args[0] }, std::string{ args[1] } };
135 reduction = parser_context->make_function(sub_parser_type, str_args, expr);
136 }
137 }
138 return data;
139 };
140
141#if HAVE_TBB
142 tbb::enumerable_thread_specific<ThreadLocalData> thread_data{ std::move(init_thread_data) };
143#else
144 struct ThreadData : public std::array<ThreadLocalData,1> {
145 ThreadLocalData& local() {
146 return this->front();
147 }
148 };
149 ThreadData thread_data{init_thread_data()};
150#endif
151
152 if (config.getSubKeys().empty()) {
153 return {};
154 }
155
156 spdlog::info("Evaluating reduce operators");
157 PDELab::forEachEntity(exec, basis.entitySetPartition(), [&thread_data, time](const auto& cell) {
158
159 thread_local LocalBasisCache<typename FEM::Traits::LocalBasisType::Traits> fe_cache = {};
160
161 auto& data = thread_data.local();
162 data.lbasis.bind(cell);
163 data.lcoeff.load(data.lbasis, std::false_type{});
164
165 const auto& geo = move_geometry(time, data.lbasis, data.lcoeff, data.leqs, fe_cache, cell.geometry(), std::identity{});
166
167 std::size_t const order = 4;
168 thread_local const auto quad_rule =
169 QuadratureRules<typename Geometry::ctype, Geometry::coorddimension>::rule(geo.type(), order);
170
171 for (std::size_t q = 0; q != quad_rule.size(); ++q) {
172 const auto [position, weight] = quad_rule[q];
173 if (not geo.affine() or not data.geojacinv_opt)
174 data.geojacinv_opt.emplace(geo.jacobianInverse(position));
175 const auto& geojacinv = *(data.geojacinv_opt);
176 data.leqs->position = geo.global(position);
177 data.leqs->integration_factor = weight * geo.integrationElement(position);
178
179 // evaluate values at quadrature point
180 forEachLeafNode(data.lbasis.tree(), [&](const auto& node) {
181 if (node.size() == 0)
182 return;
183 auto& value = data.leqs->get_value(node);
184 auto& gradient = data.leqs->get_gradient(node);
185 value = 0.;
186 gradient = 0.;
187 fe_cache.bind(node.finiteElement(), quad_rule);
188 const auto& phi = fe_cache.evaluateFunction(q);
189 const auto& jacphi = fe_cache.evaluateJacobian(q);
190 for (std::size_t dof = 0; dof != node.size(); ++dof) {
191 value += data.lcoeff(node, dof) * phi[dof];
192 gradient += data.lcoeff(node, dof) * (jacphi[dof] * geojacinv)[0];
193 }
194 });
195
196 // evaluate all functions at the current quadrature point
197 for (std::size_t i = 0; i != data.values.size(); ++i)
198 data.values[i] = data.reductions[i](data.evaluations[i](), data.values[i]);
199 }
200
201 data.lbasis.unbind();
202 data.leqs->clear();
203 });
204
205 if (basis.entitySet().grid().comm().size() > 1){
206 throw format_exception(ParallelError{}, "Transoform and Reduce is not implemented in parallel");
207 }
208
209 // gather values from every thread_data
210 std::vector<double> values;
211 for (const auto& key : config.getSubKeys())
212 values.emplace_back(config.sub(key, true).get("initial.value", 0.));
213 for (std::size_t i = 0; i != values.size(); ++i)
214 for (const auto& data : thread_data)
215 values[i] = data.reductions[i](data.values[i], values[i]);
216
217 // report results
218 std::size_t max_key_chars = 0;
219 std::size_t max_val_chars = 5;
220 for (const auto& key : config.getSubKeys())
221 max_key_chars = std::max(max_key_chars, key.size());
222
223 std::map<std::string, double> key_value;
224
225 spdlog::info(" ┌{0:─^{1}}┐", "", max_key_chars + max_val_chars + 13);
226 std::string error_msg;
227 auto values_it = values.begin();
228 for (const auto& key : config.getSubKeys()) {
229 auto& value = key_value[key] = *(values_it++);
230
231 if (config.sub(key).hasKey("transformation.expression")) {
232 auto [args, expr] = parser_context->parse_function_expression(config.sub(key)["transformation.expression"]);
233 if (args.size() != 1)
234 throw format_exception(IOError{}, "Warning function must have exactly 1 argument");
235 auto sub_parser_type = string2parser.at(config.sub(key).get("transformation.parser_type", std::string{parser2string.at(parser_type)}));
236 auto transformation = parser_context->make_function(sub_parser_type, std::array{std::string{args[0]}}, expr);
237 value = transformation(value);
238 }
239
240 bool processed = false;
241
242 if (config.sub(key).hasKey("error.expression")) {
243 auto [args, expr] = parser_context->parse_function_expression(config.sub(key)["error.expression"]);
244 if (args.size() != 1)
245 throw format_exception(IOError{}, "Error function must have exactly 1 argument");
246 auto sub_parser_type = string2parser.at(config.sub(key).get("error.parser_type", std::string{parser2string.at(parser_type)}));
247 auto error = parser_context->make_function(sub_parser_type, std::array{std::string{args[0]}}, expr);
248 if (FloatCmp::ne(error(value), 0.)) {
249 processed = true;
250 spdlog::error(" | {0:>{2}} := {1: .{3}e} |", DUNE_COPASI_FMT_STYLED_BOLD(key), value, max_key_chars, max_val_chars);
251 if (not error_msg.empty())
252 error_msg += '\n';
253 error_msg += fmt::format("Reduction on the token '{}' raised an error because the "
254 "expression '{}' with evaluates to false with '{} := {}'",
255 key,
256 expr,
257 args[0],
258 value);
259 }
260 }
261
262 if (not processed and config.sub(key).hasKey("warn.expression")) {
263 auto [args, expr] = parser_context->parse_function_expression(config.sub(key)["warn.expression"]);
264 if (args.size() != 1)
265 throw format_exception(IOError{}, "Warning function must have exactly 1 argument");
266 auto sub_parser_type = string2parser.at(config.sub(key).get("warn.parser_type", std::string{parser2string.at(parser_type)}));
267 auto warn = parser_context->make_function(sub_parser_type, std::array{std::string{args[0]}}, expr);
268 if (FloatCmp::ne(warn(value), 0.)) {
269 processed = true;
270 spdlog::warn("| {0:>{2}} := {1: .{3}e} |", DUNE_COPASI_FMT_STYLED_BOLD(key), value, max_key_chars, max_val_chars);
271 }
272 }
273
274 if (not processed and not config.sub(key).get("quiet", false)) {
275 spdlog::info(" | {0:>{2}} := {1: .{3}e} |", DUNE_COPASI_FMT_STYLED_BOLD(key), value, max_key_chars, max_val_chars);
276 }
277 }
278 spdlog::info(" └{0:─^{1}}┘", "", max_key_chars+max_val_chars+13);
279
280 if (not error_msg.empty()) {
281 throw format_exception(ReductionError{}, "{}", error_msg);
282 }
283
284 return key_value;
285}
286
287} // namespace Dune::Copasi::DiffusionReaction
288
289#endif // DUNE_COPASI_MODEL_DIFFUSION_REACTION_REDUCE_HH
static std::unique_ptr< LocalEquations > make(const LocalBasis &lbasis, const ParameterTree &eqs_cfg={}, const ParameterTree &domain_cfg={}, std::shared_ptr< const FunctorFactory< dim > > functor_factory=nullptr, std::shared_ptr< const CellData< CellDataGridView, CellDataType > > grid_cell_data=nullptr, BitFlags< FactoryFalgs > opts=BitFlags< FactoryFalgs >::no_flags())
Definition: local_equations.hh:264
Definition: compartment_tree.hh:99
#define DUNE_COPASI_FMT_STYLED_BOLD(key)
Definition: fmt_style.hh:13
Definition: factory.hh:28
auto format_exception(Exception &&e, fmt::format_string< Args... > format, Args &&... args)
Definition: exceptions.hh:23
Dune::Concept::Geometry auto move_geometry(auto time, const PDELab::Concept::LocalBasis auto &lbasis, const PDELab::Concept::LocalConstContainer auto &lcoeff, auto &local_eqs, auto &fe_cache, const Geometry &geo, auto quad_proj)
Moves a geometry based on the local deformation equations.
Definition: move_geometry.hh:37