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

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FEM/assembler/cg.hh Source File
DiFfRG
cg.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 CG
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;
40
41 ScratchData(const Mapping<dim> &mapping, const FiniteElement<dim> &fe,
42 const dealii::Quadrature<dim> &quadrature, const dealii::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 const uint n_dofs = fe.n_dofs_per_cell();
66 comp.resize(n_dofs);
67 cached_shape_values.resize(n_dofs);
68 cached_shape_grads.resize(n_dofs);
69 cached_shape_hessians.resize(n_dofs);
70 }
71
73 : n_components(scratch_data.fe_values.get_fe().n_components()),
74 fe_values(scratch_data.fe_values.get_mapping(), scratch_data.fe_values.get_fe(),
75 scratch_data.fe_values.get_quadrature(), scratch_data.fe_values.get_update_flags()),
76 fe_interface_values(scratch_data.fe_interface_values.get_mapping(),
77 scratch_data.fe_interface_values.get_fe(),
78 scratch_data.fe_interface_values.get_quadrature(),
79 scratch_data.fe_interface_values.get_update_flags())
80 {
81 const uint q_size = scratch_data.fe_values.get_quadrature().size();
82 const uint q_face_size = scratch_data.fe_interface_values.get_quadrature().size();
83
84 solution.resize(q_size, Vector<NumberType>(n_components));
85 solution_grad.resize(q_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
86 solution_hess.resize(q_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
87 solution_dot.resize(q_size, Vector<NumberType>(n_components));
88 solution_interface[0].resize(q_face_size, Vector<NumberType>(n_components));
89 solution_interface[1].resize(q_face_size, Vector<NumberType>(n_components));
90 solution_grad_interface[0].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
91 solution_grad_interface[1].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
92 solution_hess_interface[0].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
93 solution_hess_interface[1].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
94 const uint n_dofs_copy = scratch_data.comp.size();
95 comp.resize(n_dofs_copy);
96 cached_shape_values.resize(n_dofs_copy);
97 cached_shape_grads.resize(n_dofs_copy);
98 cached_shape_hessians.resize(n_dofs_copy);
99 }
100
102
103 FEValues<dim> fe_values;
104 FEInterfaceValues<dim> fe_interface_values;
105
106 std::vector<Vector<NumberType>> solution;
107 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad;
108 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess;
109 std::vector<Vector<NumberType>> solution_dot;
110 array<std::vector<Vector<NumberType>>, 2> solution_interface;
111 array<std::vector<std::vector<Tensor<1, dim, NumberType>>>, 2> solution_grad_interface;
112 array<std::vector<std::vector<Tensor<2, dim, NumberType>>>, 2> solution_hess_interface;
113
114 std::vector<uint> comp;
115
116 // Cached per-DoF shape function data for jacobian assembly
117 std::vector<double> cached_shape_values;
118 std::vector<Tensor<1, dim>> cached_shape_grads;
119 std::vector<Tensor<2, dim>> cached_shape_hessians;
120 };
121
122 template <typename NumberType> struct CopyData_R {
123 Vector<NumberType> cell_residual;
124 std::vector<types::global_dof_index> local_dof_indices;
125
126 template <class Iterator> void reinit(const Iterator &cell, uint dofs_per_cell)
127 {
128 cell_residual.reinit(dofs_per_cell);
129 local_dof_indices.resize(dofs_per_cell);
130 cell->get_dof_indices(local_dof_indices);
131 }
132 };
133
134 template <typename NumberType> struct CopyData_J {
135 FullMatrix<NumberType> cell_jacobian;
136 FullMatrix<NumberType> extractor_cell_jacobian;
137 std::vector<types::global_dof_index> local_dof_indices;
138
139 template <class Iterator> void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
140 {
141 cell_jacobian.reinit(dofs_per_cell, dofs_per_cell);
142 if (n_extractors > 0) extractor_cell_jacobian.reinit(dofs_per_cell, n_extractors);
143 local_dof_indices.resize(dofs_per_cell);
144 cell->get_dof_indices(local_dof_indices);
145 }
146 };
147
148 template <typename NumberType> struct CopyData_I {
150 std::array<uint, 2> cell_indices;
151 std::array<double, 2> values;
152 };
153 std::vector<CopyFaceData_I> face_data;
154 double value;
156 };
157 } // namespace internal
158
164 template <typename Discretization_, typename Model_> class Assembler : public FEMAssembler<Discretization_, Model_>
165 {
167
168 public:
169 using Discretization = Discretization_;
170 using Model = Model_;
173
175 static constexpr uint dim = Discretization::dim;
176
178 : Base(discretization, model, json),
179 quadrature(fe.degree + 1 + json.get_uint("/discretization/overintegration")),
180 quadrature_face(fe.degree + 1 + json.get_uint("/discretization/overintegration"))
181 {
182 static_assert(Components::count_fe_subsystems() == 1, "A CG model cannot have multiple submodels!");
183 reinit();
184 }
185
186 virtual void reinit_vector(VectorType &vec) const override { vec.reinit(dof_handler.n_dofs()); }
187
188 virtual void reinit() override
189 {
190 Timer timer;
191
192 Base::reinit();
193
194 // Mass sparsity pattern
195 {
196 DynamicSparsityPattern dsp(dof_handler.n_dofs());
197 DoFTools::make_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
198 /*keep_constrained_dofs = */ true);
199 sparsity_pattern_mass.copy_from(dsp);
201 MatrixCreator::create_mass_matrix(dof_handler, quadrature, mass_matrix, (Function<dim, NumberType> *)nullptr,
202 discretization.get_constraints());
203 }
204 // Jacobian sparsity pattern
205 {
206 DynamicSparsityPattern dsp(dof_handler.n_dofs());
207 DoFTools::make_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
208 /*keep_constrained_dofs = */ true);
209 sparsity_pattern_jacobian.copy_from(dsp);
210 }
211 timings_reinit.push_back(timer.wall_time());
212
213 // Boundary dofs
214 std::vector<IndexSet> component_boundary_dofs(Components::count_fe_functions());
215 for (uint c = 0; c < Components::count_fe_functions(); ++c) {
216 ComponentMask component_mask(Components::count_fe_functions(), false);
217 component_mask.set(c, true);
218 component_boundary_dofs[c] = DoFTools::extract_boundary_dofs(dof_handler, component_mask);
219 }
220 std::vector<std::vector<Point<dim>>> component_boundary_points(Components::count_fe_functions());
221 for (uint c = 0; c < Components::count_fe_functions(); ++c) {
222 component_boundary_points[c].resize(component_boundary_dofs[c].n_elements());
223 for (uint i = 0; i < component_boundary_dofs[c].n_elements(); ++i)
224 component_boundary_points[c][i] =
225 discretization.get_support_point(component_boundary_dofs[c].nth_index_in_set(i));
226 }
227
228 auto &constraints = discretization.get_constraints();
229 constraints.clear();
230 DoFTools::make_hanging_node_constraints(dof_handler, constraints);
231 model.affine_constraints(constraints, component_boundary_dofs, component_boundary_points);
232 constraints.close();
233 }
234
235 virtual void rebuild_jacobian_sparsity() override
236 {
237 // Jacobian sparsity pattern
238 DynamicSparsityPattern dsp(dof_handler.n_dofs());
239 DoFTools::make_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
240 /*keep_constrained_dofs = */ true);
241 for (uint row = 0; row < dsp.n_rows(); ++row)
242 for (const auto &col : extractor_dof_indices)
243 dsp.add(row, col);
244 sparsity_pattern_jacobian.copy_from(dsp);
245 }
246
247 virtual const SparsityPattern &get_sparsity_pattern_jacobian() const override
248 {
250 }
251 virtual const SparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
252
260 virtual void refinement_indicator(Vector<double> &indicator, const VectorType &solution_global) override
261 {
262 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
264 using CopyData = internal::CopyData_I<NumberType>;
265
266 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData &copy_data) {
267 scratch_data.fe_values.reinit(t_cell);
268 const auto &fe_v = scratch_data.fe_values;
269 copy_data.cell_index = t_cell->active_cell_index();
270 copy_data.value = 0;
271
272 const auto &JxW = fe_v.get_JxW_values();
273 const auto &q_points = fe_v.get_quadrature_points();
274 const auto &q_indices = fe_v.quadrature_point_indices();
275
276 auto &solution = scratch_data.solution;
277 auto &solution_grad = scratch_data.solution_grad;
278 auto &solution_hess = scratch_data.solution_hess;
279 fe_v.get_function_values(solution_global, solution);
280 fe_v.get_function_gradients(solution_global, solution_grad);
281 fe_v.get_function_hessians(solution_global, solution_hess);
282
283 double local_indicator = 0.;
284 for (const auto &q_index : q_indices) {
285 const auto &x_q = q_points[q_index];
286 model.cell_indicator(local_indicator, x_q,
287 i_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index]));
288 copy_data.value += JxW[q_index] * local_indicator;
289 }
290 };
291 const auto copier = [&](const CopyData &c) { indicator[c.cell_index] += c.value; };
292
293 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
294 CopyData copy_data;
295 MeshWorker::AssembleFlags assemble_flags = MeshWorker::assemble_own_cells;
296
297 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
298 copy_data, assemble_flags, nullptr, nullptr, threads, batch_size);
299 }
300
301 virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot,
302 NumberType weight) override
303 {
304 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
306 using CopyData = internal::CopyData_R<NumberType>;
307 const auto &constraints = discretization.get_constraints();
308
309 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
310 scratch_data.fe_values.reinit(cell);
311 const auto &fe_v = scratch_data.fe_values;
312 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
313
314 copy_data.reinit(cell, n_dofs);
315 const auto &JxW = fe_v.get_JxW_values();
316 const auto &q_points = fe_v.get_quadrature_points();
317 const auto &q_indices = fe_v.quadrature_point_indices();
318
319 auto &comp = scratch_data.comp;
320 for (uint i = 0; i < n_dofs; ++i)
321 comp[i] = fe_v.get_fe().system_to_component_index(i).first;
322
323 auto &solution = scratch_data.solution;
324 auto &solution_dot = scratch_data.solution_dot;
325 fe_v.get_function_values(solution_global, solution);
326 fe_v.get_function_values(solution_global_dot, solution_dot);
327
328 array<NumberType, Components::count_fe_functions()> mass{};
329 for (const auto &q_index : q_indices) {
330 const auto &x_q = q_points[q_index];
331 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
332
333 for (uint i = 0; i < n_dofs; ++i) {
334 copy_data.cell_residual(i) += weight * JxW[q_index] * fe_v.shape_value_component(i, q_index, comp[i]) *
335 mass[comp[i]]; // +phi_i(x_q) * mass(x_q, u_q)
336 }
337 }
338 };
339 const auto copier = [&](const CopyData &c) {
340 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, mass);
341 };
342
343 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
344 CopyData copy_data;
345 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
346
347 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
348 copy_data, flags, nullptr, nullptr, threads, batch_size);
349 }
350
351 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
352 const VectorType &solution_global_dot, NumberType weight_mass,
353 const VectorType &variables = VectorType()) override
354 {
355 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
357 using CopyData = internal::CopyData_R<NumberType>;
358 const auto &constraints = discretization.get_constraints();
359
360 // Find the EoM and extract whatever data is needed for the model.
361 std::array<NumberType, Components::count_extractors()> __extracted_data{{}};
362 if constexpr (Components::count_extractors() > 0)
363 this->extract(__extracted_data, solution_global, variables, true, false, true);
364 const auto &extracted_data = __extracted_data;
365
366 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
367 scratch_data.fe_values.reinit(cell);
368 const auto &fe_v = scratch_data.fe_values;
369 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
370
371 copy_data.reinit(cell, n_dofs);
372 const auto &JxW = fe_v.get_JxW_values();
373 const auto &q_points = fe_v.get_quadrature_points();
374 const auto &q_indices = fe_v.quadrature_point_indices();
375
376 auto &comp = scratch_data.comp;
377 for (uint i = 0; i < n_dofs; ++i)
378 comp[i] = fe.system_to_component_index(i).first;
379
380 auto &solution = scratch_data.solution;
381 auto &solution_grad = scratch_data.solution_grad;
382 auto &solution_hess = scratch_data.solution_hess;
383 auto &solution_dot = scratch_data.solution_dot;
384 fe_v.get_function_values(solution_global, solution);
385 fe_v.get_function_gradients(solution_global, solution_grad);
386 fe_v.get_function_hessians(solution_global, solution_hess);
387 fe_v.get_function_values(solution_global_dot, solution_dot);
388
389 array<Tensor<1, dim, NumberType>, Components::count_fe_functions()> flux{};
390 array<NumberType, Components::count_fe_functions()> source{};
391 array<NumberType, Components::count_fe_functions()> mass{};
392 for (const auto &q_index : q_indices) {
393 const auto &x_q = q_points[q_index];
394 model.flux(
395 flux, x_q,
396 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
397 model.source(
398 source, x_q,
399 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
400 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
401
402 for (uint i = 0; i < n_dofs; ++i) {
403 const auto &ci = comp[i];
404 copy_data.cell_residual(i) +=
405 JxW[q_index] * weight * // dx *
406 (-scalar_product(fe_v.shape_grad_component(i, q_index, ci),
407 flux[ci]) // -dphi_i(x_q) * flux(x_q, u_q)
408 + fe_v.shape_value_component(i, q_index, ci) * source[ci]); // -phi_i(x_q) * source(x_q, u_q)
409 copy_data.cell_residual(i) += weight_mass * JxW[q_index] * fe_v.shape_value_component(i, q_index, ci) *
410 mass[ci]; // +phi_i(x_q) * mass(x_q, u_q)
411 }
412 }
413 };
414 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
415 CopyData &copy_data) {
416 scratch_data.fe_interface_values.reinit(cell, face_no);
417 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
418 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
419
420 const auto &JxW = fe_fv.get_JxW_values();
421 const auto &q_points = fe_fv.get_quadrature_points();
422 const auto &q_indices = fe_fv.quadrature_point_indices();
423 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
424
425 auto &solution = scratch_data.solution_interface[0];
426 auto &solution_grad = scratch_data.solution_grad_interface[0];
427 auto &solution_hess = scratch_data.solution_hess_interface[0];
428 fe_fv.get_function_values(solution_global, solution);
429 fe_fv.get_function_gradients(solution_global, solution_grad);
430 fe_fv.get_function_hessians(solution_global, solution_hess);
431
432 auto &comp = scratch_data.comp;
433 for (uint i = 0; i < n_dofs; ++i)
434 comp[i] = fe.system_to_component_index(i).first;
435
436 array<Tensor<1, dim, NumberType>, Components::count_fe_functions()> numflux{};
437 for (const auto &q_index : q_indices) {
438 const auto &x_q = q_points[q_index];
439 model.boundary_numflux(
440 numflux, normals[q_index], x_q,
441 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
442
443 for (uint i = 0; i < n_dofs; ++i) {
444 const auto &ci = comp[i];
445 copy_data.cell_residual(i) +=
446 weight * JxW[q_index] * // dx
447 (fe_fv.shape_value_component(i, q_index, ci) *
448 scalar_product(numflux[ci], normals[q_index])); // phi_i(x_q) * numflux(x_q, u_q) * n(x_q)
449 }
450 }
451 };
452 const auto copier = [&](const CopyData &c) {
453 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, residual);
454 };
455
456 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
457 CopyData copy_data;
458 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces;
459
460 Timer timer;
461 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
462 copy_data, flags, boundary_worker, nullptr, threads, batch_size);
463 timings_residual.push_back(timer.wall_time());
464 }
465
466 virtual void jacobian_mass(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
467 const VectorType &solution_global_dot, NumberType alpha, NumberType beta) override
468 {
469 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
471 using CopyData = internal::CopyData_J<NumberType>;
472 const auto &constraints = discretization.get_constraints();
473
474 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
475 scratch_data.fe_values.reinit(cell);
476 const auto &fe_v = scratch_data.fe_values;
477 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
478
479 copy_data.reinit(cell, n_dofs, Components::count_extractors());
480 const auto &JxW = fe_v.get_JxW_values();
481 const auto &q_points = fe_v.get_quadrature_points();
482 const auto &q_indices = fe_v.quadrature_point_indices();
483
484 auto &comp = scratch_data.comp;
485 for (uint i = 0; i < n_dofs; ++i)
486 comp[i] = fe.system_to_component_index(i).first;
487
488 auto &solution = scratch_data.solution;
489 auto &solution_dot = scratch_data.solution_dot;
490
491 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass;
492 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass_dot;
493
494 fe_v.get_function_values(solution_global, solution);
495 fe_v.get_function_values(solution_global_dot, solution_dot);
496 for (const auto &q_index : q_indices) {
497 const auto &x_q = q_points[q_index];
498 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
499 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
500
501 for (uint i = 0; i < n_dofs; ++i) {
502 for (uint j = 0; j < n_dofs; ++j) {
503 copy_data.cell_jacobian(i, j) +=
504 JxW[q_index] * fe_v.shape_value_component(j, q_index, comp[j]) *
505 fe_v.shape_value_component(i, q_index, comp[i]) *
506 (alpha * j_mass_dot(comp[i], comp[j]) + beta * j_mass(comp[i], comp[j]));
507 }
508 }
509 }
510 };
511 const auto copier = [&](const CopyData &c) {
512 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
513 };
514
515 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
516 CopyData copy_data;
517 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
518
519 Timer timer;
520 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
521 copy_data, flags, nullptr, nullptr, threads, batch_size);
522 timings_jacobian.push_back(timer.wall_time());
523 }
524
525 virtual void jacobian(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global, NumberType weight,
526 const VectorType &solution_global_dot, NumberType alpha, NumberType beta,
527 const VectorType &variables = VectorType()) override
528 {
529 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
531 using CopyData = internal::CopyData_J<NumberType>;
532 const auto &constraints = discretization.get_constraints();
533
534 // Find the EoM and extract whatever data is needed for the model.
535 std::array<NumberType, Components::count_extractors()> extracted_data{{}};
536 if constexpr (Components::count_extractors() > 0) {
537 this->extract(extracted_data, solution_global, variables, true, true, true);
538 if (this->jacobian_extractors(this->extractor_jacobian, solution_global, variables))
540 }
541
542 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
543 scratch_data.fe_values.reinit(cell);
544 const auto &fe_v = scratch_data.fe_values;
545 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
546
547 copy_data.reinit(cell, n_dofs, Components::count_extractors());
548 const auto &JxW = fe_v.get_JxW_values();
549 const auto &q_points = fe_v.get_quadrature_points();
550 const auto &q_indices = fe_v.quadrature_point_indices();
551
552 auto &comp = scratch_data.comp;
553 for (uint i = 0; i < n_dofs; ++i)
554 comp[i] = fe.system_to_component_index(i).first;
555
556 auto &solution = scratch_data.solution;
557 auto &solution_dot = scratch_data.solution_dot;
558 auto &solution_grad = scratch_data.solution_grad;
559 auto &solution_hess = scratch_data.solution_hess;
560
561 fe_v.get_function_values(solution_global, solution);
562 fe_v.get_function_gradients(solution_global, solution_grad);
563 fe_v.get_function_hessians(solution_global, solution_hess);
564 fe_v.get_function_values(solution_global_dot, solution_dot);
565
566 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions()> j_flux;
567 SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()> j_grad_flux;
568 SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()> j_hess_flux;
569 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
570 j_extr_flux;
571 SimpleMatrix<NumberType, Components::count_fe_functions()> j_source;
572 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions()> j_grad_source;
573 SimpleMatrix<Tensor<2, dim, NumberType>, Components::count_fe_functions()> j_hess_source;
574 SimpleMatrix<NumberType, Components::count_fe_functions(), Components::count_extractors()> j_extr_source;
575 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass;
576 SimpleMatrix<NumberType, Components::count_fe_functions()> j_mass_dot;
577
578 for (const auto &q_index : q_indices) {
579 const auto &x_q = q_points[q_index];
580 model.template jacobian_flux_source<0, 0>(
581 j_flux, j_source, x_q,
582 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
583 model.template jacobian_flux_source_grad<1>(
584 j_grad_flux, j_grad_source, x_q,
585 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
586 model.template jacobian_flux_source_hess<2>(
587 j_hess_flux, j_hess_source, x_q,
588 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
589 if constexpr (Components::count_extractors() > 0) {
590 model.template jacobian_flux_source_extr<3>(
591 j_extr_flux, j_extr_source, x_q,
592 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
593 }
594 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
595 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
596
597 // Cache per-DoF shape function data for this quadrature point
598 auto &sv = scratch_data.cached_shape_values;
599 auto &sg = scratch_data.cached_shape_grads;
600 auto &sh = scratch_data.cached_shape_hessians;
601 for (uint k = 0; k < n_dofs; ++k) {
602 sv[k] = fe_v.shape_value_component(k, q_index, comp[k]);
603 sg[k] = fe_v.shape_grad_component(k, q_index, comp[k]);
604 sh[k] = fe_v.shape_hessian_component(k, q_index, comp[k]);
605 }
606
607 const auto do_work = [&](uint i_begin, uint i_end, uint j_begin, uint j_end) {
608 for (uint i = i_begin; i < i_end; ++i) {
609 const auto &ci = comp[i];
610 const auto &sv_i = sv[i];
611 const auto &sg_i = sg[i];
612 for (uint j = j_begin; j < j_end; ++j) {
613 const auto &cj = comp[j];
614 NumberType contribution = weight * JxW[q_index] *
615 (sv[j] * (-scalar_product(sg_i, j_flux(ci, cj)) + sv_i * j_source(ci, cj)) +
616 scalar_product(sg[j], -scalar_product(sg_i, j_grad_flux(ci, cj)) +
617 sv_i * j_grad_source(ci, cj)) +
618 scalar_product(sh[j], -scalar_product(sg_i, j_hess_flux(ci, cj)) +
619 sv_i * j_hess_source(ci, cj)));
620 contribution += JxW[q_index] * sv[j] * sv_i * (alpha * j_mass_dot(ci, cj) + beta * j_mass(ci, cj));
621 copy_data.cell_jacobian(i, j) += contribution;
622 }
623 }
624 };
625
626 if (n_dofs * n_dofs < 64)
627 do_work(0, n_dofs, 0, n_dofs);
628 else
629 tbb::parallel_for(
630 tbb::blocked_range2d<uint>(0, n_dofs, 0, n_dofs), [&](const tbb::blocked_range2d<uint> &range) {
631 do_work(range.rows().begin(), range.rows().end(), range.cols().begin(), range.cols().end());
632 });
633
634 // extractor contribution
635 if constexpr (Components::count_extractors() > 0) {
636 for (uint i = 0; i < n_dofs; ++i) {
637 for (uint e = 0; e < Components::count_extractors(); ++e)
638 copy_data.extractor_cell_jacobian(i, e) +=
639 weight * JxW[q_index] * // dx * phi_j * (
640 (-scalar_product(fe_v.shape_grad_component(i, q_index, comp[i]),
641 j_extr_flux(comp[i], e)) // -dphi_i * jflux
642 + fe_v.shape_value_component(i, q_index, comp[i]) *
643 j_extr_source(comp[i], e)); // -phi_i * jsource)
644 }
645 }
646 }
647 };
648
649 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
650 CopyData &copy_data) {
651 scratch_data.fe_interface_values.reinit(cell, face_no);
652 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
653 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
654
655 const auto &JxW = fe_fv.get_JxW_values();
656 const auto &q_points = fe_fv.get_quadrature_points();
657 const auto &q_indices = fe_fv.quadrature_point_indices();
658 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
659
660 auto &solution = scratch_data.solution_interface[0];
661 auto &solution_grad = scratch_data.solution_grad_interface[0];
662 auto &solution_hess = scratch_data.solution_hess;
663 fe_fv.get_function_values(solution_global, solution);
664 fe_fv.get_function_gradients(solution_global, solution_grad);
665 fe_fv.get_function_hessians(solution_global, solution_hess);
666
667 auto &comp = scratch_data.comp;
668 for (uint i = 0; i < n_dofs; ++i)
669 comp[i] = fe.system_to_component_index(i).first;
670
671 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions()> j_boundary_numflux;
672 SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()>
673 j_grad_boundary_numflux;
674 SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()>
675 j_hess_boundary_numflux;
676 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
677 j_extr_boundary_numflux;
678 for (const auto &q_index : q_indices) {
679 const auto &x_q = q_points[q_index];
680 model.template jacobian_boundary_numflux<0, 0>(
681 j_boundary_numflux, normals[q_index], x_q,
682 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
683 model.template jacobian_boundary_numflux_grad<1>(
684 j_grad_boundary_numflux, normals[q_index], x_q,
685 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
686 model.template jacobian_boundary_numflux_hess<2>(
687 j_hess_boundary_numflux, normals[q_index], x_q,
688 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
689 if constexpr (Components::count_extractors() > 0) {
690 model.template jacobian_boundary_numflux_extr<3>(
691 j_extr_boundary_numflux, normals[q_index], x_q,
692 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
693 }
694
695 // Cache per-DoF shape function data for boundary
696 auto &bsv = scratch_data.cached_shape_values;
697 auto &bsg = scratch_data.cached_shape_grads;
698 auto &bsh = scratch_data.cached_shape_hessians;
699 for (uint k = 0; k < n_dofs; ++k) {
700 bsv[k] = fe_fv.shape_value_component(k, q_index, comp[k]);
701 bsg[k] = fe_fv.shape_grad_component(k, q_index, comp[k]);
702 bsh[k] = fe_fv.shape_hessian_component(k, q_index, comp[k]);
703 }
704
705 const auto do_bnd_work = [&](uint i_begin, uint i_end, uint j_begin, uint j_end) {
706 for (uint i = i_begin; i < i_end; ++i) {
707 const auto &ci = comp[i];
708 const auto &sv_i = bsv[i];
709 for (uint j = j_begin; j < j_end; ++j) {
710 const auto &cj = comp[j];
711 const auto n_dot_jnf = scalar_product(j_boundary_numflux(ci, cj), normals[q_index]);
712 const auto n_dot_jgnf = scalar_product(j_grad_boundary_numflux(ci, cj), normals[q_index]);
713 const auto n_dot_jhnf = scalar_product(j_hess_boundary_numflux(ci, cj), normals[q_index]);
714 copy_data.cell_jacobian(i, j) +=
715 weight * JxW[q_index] *
716 (bsv[j] * sv_i * n_dot_jnf + scalar_product(bsg[j], sv_i * n_dot_jgnf) +
717 scalar_product(bsh[j], sv_i * n_dot_jhnf));
718 }
719 }
720 };
721
722 if (n_dofs * n_dofs < 64)
723 do_bnd_work(0, n_dofs, 0, n_dofs);
724 else
725 tbb::parallel_for(
726 tbb::blocked_range2d<uint>(0, n_dofs, 0, n_dofs), [&](const tbb::blocked_range2d<uint> &range) {
727 do_bnd_work(range.rows().begin(), range.rows().end(), range.cols().begin(), range.cols().end());
728 });
729
730 // extractor contribution
731 if constexpr (Components::count_extractors() > 0) {
732 for (uint i = 0; i < n_dofs; ++i) {
733 for (uint e = 0; e < Components::count_extractors(); ++e)
734 copy_data.extractor_cell_jacobian(i, e) +=
735 weight * JxW[q_index] * // dx * phi_j(x_q)
736 (fe_fv.shape_value_component(i, q_index, comp[i]) *
737 scalar_product(j_extr_boundary_numflux(comp[i], e),
738 normals[q_index])); // phi_i(x_q) * j_numflux(x_q, u_q) * n(x_q)
739 }
740 }
741 }
742 };
743
744 const auto copier = [&](const CopyData &c) {
745 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
746 if constexpr (Components::count_extractors() > 0) {
747 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices.size(), extractor_dof_indices.size());
748 c.extractor_cell_jacobian.mmult(extractor_dependence, this->extractor_jacobian);
749 constraints.distribute_local_to_global(extractor_dependence, c.local_dof_indices, extractor_dof_indices,
750 jacobian);
751 }
752 };
753
754 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
755 CopyData copy_data;
756 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces;
757
758 Timer timer;
759 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
760 copy_data, flags, boundary_worker, nullptr, threads, batch_size);
761 timings_jacobian.push_back(timer.wall_time());
762 }
763
764 void log(const std::string logger)
765 {
766 std::stringstream ss;
767 ss << "CG Assembler: " << std::endl;
768 ss << " Reinit: " << average_time_reinit() * 1000 << "ms (" << num_reinits() << ")" << std::endl;
769 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
770 << std::endl;
771 ss << " Jacobian: " << average_time_jacobian_assembly() * 1000 << "ms (" << num_jacobians() << ")"
772 << std::endl;
773 spdlog::get(logger)->info(ss.str());
774 }
775
776 double average_time_reinit() const
777 {
778 double t = 0.;
779 double n = timings_reinit.size();
780 for (const auto &t_ : timings_reinit)
781 t += t_ / n;
782 return t;
783 }
784 uint num_reinits() const { return timings_reinit.size(); }
785
787 {
788 double t = 0.;
789 double n = timings_residual.size();
790 for (const auto &t_ : timings_residual)
791 t += t_ / n;
792 return t;
793 }
794 uint num_residuals() const { return timings_residual.size(); }
795
797 {
798 double t = 0.;
799 double n = timings_jacobian.size();
800 for (const auto &t_ : timings_jacobian)
801 t += t_ / n;
802 return t;
803 }
804 uint num_jacobians() const { return timings_jacobian.size(); }
805
806 protected:
808 using Base::dof_handler;
809 using Base::fe;
810 using Base::mapping;
811 using Base::model;
812
813 QGauss<dim> quadrature;
814 QGauss<dim - 1> quadrature_face;
815 using Base::batch_size;
816 using Base::threads;
817
818 SparsityPattern sparsity_pattern_mass;
820 SparseMatrix<NumberType> mass_matrix;
821
822 std::vector<double> timings_reinit;
823 std::vector<double> timings_residual;
824 std::vector<double> timings_jacobian;
825
827 };
828 } // namespace CG
829} // namespace DiFfRG
The basic assembler that can be used for any standard CG scheme with flux and source.
Definition cg.hh:165
double average_time_residual_assembly()
Definition cg.hh:786
std::vector< double > timings_jacobian
Definition cg.hh:824
Model & model
Definition common.hh:338
double average_time_jacobian_assembly()
Definition cg.hh:796
const Mapping< dim > & mapping
Definition common.hh:341
uint num_reinits() const
Definition cg.hh:784
virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight) override
Definition cg.hh:301
const FiniteElement< dim > & fe
Definition common.hh:339
double average_time_reinit() const
Definition cg.hh:776
virtual void refinement_indicator(Vector< double > &indicator, const VectorType &solution_global) override
refinement indicator for adaptivity. Only calls the model's cell_indicator function,...
Definition cg.hh:260
static constexpr uint dim
Definition cg.hh:175
Model_ Model
Definition cg.hh:170
uint threads
Definition common.hh:343
typename Discretization::Components Components
Definition cg.hh:174
typename Discretization::NumberType NumberType
Definition cg.hh:171
std::vector< double > timings_residual
Definition cg.hh:823
virtual void rebuild_jacobian_sparsity() override
Definition cg.hh:235
std::vector< double > timings_reinit
Definition cg.hh:822
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition cg.hh:188
SparsityPattern sparsity_pattern_mass
Definition cg.hh:818
QGauss< dim > quadrature
Definition cg.hh:813
Assembler(Discretization &discretization, Model &model, const JSONValue &json)
Definition cg.hh:177
typename Discretization::VectorType VectorType
Definition cg.hh:172
const DoFHandler< dim > & dof_handler
Definition common.hh:340
virtual const SparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition cg.hh:247
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 cg.hh:351
SparsityPattern sparsity_pattern_jacobian
Definition cg.hh:819
QGauss< dim - 1 > quadrature_face
Definition cg.hh:814
virtual const SparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition cg.hh:251
Discretization & discretization
Definition common.hh:337
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 cg.hh:525
SparseMatrix< NumberType > mass_matrix
Definition cg.hh:820
uint batch_size
Definition common.hh:344
Discretization_ Discretization
Definition cg.hh:169
std::vector< types::global_dof_index > extractor_dof_indices
Definition common.hh:355
void log(const std::string logger)
Definition cg.hh:764
uint num_residuals() const
Definition cg.hh:794
virtual void jacobian_mass(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha, NumberType beta) override
Definition cg.hh:466
virtual void reinit_vector(VectorType &vec) const override
Definition cg.hh:186
uint num_jacobians() const
Definition cg.hh:804
NumberType_ NumberType
Definition cg.hh:34
Components_ Components
Definition cg.hh:33
Vector< NumberType > VectorType
Definition cg.hh:35
static constexpr uint dim
Definition cg.hh:38
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
A simple NxM-matrix class, which is used for cell-wise Jacobians.
Definition tuples.hh:156
auto fe_tie(T &&...t)
Definition cg.hh:16
auto i_tie(T &&...t)
Definition cg.hh:23
Definition complex_math.hh:10
unsigned int uint
Definition utils.hh:24
std::array< double, 2 > values
Definition cg.hh:151
std::array< uint, 2 > cell_indices
Definition cg.hh:150
double value
Definition cg.hh:154
uint cell_index
Definition cg.hh:155
std::vector< CopyFaceData_I > face_data
Definition cg.hh:153
FullMatrix< NumberType > extractor_cell_jacobian
Definition cg.hh:136
std::vector< types::global_dof_index > local_dof_indices
Definition cg.hh:137
FullMatrix< NumberType > cell_jacobian
Definition cg.hh:135
void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
Definition cg.hh:139
void reinit(const Iterator &cell, uint dofs_per_cell)
Definition cg.hh:126
Vector< NumberType > cell_residual
Definition cg.hh:123
std::vector< types::global_dof_index > local_dof_indices
Definition cg.hh:124
Class to hold data for each assembly thread, i.e. FEValues for cells, interfaces, as well as pre-allo...
Definition cg.hh:37
std::vector< Vector< NumberType > > solution
Definition cg.hh:106
std::vector< uint > comp
Definition cg.hh:114
array< std::vector< std::vector< Tensor< 1, dim, NumberType > > >, 2 > solution_grad_interface
Definition cg.hh:111
FEValues< dim > fe_values
Definition cg.hh:103
static constexpr uint dim
Definition cg.hh:38
array< std::vector< std::vector< Tensor< 2, dim, NumberType > > >, 2 > solution_hess_interface
Definition cg.hh:112
ScratchData(const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const dealii::Quadrature< dim > &quadrature, const dealii::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 cg.hh:41
std::vector< Tensor< 1, dim > > cached_shape_grads
Definition cg.hh:118
std::vector< std::vector< Tensor< 1, dim, NumberType > > > solution_grad
Definition cg.hh:107
std::vector< double > cached_shape_values
Definition cg.hh:117
FEInterfaceValues< dim > fe_interface_values
Definition cg.hh:104
ScratchData(const ScratchData< Discretization > &scratch_data)
Definition cg.hh:72
const uint n_components
Definition cg.hh:101
std::vector< Tensor< 2, dim > > cached_shape_hessians
Definition cg.hh:119
array< std::vector< Vector< NumberType > >, 2 > solution_interface
Definition cg.hh:110
std::vector< std::vector< Tensor< 2, dim, NumberType > > > solution_hess
Definition cg.hh:108
std::vector< Vector< NumberType > > solution_dot
Definition cg.hh:109
typename Discretization::NumberType NumberType
Definition cg.hh:39
Definition tuples.hh:34
A class to store a tuple with elements that can be accessed by name. The names are stored as FixedStr...
Definition tuples.hh:56