/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FEM/assembler/ddg.hh Source File#

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FEM/assembler/ddg.hh Source File
DiFfRG
ddg.hh
Go to the documentation of this file.
1#pragma once
2
3// standard library
4#include <sstream>
5
6// DiFfRG
8
9namespace DiFfRG
10{
11 namespace dDG
12 {
13 using namespace dealii;
14 using std::array;
15
16 template <typename... T> auto fe_tie(T &&...t)
17 {
18 return named_tuple<std::tuple<T &...>,
19 StringSet<"fe_functions", "fe_derivatives", "fe_hessians", "extractors", "variables">>(
20 std::tie(t...));
21 }
22
23 template <typename... T> auto i_tie(T &&...t)
24 {
25 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "fe_derivatives", "fe_hessians">>(
26 std::tie(t...));
27 }
28
29 namespace internal
30 {
37 template <typename Discretization> struct ScratchData {
38 static constexpr uint dim = Discretization::dim;
39 using NumberType = typename Discretization::NumberType;
40
41 ScratchData(const Mapping<dim> &mapping, const FiniteElement<dim> &fe, const Quadrature<dim> &quadrature,
42 const Quadrature<dim - 1> &quadrature_face,
43 const UpdateFlags update_flags = update_values | update_gradients | update_quadrature_points |
44 update_JxW_values | update_hessians,
45 const UpdateFlags interface_update_flags = update_values | update_gradients |
46 update_quadrature_points | update_JxW_values |
47 update_normal_vectors | update_hessians)
48 : n_components(fe.n_components()), fe_values(mapping, fe, quadrature, update_flags),
49 fe_interface_values(mapping, fe, quadrature_face, interface_update_flags)
50 {
51 solution.resize(quadrature.size(), Vector<NumberType>(n_components));
52 solution_grad.resize(quadrature.size(), std::vector<Tensor<1, dim, NumberType>>(n_components));
53 solution_hess.resize(quadrature.size(), std::vector<Tensor<2, dim, NumberType>>(n_components));
54 solution_dot.resize(quadrature.size(), Vector<NumberType>(n_components));
55 solution_interface[0].resize(quadrature_face.size(), Vector<NumberType>(n_components));
56 solution_interface[1].resize(quadrature_face.size(), Vector<NumberType>(n_components));
57 solution_grad_interface[0].resize(quadrature_face.size(),
58 std::vector<Tensor<1, dim, NumberType>>(n_components));
59 solution_grad_interface[1].resize(quadrature_face.size(),
60 std::vector<Tensor<1, dim, NumberType>>(n_components));
61 solution_hess_interface[0].resize(quadrature_face.size(),
62 std::vector<Tensor<2, dim, NumberType>>(n_components));
63 solution_hess_interface[1].resize(quadrature_face.size(),
64 std::vector<Tensor<2, dim, NumberType>>(n_components));
65 comp.resize(fe.n_dofs_per_cell());
66 cached_shape_values.resize(fe.n_dofs_per_cell());
67 cached_shape_grads.resize(fe.n_dofs_per_cell());
68 cached_shape_hessians.resize(fe.n_dofs_per_cell());
69 }
70
72 : n_components(scratch_data.fe_values.get_fe().n_components()),
73 fe_values(scratch_data.fe_values.get_mapping(), scratch_data.fe_values.get_fe(),
74 scratch_data.fe_values.get_quadrature(), scratch_data.fe_values.get_update_flags()),
75 fe_interface_values(scratch_data.fe_interface_values.get_mapping(),
76 scratch_data.fe_interface_values.get_fe(),
77 scratch_data.fe_interface_values.get_quadrature(),
78 scratch_data.fe_interface_values.get_update_flags())
79 {
80 const uint q_size = scratch_data.fe_values.get_quadrature().size();
81 const uint q_face_size = scratch_data.fe_interface_values.get_quadrature().size();
82
83 solution.resize(q_size, Vector<NumberType>(n_components));
84 solution_grad.resize(q_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
85 solution_hess.resize(q_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
86 solution_dot.resize(q_size, Vector<NumberType>(n_components));
87 solution_interface[0].resize(q_face_size, Vector<NumberType>(n_components));
88 solution_interface[1].resize(q_face_size, Vector<NumberType>(n_components));
89 solution_grad_interface[0].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
90 solution_grad_interface[1].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
91 solution_hess_interface[0].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
92 solution_hess_interface[1].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
93 comp.resize(scratch_data.fe_values.get_fe().n_dofs_per_cell());
94 cached_shape_values.resize(scratch_data.fe_values.get_fe().n_dofs_per_cell());
95 cached_shape_grads.resize(scratch_data.fe_values.get_fe().n_dofs_per_cell());
96 cached_shape_hessians.resize(scratch_data.fe_values.get_fe().n_dofs_per_cell());
97 }
98
100
101 FEValues<dim> fe_values;
102 FEInterfaceValues<dim> fe_interface_values;
103
104 std::vector<Vector<NumberType>> solution;
105 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad;
106 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess;
107 std::vector<Vector<NumberType>> solution_dot;
108 array<std::vector<Vector<NumberType>>, 2> solution_interface;
109 array<std::vector<std::vector<Tensor<1, dim, NumberType>>>, 2> solution_grad_interface;
110 array<std::vector<std::vector<Tensor<2, dim, NumberType>>>, 2> solution_hess_interface;
111
112 std::vector<uint> comp;
113 std::vector<double> cached_shape_values;
114 std::vector<Tensor<1, dim>> cached_shape_grads;
115 std::vector<Tensor<2, dim>> cached_shape_hessians;
116 };
117
118 template <typename NumberType> struct CopyData_R {
120 Vector<NumberType> cell_residual;
121 std::vector<types::global_dof_index> joint_dof_indices;
122
123 template <int dim> void reinit(const FEInterfaceValues<dim> &fe_iv)
124 {
125 cell_residual.reinit(fe_iv.n_current_interface_dofs());
126 joint_dof_indices = fe_iv.get_interface_dof_indices();
127 }
128 };
129
130 Vector<NumberType> cell_residual;
131 std::vector<types::global_dof_index> local_dof_indices;
132 std::vector<CopyDataFace_R> face_data;
133
134 template <class Iterator> void reinit(const Iterator &cell, uint dofs_per_cell)
135 {
136 cell_residual.reinit(dofs_per_cell);
137 local_dof_indices.resize(dofs_per_cell);
138 cell->get_dof_indices(local_dof_indices);
139 face_data.clear();
140 face_data.reserve(6);
141 }
142 };
143
144 template <typename NumberType> struct CopyData_J {
146 FullMatrix<NumberType> cell_jacobian;
147 FullMatrix<NumberType> extractor_cell_jacobian;
148 std::vector<types::global_dof_index> joint_dof_indices;
149
150 template <int dim> void reinit(const FEInterfaceValues<dim> &fe_iv, uint n_extractors)
151 {
152 uint dofs_per_cell = fe_iv.n_current_interface_dofs();
153 cell_jacobian.reinit(dofs_per_cell, dofs_per_cell);
154 if (n_extractors > 0) extractor_cell_jacobian.reinit(dofs_per_cell, n_extractors);
155 joint_dof_indices = fe_iv.get_interface_dof_indices();
156 }
157 };
158
159 FullMatrix<NumberType> cell_jacobian;
160 FullMatrix<NumberType> extractor_cell_jacobian;
161 std::vector<types::global_dof_index> local_dof_indices;
162 std::vector<CopyDataFace_J> face_data;
163
164 template <class Iterator> void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
165 {
166 cell_jacobian.reinit(dofs_per_cell, dofs_per_cell);
167 if (n_extractors > 0) extractor_cell_jacobian.reinit(dofs_per_cell, n_extractors);
168 local_dof_indices.resize(dofs_per_cell);
169 cell->get_dof_indices(local_dof_indices);
170 face_data.clear();
171 face_data.reserve(6);
172 }
173 };
174
175 template <typename NumberType> struct CopyData_I {
177 std::array<uint, 2> cell_indices;
178 std::array<double, 2> values;
179 };
180 std::vector<CopyFaceData_I> face_data;
181 double value;
183 };
184 } // namespace internal
185
191 template <typename Discretization_, typename Model_> class Assembler : public FEMAssembler<Discretization_, Model_>
192 {
194
195 public:
196 using Discretization = Discretization_;
197 using Model = Model_;
198 using NumberType = typename Discretization::NumberType;
199 using VectorType = typename Discretization::VectorType;
200
201 using Components = typename Discretization::Components;
202 static constexpr uint dim = Discretization::dim;
203
205 : Base(discretization, model, json),
206 quadrature(fe.degree + 1 + json.get_uint("/discretization/overintegration")),
207 quadrature_face(fe.degree + 1 + json.get_uint("/discretization/overintegration"))
208 {
209 static_assert(Components::count_fe_subsystems() == 1, "A dDG model cannot have multiple submodels!");
210 reinit();
211 }
212
213 virtual void reinit_vector(VectorType &vec) const override { vec.reinit(dof_handler.n_dofs()); }
214
215 virtual void reinit() override
216 {
217 Timer timer;
218
219 Base::reinit();
220
221 // Mass sparsity pattern
222 {
223 DynamicSparsityPattern dsp(dof_handler.n_dofs());
224 DoFTools::make_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
225 /*keep_constrained_dofs = */ true);
226 sparsity_pattern_mass.copy_from(dsp);
228 MatrixCreator::create_mass_matrix(dof_handler, quadrature, mass_matrix, (Function<dim, NumberType> *)nullptr,
229 discretization.get_constraints());
230 }
231 // Jacobian sparsity pattern
232 {
233 DynamicSparsityPattern dsp(dof_handler.n_dofs());
234 DoFTools::make_flux_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
235 /*keep_constrained_dofs = */ true);
236 sparsity_pattern_jacobian.copy_from(dsp);
237 }
238
239 timings_reinit.push_back(timer.wall_time());
240 }
241
242 virtual void rebuild_jacobian_sparsity() override
243 {
244 // Jacobian sparsity pattern
245 DynamicSparsityPattern dsp(dof_handler.n_dofs());
246 DoFTools::make_flux_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
247 /*keep_constrained_dofs = */ true);
248 for (uint row = 0; row < dsp.n_rows(); ++row)
249 for (const auto &col : extractor_dof_indices)
250 dsp.add(row, col);
251 sparsity_pattern_jacobian.copy_from(dsp);
252 }
253
254 virtual const SparsityPattern &get_sparsity_pattern_jacobian() const override
255 {
257 }
258 virtual const SparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
259
266 virtual void refinement_indicator(Vector<double> &indicator, const VectorType &solution_global) override
267 {
268 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
270 using CopyData = internal::CopyData_I<NumberType>;
271
272 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData &copy_data) {
273 scratch_data.fe_values.reinit(t_cell);
274 const auto &fe_v = scratch_data.fe_values;
275 copy_data.cell_index = t_cell->active_cell_index();
276 copy_data.value = 0;
277
278 const auto &JxW = fe_v.get_JxW_values();
279 const auto &q_points = fe_v.get_quadrature_points();
280 const auto &q_indices = fe_v.quadrature_point_indices();
281
282 auto &solution = scratch_data.solution;
283 auto &solution_grad = scratch_data.solution_grad;
284 auto &solution_hess = scratch_data.solution_hess;
285 fe_v.get_function_values(solution_global, solution);
286 fe_v.get_function_gradients(solution_global, solution_grad);
287 fe_v.get_function_hessians(solution_global, solution_hess);
288
289 double local_indicator = 0.;
290 for (const auto &q_index : q_indices) {
291 const auto &x_q = q_points[q_index];
292 model.cell_indicator(local_indicator, x_q,
293 i_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index]));
294 copy_data.value += JxW[q_index] * local_indicator;
295 }
296 };
297 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
298 const uint &nf, const unsigned int &nsf, Scratch &scratch_data,
299 CopyData &copy_data) {
300 scratch_data.fe_interface_values.reinit(t_cell, f, sf, t_ncell, nf, nsf);
301 const auto &fe_iv = scratch_data.fe_interface_values;
302 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
303 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
304
305 auto &copy_data_face = copy_data.face_data.emplace_back();
306 copy_data_face.cell_indices[0] = t_cell->active_cell_index();
307 copy_data_face.cell_indices[1] = t_ncell->active_cell_index();
308 copy_data_face.values[0] = 0;
309 copy_data_face.values[1] = 0;
310
311 const auto &JxW = fe_iv.get_JxW_values();
312 const auto &q_points = fe_iv.get_quadrature_points();
313 const auto &q_indices = fe_iv.quadrature_point_indices();
314 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
315
316 auto &solution_s = scratch_data.solution_interface[0];
317 auto &solution_n = scratch_data.solution_interface[1];
318 auto &solution_grad_s = scratch_data.solution_grad_interface[0];
319 auto &solution_grad_n = scratch_data.solution_grad_interface[1];
320 auto &solution_hess_s = scratch_data.solution_hess_interface[0];
321 auto &solution_hess_n = scratch_data.solution_hess_interface[1];
322 fe_iv_s.get_function_values(solution_global, solution_s);
323 fe_iv_n.get_function_values(solution_global, solution_n);
324 fe_iv_s.get_function_gradients(solution_global, solution_grad_s);
325 fe_iv_n.get_function_gradients(solution_global, solution_grad_n);
326 fe_iv_s.get_function_hessians(solution_global, solution_hess_s);
327 fe_iv_n.get_function_hessians(solution_global, solution_hess_n);
328
329 array<double, 2> local_indicator{};
330 for (const auto &q_index : q_indices) {
331 const auto &x_q = q_points[q_index];
332 model.face_indicator(local_indicator, normals[q_index], x_q,
333 i_tie(solution_s[q_index], solution_grad_s[q_index], solution_hess_s[q_index]),
334 i_tie(solution_n[q_index], solution_grad_n[q_index], solution_hess_n[q_index]));
335
336 copy_data_face.values[0] += JxW[q_index] * local_indicator[0] * (1. + t_cell->at_boundary());
337 copy_data_face.values[1] += JxW[q_index] * local_indicator[1] * (1. + t_ncell->at_boundary());
338 }
339 };
340 const auto copier = [&](const CopyData &c) {
341 for (auto &cdf : c.face_data)
342 for (uint j = 0; j < 2; ++j)
343 indicator[cdf.cell_indices[j]] += cdf.values[j];
344 indicator[c.cell_index] += c.value;
345 };
346
347 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
348 CopyData copy_data;
349 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_own_interior_faces_once;
350
351 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
352 copy_data, flags, nullptr, face_worker, threads, batch_size);
353 }
354
355 virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot,
356 NumberType weight) override
357 {
358 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
360 using CopyData = internal::CopyData_R<NumberType>;
361 const auto &constraints = discretization.get_constraints();
362
363 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
364 scratch_data.fe_values.reinit(cell);
365 const auto &fe_v = scratch_data.fe_values;
366 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
367
368 copy_data.reinit(cell, n_dofs);
369 const auto &JxW = fe_v.get_JxW_values();
370 const auto &q_points = fe_v.get_quadrature_points();
371 const auto &q_indices = fe_v.quadrature_point_indices();
372
373 auto &solution = scratch_data.solution;
374 auto &solution_dot = scratch_data.solution_dot;
375 fe_v.get_function_values(solution_global, solution);
376 fe_v.get_function_values(solution_global_dot, solution_dot);
377
378 auto &comp = scratch_data.comp;
379 for (uint i = 0; i < n_dofs; ++i)
380 comp[i] = fe.system_to_component_index(i).first;
381
382 array<NumberType, Components::count_fe_functions()> mass{};
383 for (const auto &q_index : q_indices) {
384 const auto &x_q = q_points[q_index];
385 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
386
387 for (uint i = 0; i < n_dofs; ++i) {
388 const auto component_i = comp[i];
389 copy_data.cell_residual(i) += weight * JxW[q_index] *
390 fe_v.shape_value_component(i, q_index, component_i) *
391 mass[component_i]; // +phi_i(x_q) * mass(x_q, u_q)
392 }
393 }
394 };
395 const auto copier = [&](const CopyData &c) {
396 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, mass);
397 };
398
399 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
400 CopyData copy_data;
401 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
402
403 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
404 copy_data, flags, nullptr, nullptr, threads, batch_size);
405 }
406
407 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
408 const VectorType &solution_global_dot, NumberType weight_mass,
409 const VectorType &variables = VectorType()) override
410 {
411 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
413 using CopyData = internal::CopyData_R<NumberType>;
414 const auto &constraints = discretization.get_constraints();
415
416 // Find the EoM and extract whatever data is needed for the model.
417 std::array<NumberType, Components::count_extractors()> __extracted_data{{}};
418 if constexpr (Components::count_extractors() > 0)
419 this->extract(__extracted_data, solution_global, variables, true, false, true);
420 const auto &extracted_data = __extracted_data;
421
422 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
423 scratch_data.fe_values.reinit(cell);
424 const auto &fe_v = scratch_data.fe_values;
425 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
426
427 copy_data.reinit(cell, n_dofs);
428 const auto &JxW = fe_v.get_JxW_values();
429 const auto &q_points = fe_v.get_quadrature_points();
430 const auto &q_indices = fe_v.quadrature_point_indices();
431
432 auto &solution = scratch_data.solution;
433 auto &solution_grad = scratch_data.solution_grad;
434 auto &solution_hess = scratch_data.solution_hess;
435 auto &solution_dot = scratch_data.solution_dot;
436 fe_v.get_function_values(solution_global, solution);
437 fe_v.get_function_gradients(solution_global, solution_grad);
438 fe_v.get_function_hessians(solution_global, solution_hess);
439 fe_v.get_function_values(solution_global_dot, solution_dot);
440
441 auto &comp = scratch_data.comp;
442 for (uint i = 0; i < n_dofs; ++i)
443 comp[i] = fe.system_to_component_index(i).first;
444
445 array<Tensor<1, dim, NumberType>, Components::count_fe_functions()> flux{};
446 array<NumberType, Components::count_fe_functions()> source{};
447 array<NumberType, Components::count_fe_functions()> mass{};
448 for (const auto &q_index : q_indices) {
449 const auto &x_q = q_points[q_index];
450 model.flux(
451 flux, x_q,
452 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
453 model.source(
454 source, x_q,
455 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
456 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
457
458 for (uint i = 0; i < n_dofs; ++i) {
459 const auto component_i = comp[i];
460 copy_data.cell_residual(i) += JxW[q_index] * weight * // dx *
461 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
462 flux[component_i]) // -dphi_i(x_q) * flux(x_q, u_q)
463 + fe_v.shape_value_component(i, q_index, component_i) *
464 source[component_i]); // -phi_i(x_q) * source(x_q, u_q)
465 copy_data.cell_residual(i) += weight_mass * JxW[q_index] *
466 fe_v.shape_value_component(i, q_index, component_i) *
467 mass[component_i]; // +phi_i(x_q) * mass(x_q, u_q)
468 }
469 }
470 };
471 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
472 CopyData &copy_data) {
473 scratch_data.fe_interface_values.reinit(cell, face_no);
474 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
475 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
476
477 const auto &JxW = fe_fv.get_JxW_values();
478 const auto &q_points = fe_fv.get_quadrature_points();
479 const auto &q_indices = fe_fv.quadrature_point_indices();
480 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
481
482 auto &solution = scratch_data.solution_interface[0];
483 auto &solution_grad = scratch_data.solution_grad_interface[0];
484 auto &solution_hess = scratch_data.solution_hess_interface[0];
485 fe_fv.get_function_values(solution_global, solution);
486 fe_fv.get_function_gradients(solution_global, solution_grad);
487 fe_fv.get_function_hessians(solution_global, solution_hess);
488
489 auto &comp = scratch_data.comp;
490 for (uint i = 0; i < n_dofs; ++i)
491 comp[i] = fe.system_to_component_index(i).first;
492
493 array<Tensor<1, dim, NumberType>, Components::count_fe_functions()> numflux{};
494 for (const auto &q_index : q_indices) {
495 const auto &x_q = q_points[q_index];
496 model.boundary_numflux(
497 numflux, normals[q_index], x_q,
498 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
499
500 for (uint i = 0; i < n_dofs; ++i) {
501 const auto component_i = comp[i];
502 copy_data.cell_residual(i) +=
503 weight * JxW[q_index] * // dx
504 (fe_fv.shape_value_component(i, q_index, component_i) *
505 scalar_product(numflux[component_i], normals[q_index])); // phi_i(x_q) * numflux(x_q, u_q) * n(x_q)
506 }
507 }
508 };
509 const auto face_worker = [&](const Iterator &cell, const uint &f, const unsigned int &sf, const Iterator &ncell,
510 const unsigned int &nf, const unsigned int &nsf, Scratch &scratch_data,
511 CopyData &copy_data) {
512 scratch_data.fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
513 const auto &fe_iv = scratch_data.fe_interface_values;
514 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
515 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
516 const uint n_dofs = fe_iv.n_current_interface_dofs();
517
518 copy_data.face_data.emplace_back();
519 auto &copy_data_face = copy_data.face_data.back();
520 copy_data_face.reinit(fe_iv);
521
522 const auto &JxW = fe_iv.get_JxW_values();
523 const auto &q_points = fe_iv.get_quadrature_points();
524 const auto &q_indices = fe_iv.quadrature_point_indices();
525 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
526
527 auto &solution_s = scratch_data.solution_interface[0];
528 auto &solution_n = scratch_data.solution_interface[1];
529 auto &solution_grad_s = scratch_data.solution_grad_interface[0];
530 auto &solution_grad_n = scratch_data.solution_grad_interface[1];
531 auto &solution_hess_s = scratch_data.solution_hess_interface[0];
532 auto &solution_hess_n = scratch_data.solution_hess_interface[1];
533 fe_iv_s.get_function_values(solution_global, solution_s);
534 fe_iv_n.get_function_values(solution_global, solution_n);
535 fe_iv_s.get_function_gradients(solution_global, solution_grad_s);
536 fe_iv_n.get_function_gradients(solution_global, solution_grad_n);
537 fe_iv_s.get_function_hessians(solution_global, solution_hess_s);
538 fe_iv_n.get_function_hessians(solution_global, solution_hess_n);
539
540 auto &comp = scratch_data.comp;
541 comp.resize(n_dofs);
542 for (uint i = 0; i < n_dofs; ++i) {
543 const auto &cd_i = fe_iv.interface_dof_to_dof_indices(i);
544 comp[i] = cd_i[0] == numbers::invalid_unsigned_int ? fe.system_to_component_index(cd_i[1]).first
545 : fe.system_to_component_index(cd_i[0]).first;
546 }
547
548 array<Tensor<1, dim, NumberType>, Components::count_fe_functions()> numflux{};
549 for (const auto &q_index : q_indices) {
550 const auto &x_q = q_points[q_index];
551 model.numflux(numflux, normals[q_index], x_q,
552 fe_tie(solution_s[q_index], solution_grad_s[q_index], solution_hess_s[q_index],
553 extracted_data, variables),
554 fe_tie(solution_n[q_index], solution_grad_n[q_index], solution_hess_n[q_index],
555 extracted_data, variables));
556
557 for (uint i = 0; i < n_dofs; ++i) {
558 const auto component_i = comp[i];
559 copy_data_face.cell_residual(i) +=
560 weight * JxW[q_index] * // dx
561 (fe_iv.jump_in_shape_values(i, q_index, component_i) *
562 scalar_product(numflux[component_i],
563 normals[q_index])); // [[phi_i(x_q)]] * numflux(x_q, u_q) * n(x_q)
564 }
565 }
566 };
567 const auto copier = [&](const CopyData &c) {
568 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, residual);
569 for (auto &cdf : c.face_data)
570 constraints.distribute_local_to_global(cdf.cell_residual, cdf.joint_dof_indices, residual);
571 };
572
573 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
574 CopyData copy_data;
575 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
576 MeshWorker::assemble_own_interior_faces_once;
577
578 Timer timer;
579 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
580 copy_data, flags, boundary_worker, face_worker, threads, batch_size);
581 timings_residual.push_back(timer.wall_time());
582 }
583
584 virtual void jacobian_mass(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
585 const VectorType &solution_global_dot, NumberType alpha, NumberType beta) override
586 {
587 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
589 using CopyData = internal::CopyData_J<NumberType>;
590 const auto &constraints = discretization.get_constraints();
591
592 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
593 scratch_data.fe_values.reinit(cell);
594 const auto &fe_v = scratch_data.fe_values;
595 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
596
597 copy_data.reinit(cell, n_dofs, Components::count_extractors());
598 const auto &JxW = fe_v.get_JxW_values();
599 const auto &q_points = fe_v.get_quadrature_points();
600 const auto &q_indices = fe_v.quadrature_point_indices();
601 auto &solution = scratch_data.solution;
602 auto &solution_dot = scratch_data.solution_dot;
603
604 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass;
605 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass_dot;
606
607 fe_v.get_function_values(solution_global, solution);
608 fe_v.get_function_values(solution_global_dot, solution_dot);
609
610 auto &comp = scratch_data.comp;
611 for (uint i = 0; i < n_dofs; ++i)
612 comp[i] = fe.system_to_component_index(i).first;
613
614 for (const auto &q_index : q_indices) {
615 const auto &x_q = q_points[q_index];
616 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
617 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
618
619 for (uint i = 0; i < n_dofs; ++i) {
620 const auto component_i = comp[i];
621 for (uint j = 0; j < n_dofs; ++j) {
622 const auto component_j = comp[j];
623 copy_data.cell_jacobian(i, j) +=
624 JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) *
625 fe_v.shape_value_component(i, q_index, component_i) *
626 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
627 }
628 }
629 }
630 };
631 const auto copier = [&](const CopyData &c) {
632 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
633 };
634
635 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
636 CopyData copy_data;
637 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
638
639 Timer timer;
640 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
641 copy_data, flags, nullptr, nullptr, threads, batch_size);
642 timings_jacobian.push_back(timer.wall_time());
643 }
644
645 virtual void jacobian(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global, NumberType weight,
646 const VectorType &solution_global_dot, NumberType alpha, NumberType beta,
647 const VectorType &variables = VectorType()) override
648 {
649 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
651 using CopyData = internal::CopyData_J<NumberType>;
652 const auto &constraints = discretization.get_constraints();
653
654 // Find the EoM and extract whatever data is needed for the model.
655 std::array<NumberType, Components::count_extractors()> __extracted_data{{}};
656 if constexpr (Components::count_extractors() > 0) {
657 this->extract(__extracted_data, solution_global, variables, true, true, true);
658 if (this->jacobian_extractors(this->extractor_jacobian, solution_global, variables))
660 }
661 const auto &extracted_data = __extracted_data;
662
663 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
664 scratch_data.fe_values.reinit(cell);
665 const auto &fe_v = scratch_data.fe_values;
666 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
667
668 copy_data.reinit(cell, n_dofs, Components::count_extractors());
669 const auto &JxW = fe_v.get_JxW_values();
670 const auto &q_points = fe_v.get_quadrature_points();
671 const auto &q_indices = fe_v.quadrature_point_indices();
672
673 auto &solution = scratch_data.solution;
674 auto &solution_dot = scratch_data.solution_dot;
675 auto &solution_grad = scratch_data.solution_grad;
676 auto &solution_hess = scratch_data.solution_hess;
677
678 fe_v.get_function_values(solution_global, solution);
679 fe_v.get_function_values(solution_global_dot, solution_dot);
680 fe_v.get_function_gradients(solution_global, solution_grad);
681 fe_v.get_function_hessians(solution_global, solution_hess);
682
683 auto &comp = scratch_data.comp;
684 for (uint i = 0; i < n_dofs; ++i)
685 comp[i] = fe.system_to_component_index(i).first;
686
687 auto &cached_shape_values = scratch_data.cached_shape_values;
688 auto &cached_shape_grads = scratch_data.cached_shape_grads;
689 auto &cached_shape_hessians = scratch_data.cached_shape_hessians;
690
691 SimpleMatrix<Tensor<1, dim>, Components::count_fe_functions()> j_flux;
692 SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()> j_grad_flux;
693 SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()> j_hess_flux;
694 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
695 j_extr_flux;
696 SimpleMatrix<NumberType, Components::count_fe_functions()> j_source;
697 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions()> j_grad_source;
698 SimpleMatrix<Tensor<2, dim, NumberType>, Components::count_fe_functions()> j_hess_source;
699 SimpleMatrix<NumberType, Components::count_fe_functions(), Components::count_extractors()> j_extr_source;
700 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass;
701 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass_dot;
702
703 for (const auto &q_index : q_indices) {
704 const auto &x_q = q_points[q_index];
705 model.template jacobian_flux_source<0, 0>(
706 j_flux, j_source, x_q,
707 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
708 model.template jacobian_flux_source_grad<1>(
709 j_grad_flux, j_grad_source, x_q,
710 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
711 model.template jacobian_flux_source_hess<2>(
712 j_hess_flux, j_hess_source, x_q,
713 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
714 if constexpr (Components::count_extractors() > 0) {
715 model.template jacobian_flux_source_extr<3>(
716 j_extr_flux, j_extr_source, x_q,
717 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
718 }
719 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
720 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
721
722 // Cache shape values, gradients, and hessians for all DoFs at this quadrature point
723 for (uint k = 0; k < n_dofs; ++k) {
724 cached_shape_values[k] = fe_v.shape_value_component(k, q_index, comp[k]);
725 cached_shape_grads[k] = fe_v.shape_grad_component(k, q_index, comp[k]);
726 cached_shape_hessians[k] = fe_v.shape_hessian_component(k, q_index, comp[k]);
727 }
728
729 for (uint i = 0; i < n_dofs; ++i) {
730 const auto component_i = comp[i];
731 const auto &shape_value_i = cached_shape_values[i];
732 const auto &shape_grad_i = cached_shape_grads[i];
733 for (uint j = 0; j < n_dofs; ++j) {
734 const auto component_j = comp[j];
735 const auto &shape_value_j = cached_shape_values[j];
736 const auto &shape_grad_j = cached_shape_grads[j];
737 const auto &shape_hessian_j = cached_shape_hessians[j];
738 // consolidated contribution: scalar + gradient + hessian + mass
739 copy_data.cell_jacobian(i, j) +=
740 weight * JxW[q_index] *
741 (shape_value_j * // dx * phi_j * (
742 (-scalar_product(shape_grad_i,
743 j_flux(component_i, component_j)) // -dphi_i * jflux
744 + shape_value_i * j_source(component_i, component_j)) // -phi_i * jsource)
745 + scalar_product(shape_grad_j, // gradient contribution
746 -scalar_product(shape_grad_i,
747 j_grad_flux(component_i, component_j)) // -dphi_i * jflux
748 + shape_value_i *
749 j_grad_source(component_i, component_j)) // -phi_i * jsource
750 + scalar_product(shape_hessian_j, // hessian contribution
751 -scalar_product(shape_grad_i,
752 j_hess_flux(component_i, component_j)) +
753 shape_value_i * j_hess_source(component_i, component_j))) +
754 JxW[q_index] * shape_value_j * shape_value_i * // mass contribution
755 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
756 }
757 // extractor contribution
758 if constexpr (Components::count_extractors() > 0)
759 for (uint e = 0; e < Components::count_extractors(); ++e)
760 copy_data.extractor_cell_jacobian(i, e) +=
761 weight * JxW[q_index] * // dx * phi_j * (
762 (-scalar_product(shape_grad_i,
763 j_extr_flux(component_i, e)) // -dphi_i * jflux
764 + shape_value_i *
765 j_extr_source(component_i, e)); // -phi_i * jsource)
766 }
767 }
768 };
769 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
770 CopyData &copy_data) {
771 scratch_data.fe_interface_values.reinit(cell, face_no);
772 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
773 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
774
775 const auto &JxW = fe_fv.get_JxW_values();
776 const auto &q_points = fe_fv.get_quadrature_points();
777 const auto &q_indices = fe_fv.quadrature_point_indices();
778 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
779
780 auto &solution = scratch_data.solution_interface[0];
781 auto &solution_grad = scratch_data.solution_grad_interface[0];
782 auto &solution_hess = scratch_data.solution_hess_interface[0];
783 fe_fv.get_function_values(solution_global, solution);
784 fe_fv.get_function_gradients(solution_global, solution_grad);
785 fe_fv.get_function_hessians(solution_global, solution_hess);
786
787 auto &comp = scratch_data.comp;
788 for (uint i = 0; i < n_dofs; ++i)
789 comp[i] = fe.system_to_component_index(i).first;
790
791 auto &cached_shape_values = scratch_data.cached_shape_values;
792 auto &cached_shape_grads = scratch_data.cached_shape_grads;
793 auto &cached_shape_hessians = scratch_data.cached_shape_hessians;
794
795 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions()> j_boundary_numflux;
796 SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()>
797 j_grad_boundary_numflux;
798 SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()>
799 j_hess_boundary_numflux;
800 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
801 j_extr_boundary_numflux;
802 for (const auto &q_index : q_indices) {
803 const auto &x_q = q_points[q_index];
804 model.template jacobian_boundary_numflux<0, 0>(
805 j_boundary_numflux, normals[q_index], x_q,
806 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
807 model.template jacobian_boundary_numflux_grad<1>(
808 j_grad_boundary_numflux, normals[q_index], x_q,
809 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
810 model.template jacobian_boundary_numflux_hess<2>(
811 j_hess_boundary_numflux, normals[q_index], x_q,
812 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
813 if constexpr (Components::count_extractors() > 0) {
814 model.template jacobian_boundary_numflux_extr<3>(
815 j_extr_boundary_numflux, normals[q_index], x_q,
816 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
817 }
818
819 // Cache shape values, gradients, and hessians for all DoFs at this quadrature point
820 for (uint k = 0; k < n_dofs; ++k) {
821 cached_shape_values[k] = fe_fv.shape_value_component(k, q_index, comp[k]);
822 cached_shape_grads[k] = fe_fv.shape_grad_component(k, q_index, comp[k]);
823 cached_shape_hessians[k] = fe_fv.shape_hessian_component(k, q_index, comp[k]);
824 }
825
826 for (uint i = 0; i < n_dofs; ++i) {
827 const auto component_i = comp[i];
828 const auto &shape_value_i = cached_shape_values[i];
829 for (uint j = 0; j < n_dofs; ++j) {
830 const auto component_j = comp[j];
831 const auto &shape_value_j = cached_shape_values[j];
832 const auto &shape_grad_j = cached_shape_grads[j];
833 const auto &shape_hessian_j = cached_shape_hessians[j];
834 // consolidated contribution: scalar + gradient + hessian
835 copy_data.cell_jacobian(i, j) +=
836 weight * JxW[q_index] *
837 (shape_value_j * // dx * phi_j(x_q)
838 (shape_value_i *
839 scalar_product(j_boundary_numflux(component_i, component_j),
840 normals[q_index])) // phi_i(x_q) * j_numflux(x_q, u_q) * n(x_q)
841 + scalar_product(shape_grad_j, // gradient contribution
842 shape_value_i *
843 scalar_product(j_grad_boundary_numflux(component_i, component_j),
844 normals[q_index]))
845 + scalar_product(shape_hessian_j, // hessian contribution
846 shape_value_i *
847 scalar_product(j_hess_boundary_numflux(component_i, component_j),
848 normals[q_index])));
849 }
850 // extractor contribution
851 if constexpr (Components::count_extractors() > 0)
852 for (uint e = 0; e < Components::count_extractors(); ++e)
853 copy_data.extractor_cell_jacobian(i, e) +=
854 weight * JxW[q_index] * // dx * phi_j(x_q)
855 (shape_value_i *
856 scalar_product(j_extr_boundary_numflux(component_i, e),
857 normals[q_index])); // phi_i(x_q) * j_numflux(x_q, u_q) * n(x_q)
858 }
859 }
860 };
861 const auto face_worker = [&](const Iterator &cell, const uint &f, const uint &sf, const Iterator &ncell,
862 const uint &nf, const uint &nsf, Scratch &scratch_data, CopyData &copy_data) {
863 scratch_data.fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
864 const auto &fe_iv = scratch_data.fe_interface_values;
865 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
866 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
867 const uint n_dofs = fe_iv.n_current_interface_dofs();
868
869 copy_data.face_data.emplace_back();
870 auto &copy_data_face = copy_data.face_data.back();
871 copy_data_face.reinit(fe_iv, Components::count_extractors());
872
873 const auto &JxW = fe_iv.get_JxW_values();
874 const auto &q_points = fe_iv.get_quadrature_points();
875 const auto &q_indices = fe_iv.quadrature_point_indices();
876 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
877
878 auto &solution_s = scratch_data.solution_interface[0];
879 auto &solution_n = scratch_data.solution_interface[1];
880 auto &solution_grad_s = scratch_data.solution_grad_interface[0];
881 auto &solution_grad_n = scratch_data.solution_grad_interface[1];
882 auto &solution_hess_s = scratch_data.solution_hess_interface[0];
883 auto &solution_hess_n = scratch_data.solution_hess_interface[1];
884 fe_iv_s.get_function_values(solution_global, solution_s);
885 fe_iv_n.get_function_values(solution_global, solution_n);
886 fe_iv_s.get_function_gradients(solution_global, solution_grad_s);
887 fe_iv_n.get_function_gradients(solution_global, solution_grad_n);
888 fe_iv_s.get_function_hessians(solution_global, solution_hess_s);
889 fe_iv_n.get_function_hessians(solution_global, solution_hess_n);
890
891 // Pre-compute interface DoF component indices and face numbers
892 auto &comp = scratch_data.comp;
893 comp.resize(n_dofs);
894 std::vector<uint> face_no_arr(n_dofs);
895 std::vector<uint> local_dof_arr(n_dofs);
896 for (uint i = 0; i < n_dofs; ++i) {
897 const auto &cd_i = fe_iv.interface_dof_to_dof_indices(i);
898 face_no_arr[i] = cd_i[0] == numbers::invalid_unsigned_int ? 1 : 0;
899 local_dof_arr[i] = cd_i[face_no_arr[i]];
900 comp[i] = fe.system_to_component_index(local_dof_arr[i]).first;
901 }
902
903 auto &cached_shape_values = scratch_data.cached_shape_values;
904 auto &cached_shape_grads = scratch_data.cached_shape_grads;
905 auto &cached_shape_hessians = scratch_data.cached_shape_hessians;
906 cached_shape_values.resize(n_dofs);
907 cached_shape_grads.resize(n_dofs);
908 cached_shape_hessians.resize(n_dofs);
909
910 array<SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions()>, 2> j_numflux;
911 array<SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()>, 2>
912 j_grad_numflux;
913 array<SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()>, 2>
914 j_hess_numflux;
915 array<SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(),
916 Components::count_extractors()>,
917 2>
918 j_extr_numflux;
919 for (const auto &q_index : q_indices) {
920 const auto &x_q = q_points[q_index];
921 model.template jacobian_numflux<0, 0>(j_numflux, normals[q_index], x_q,
922 fe_tie(solution_s[q_index], solution_grad_s[q_index],
923 solution_hess_s[q_index], extracted_data, variables),
924 fe_tie(solution_n[q_index], solution_grad_n[q_index],
925 solution_hess_n[q_index], extracted_data, variables));
926 model.template jacobian_numflux_grad<1>(j_grad_numflux, normals[q_index], x_q,
927 fe_tie(solution_s[q_index], solution_grad_s[q_index],
928 solution_hess_s[q_index], extracted_data, variables),
929 fe_tie(solution_n[q_index], solution_grad_n[q_index],
930 solution_hess_n[q_index], extracted_data, variables));
931 model.template jacobian_numflux_hess<2>(j_hess_numflux, normals[q_index], x_q,
932 fe_tie(solution_s[q_index], solution_grad_s[q_index],
933 solution_hess_s[q_index], extracted_data, variables),
934 fe_tie(solution_n[q_index], solution_grad_n[q_index],
935 solution_hess_n[q_index], extracted_data, variables));
936 if constexpr (Components::count_extractors() > 0) {
937 model.template jacobian_numflux_extr<3>(j_extr_numflux, normals[q_index], x_q,
938 fe_tie(solution_s[q_index], solution_grad_s[q_index],
939 solution_hess_s[q_index], extracted_data, variables),
940 fe_tie(solution_n[q_index], solution_grad_n[q_index],
941 solution_hess_n[q_index], extracted_data, variables));
942 }
943
944 // Cache shape values, gradients, and hessians for all interface DoFs at this quadrature point
945 for (uint k = 0; k < n_dofs; ++k) {
946 cached_shape_values[k] =
947 fe_iv.get_fe_face_values(face_no_arr[k]).shape_value_component(local_dof_arr[k], q_index, comp[k]);
948 cached_shape_grads[k] =
949 fe_iv.get_fe_face_values(face_no_arr[k]).shape_grad_component(local_dof_arr[k], q_index, comp[k]);
950 cached_shape_hessians[k] =
951 fe_iv.get_fe_face_values(face_no_arr[k]).shape_hessian_component(local_dof_arr[k], q_index, comp[k]);
952 }
953
954 for (uint i = 0; i < n_dofs; ++i) {
955 const auto component_i = comp[i];
956 const auto face_no_i = face_no_arr[i];
957 const auto jump_i = fe_iv.jump_in_shape_values(i, q_index, component_i);
958 for (uint j = 0; j < n_dofs; ++j) {
959 const auto component_j = comp[j];
960 const auto face_no_j = face_no_arr[j];
961 const auto &shape_value_j = cached_shape_values[j];
962 const auto &shape_grad_j = cached_shape_grads[j];
963 const auto &shape_hessian_j = cached_shape_hessians[j];
964 // consolidated contribution: scalar + gradient + hessian
965 copy_data_face.cell_jacobian(i, j) +=
966 weight * JxW[q_index] *
967 (shape_value_j * // dx * phi_j(x_q)
968 (jump_i *
969 scalar_product(j_numflux[face_no_j](component_i, component_j),
970 normals[q_index])) // [[phi_i(x_q)]] * j_numflux(x_q, u_q)
971 + scalar_product(shape_grad_j, // gradient contribution
972 jump_i *
973 scalar_product(j_grad_numflux[face_no_j](component_i, component_j),
974 normals[q_index]))
975 + scalar_product(shape_hessian_j, // hessian contribution
976 jump_i *
977 scalar_product(j_hess_numflux[face_no_j](component_i, component_j),
978 normals[q_index])));
979 }
980 // extractor contribution
981 if constexpr (Components::count_extractors() > 0)
982 for (uint e = 0; e < Components::count_extractors(); ++e)
983 copy_data_face.extractor_cell_jacobian(i, e) +=
984 weight * JxW[q_index] * // dx * phi_j(x_q)
985 (jump_i *
986 scalar_product(j_extr_numflux[face_no_i](component_i, e),
987 normals[q_index])); // [[phi_i(x_q)]] * j_numflux(x_q, u_q)
988 }
989 }
990 };
991 const auto copier = [&](const CopyData &c) {
992 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
993 for (auto &cdf : c.face_data) {
994 constraints.distribute_local_to_global(cdf.cell_jacobian, cdf.joint_dof_indices, jacobian);
995 if constexpr (Components::count_extractors() > 0) {
996 FullMatrix<NumberType> extractor_dependence(cdf.joint_dof_indices.size(), extractor_dof_indices.size());
997 cdf.extractor_cell_jacobian.mmult(extractor_dependence, this->extractor_jacobian);
998 constraints.distribute_local_to_global(extractor_dependence, cdf.joint_dof_indices, extractor_dof_indices,
999 jacobian);
1000 }
1001 }
1002 if constexpr (Components::count_extractors() > 0) {
1003 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices.size(), extractor_dof_indices.size());
1004 c.extractor_cell_jacobian.mmult(extractor_dependence, this->extractor_jacobian);
1005 constraints.distribute_local_to_global(extractor_dependence, c.local_dof_indices, extractor_dof_indices,
1006 jacobian);
1007 }
1008 };
1009
1010 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
1011 CopyData copy_data;
1012 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
1013 MeshWorker::assemble_own_interior_faces_once;
1014
1015 Timer timer;
1016 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
1017 copy_data, flags, boundary_worker, face_worker, threads, batch_size);
1018 timings_jacobian.push_back(timer.wall_time());
1019 }
1020
1021 void log(const std::string logger)
1022 {
1023 std::stringstream ss;
1024 ss << "dDG Assembler: " << std::endl;
1025 ss << " Reinit: " << average_time_reinit() * 1000 << "ms (" << num_reinits() << ")" << std::endl;
1026 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
1027 << std::endl;
1028 ss << " Jacobian: " << average_time_jacobian_assembly() * 1000 << "ms (" << num_jacobians() << ")"
1029 << std::endl;
1030 spdlog::get(logger)->info(ss.str());
1031 }
1032
1033 double average_time_reinit() const
1034 {
1035 double t = 0.;
1036 double n = timings_reinit.size();
1037 for (const auto &t_ : timings_reinit)
1038 t += t_ / n;
1039 return t;
1040 }
1041 uint num_reinits() const { return timings_reinit.size(); }
1042
1044 {
1045 double t = 0.;
1046 double n = timings_residual.size();
1047 for (const auto &t_ : timings_residual)
1048 t += t_ / n;
1049 return t;
1050 }
1051 uint num_residuals() const { return timings_residual.size(); }
1052
1054 {
1055 double t = 0.;
1056 double n = timings_jacobian.size();
1057 for (const auto &t_ : timings_jacobian)
1058 t += t_ / n;
1059 return t;
1060 }
1061 uint num_jacobians() const { return timings_jacobian.size(); }
1062
1063 protected:
1065 using Base::dof_handler;
1066 using Base::fe;
1067 using Base::mapping;
1068 using Base::model;
1069
1070 QGauss<dim> quadrature;
1072 using Base::batch_size;
1073 using Base::threads;
1074
1075 SparsityPattern sparsity_pattern_mass;
1077 SparseMatrix<NumberType> mass_matrix;
1078
1079 std::vector<double> timings_reinit;
1080 std::vector<double> timings_residual;
1081 std::vector<double> timings_jacobian;
1082
1084 };
1085 } // namespace dDG
1086} // namespace DiFfRG
The basic assembler that can be used for any standard CG scheme with flux and source.
Definition common.hh:39
Model & model
Definition common.hh:338
const Mapping< dim > & mapping
Definition common.hh:341
FullMatrix< NumberType > extractor_jacobian
Definition common.hh:351
const FiniteElement< dim > & fe
Definition common.hh:339
void extract(std::array< NumberType, Components::count_extractors()> &data, const VectorType &solution_global, const VectorType &variables, bool search_EoM, bool set_EoM, bool postprocess) const
Definition common.hh:200
bool jacobian_extractors(FullMatrix< NumberType > &extractor_jacobian, const VectorType &solution_global, const VectorType &variables)
Definition common.hh:237
uint threads
Definition common.hh:343
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition common.hh:108
const DoFHandler< dim > & dof_handler
Definition common.hh:340
Discretization & discretization
Definition common.hh:337
uint batch_size
Definition common.hh:344
std::vector< types::global_dof_index > extractor_dof_indices
Definition common.hh:355
A wrapper around the boost json value class.
Definition json.hh:19
Definition quadrature.hh:56
size_t size() const
A simple NxM-matrix class, which is used for cell-wise Jacobians.
Definition tuples.hh:156
The basic assembler that can be used for any standard DG scheme with flux and source.
Definition ddg.hh:192
Model & model
Definition common.hh:338
uint num_residuals() const
Definition ddg.hh:1051
SparsityPattern sparsity_pattern_jacobian
Definition ddg.hh:1076
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition ddg.hh:215
const Mapping< dim > & mapping
Definition common.hh:341
Discretization_ Discretization
Definition ddg.hh:196
virtual void refinement_indicator(Vector< double > &indicator, const VectorType &solution_global) override
refinement indicator for adaptivity. Calls the model's cell_indicator and face_indicator functions.
Definition ddg.hh:266
double average_time_reinit() const
Definition ddg.hh:1033
virtual void jacobian(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType alpha, NumberType beta, const VectorType &variables=VectorType()) override
Definition ddg.hh:645
QGauss< dim > quadrature
Definition ddg.hh:1070
double average_time_jacobian_assembly()
Definition ddg.hh:1053
uint num_reinits() const
Definition ddg.hh:1041
const FiniteElement< dim > & fe
Definition common.hh:339
Model_ Model
Definition ddg.hh:197
virtual const SparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition ddg.hh:254
QGauss< dim - 1 > quadrature_face
Definition ddg.hh:1071
std::vector< double > timings_jacobian
Definition ddg.hh:1081
SparseMatrix< NumberType > mass_matrix
Definition ddg.hh:1077
uint threads
Definition common.hh:343
uint num_jacobians() const
Definition ddg.hh:1061
typename Discretization::NumberType NumberType
Definition ddg.hh:198
static constexpr uint dim
Definition ddg.hh:202
std::vector< double > timings_residual
Definition ddg.hh:1080
const DoFHandler< dim > & dof_handler
Definition common.hh:340
Discretization & discretization
Definition common.hh:337
virtual const SparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition ddg.hh:258
virtual void jacobian_mass(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha, NumberType beta) override
Definition ddg.hh:584
void log(const std::string logger)
Definition ddg.hh:1021
typename Discretization::VectorType VectorType
Definition ddg.hh:199
uint batch_size
Definition common.hh:344
virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight) override
Definition ddg.hh:355
virtual void rebuild_jacobian_sparsity() override
Definition ddg.hh:242
std::vector< types::global_dof_index > extractor_dof_indices
Definition common.hh:355
std::vector< double > timings_reinit
Definition ddg.hh:1079
typename Discretization::Components Components
Definition ddg.hh:201
virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType weight_mass, const VectorType &variables=VectorType()) override
Definition ddg.hh:407
virtual void reinit_vector(VectorType &vec) const override
Definition ddg.hh:213
SparsityPattern sparsity_pattern_mass
Definition ddg.hh:1075
double average_time_residual_assembly()
Definition ddg.hh:1043
Assembler(Discretization &discretization, Model &model, const JSONValue &json)
Definition ddg.hh:204
auto i_tie(T &&...t)
Definition ddg.hh:23
auto fe_tie(T &&...t)
Definition ddg.hh:16
Definition complex_math.hh:10
unsigned int uint
Definition utils.hh:24
Definition tuples.hh:34
std::array< uint, 2 > cell_indices
Definition ddg.hh:177
std::array< double, 2 > values
Definition ddg.hh:178
std::vector< CopyFaceData_I > face_data
Definition ddg.hh:180
uint cell_index
Definition ddg.hh:182
double value
Definition ddg.hh:181
FullMatrix< NumberType > extractor_cell_jacobian
Definition ddg.hh:147
std::vector< types::global_dof_index > joint_dof_indices
Definition ddg.hh:148
FullMatrix< NumberType > cell_jacobian
Definition ddg.hh:146
void reinit(const FEInterfaceValues< dim > &fe_iv, uint n_extractors)
Definition ddg.hh:150
void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
Definition ddg.hh:164
FullMatrix< NumberType > extractor_cell_jacobian
Definition ddg.hh:160
std::vector< CopyDataFace_J > face_data
Definition ddg.hh:162
FullMatrix< NumberType > cell_jacobian
Definition ddg.hh:159
std::vector< types::global_dof_index > local_dof_indices
Definition ddg.hh:161
std::vector< types::global_dof_index > joint_dof_indices
Definition ddg.hh:121
void reinit(const FEInterfaceValues< dim > &fe_iv)
Definition ddg.hh:123
Vector< NumberType > cell_residual
Definition ddg.hh:120
Vector< NumberType > cell_residual
Definition ddg.hh:130
std::vector< CopyDataFace_R > face_data
Definition ddg.hh:132
void reinit(const Iterator &cell, uint dofs_per_cell)
Definition ddg.hh:134
std::vector< types::global_dof_index > local_dof_indices
Definition ddg.hh:131
Class to hold data for each assembly thread, i.e. FEValues for cells, interfaces, as well as pre-allo...
Definition ddg.hh:37
std::vector< Vector< NumberType > > solution_dot
Definition ddg.hh:107
std::vector< Tensor< 1, dim > > cached_shape_grads
Definition ddg.hh:114
std::vector< double > cached_shape_values
Definition ddg.hh:113
std::vector< Vector< NumberType > > solution
Definition ddg.hh:104
std::vector< uint > comp
Definition ddg.hh:112
typename Discretization::NumberType NumberType
Definition ddg.hh:39
std::vector< std::vector< Tensor< 1, dim, NumberType > > > solution_grad
Definition ddg.hh:105
const uint n_components
Definition ddg.hh:99
array< std::vector< std::vector< Tensor< 1, dim, NumberType > > >, 2 > solution_grad_interface
Definition ddg.hh:109
ScratchData(const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< dim > &quadrature, const Quadrature< dim - 1 > &quadrature_face, const UpdateFlags update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values|update_hessians, const UpdateFlags interface_update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values|update_normal_vectors|update_hessians)
Definition ddg.hh:41
ScratchData(const ScratchData< Discretization > &scratch_data)
Definition ddg.hh:71
array< std::vector< Vector< NumberType > >, 2 > solution_interface
Definition ddg.hh:108
FEValues< dim > fe_values
Definition ddg.hh:101
std::vector< Tensor< 2, dim > > cached_shape_hessians
Definition ddg.hh:115
FEInterfaceValues< dim > fe_interface_values
Definition ddg.hh:102
static constexpr uint dim
Definition ddg.hh:38
std::vector< std::vector< Tensor< 2, dim, NumberType > > > solution_hess
Definition ddg.hh:106
array< std::vector< std::vector< Tensor< 2, dim, NumberType > > >, 2 > solution_hess_interface
Definition ddg.hh:110
A class to store a tuple with elements that can be accessed by name. The names are stored as FixedStr...
Definition tuples.hh:56