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

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FEM/assembler/dg.hh Source File
DiFfRG
dg.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 DG
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 &...>, StringSet<"fe_functions", "extractors", "variables">>(std::tie(t...));
19 }
20
21 template <typename... T> auto i_tie(T &&...t)
22 {
23 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "fe_derivatives", "fe_hessians">>(
24 std::tie(t...));
25 }
26
27 namespace internal
28 {
33 template <typename Discretization> struct ScratchData {
34 static constexpr int dim = Discretization::dim;
36 using VectorType = Vector<NumberType>;
37
38 ScratchData(const Mapping<dim> &mapping, const FiniteElement<dim> &fe, const Quadrature<dim> &quadrature,
39 const Quadrature<dim - 1> &quadrature_face,
40 const UpdateFlags update_flags = update_values | update_gradients | update_quadrature_points |
41 update_JxW_values,
42 const UpdateFlags interface_update_flags = update_values | update_gradients |
43 update_quadrature_points | update_JxW_values |
44 update_normal_vectors)
45 : n_components(fe.n_components()), fe_values(mapping, fe, quadrature, update_flags),
46 fe_interface_values(mapping, fe, quadrature_face, interface_update_flags)
47 {
48 solution.resize(quadrature.size(), VectorType(n_components));
49 solution_dot.resize(quadrature.size(), VectorType(n_components));
50 solution_interface[0].resize(quadrature_face.size(), VectorType(n_components));
51 solution_interface[1].resize(quadrature_face.size(), VectorType(n_components));
52 comp.resize(fe.n_dofs_per_cell());
53 for (uint i = 0; i < comp.size(); ++i)
54 comp[i] = fe.system_to_component_index(i).first;
55 }
56
58 : n_components(scratch_data.fe_values.get_fe().n_components()),
59 fe_values(scratch_data.fe_values.get_mapping(), scratch_data.fe_values.get_fe(),
60 scratch_data.fe_values.get_quadrature(), scratch_data.fe_values.get_update_flags()),
61 fe_interface_values(scratch_data.fe_interface_values.get_mapping(),
62 scratch_data.fe_interface_values.get_fe(),
63 scratch_data.fe_interface_values.get_quadrature(),
64 scratch_data.fe_interface_values.get_update_flags())
65 {
66 const auto &fe = scratch_data.fe_values.get_fe();
67 solution.resize(scratch_data.fe_values.get_quadrature().size(), VectorType(n_components));
68 solution_dot.resize(scratch_data.fe_values.get_quadrature().size(), VectorType(n_components));
69 solution_interface[0].resize(scratch_data.fe_interface_values.get_quadrature().size(),
71 solution_interface[1].resize(scratch_data.fe_interface_values.get_quadrature().size(),
73 comp.resize(fe.n_dofs_per_cell());
74 for (uint i = 0; i < comp.size(); ++i)
75 comp[i] = fe.system_to_component_index(i).first;
76 }
77
79
80 FEValues<dim> fe_values;
81 FEInterfaceValues<dim> fe_interface_values;
82
83 std::vector<VectorType> solution;
84 std::vector<VectorType> solution_dot;
85 array<std::vector<VectorType>, 2> solution_interface;
86 std::vector<uint> comp;
87 };
88
89 // TODO fewer memory allocations
90 template <typename NumberType> struct CopyData_R {
92 Vector<NumberType> cell_residual;
93 std::vector<types::global_dof_index> joint_dof_indices;
94
95 template <int dim> void reinit(const FEInterfaceValues<dim> &fe_iv)
96 {
97 cell_residual.reinit(fe_iv.n_current_interface_dofs());
98 joint_dof_indices = fe_iv.get_interface_dof_indices();
99 }
100 };
101
102 Vector<NumberType> cell_residual;
103 Vector<NumberType> cell_mass;
104 std::vector<types::global_dof_index> local_dof_indices;
105 std::vector<CopyDataFace_R> face_data;
106
107 template <class Iterator> void reinit(const Iterator &cell, uint dofs_per_cell)
108 {
109 cell_residual.reinit(dofs_per_cell);
110 cell_mass.reinit(dofs_per_cell);
111 local_dof_indices.resize(dofs_per_cell);
112 cell->get_dof_indices(local_dof_indices);
113 face_data.reserve(cell->n_faces());
114 }
115 };
116
117 // TODO fewer memory allocations
118 template <typename NumberType> struct CopyData_J {
120 FullMatrix<NumberType> cell_jacobian;
121 FullMatrix<NumberType> extractor_cell_jacobian;
122 std::vector<types::global_dof_index> joint_dof_indices;
123
124 template <int dim> void reinit(const FEInterfaceValues<dim> &fe_iv, uint n_extractors)
125 {
126 uint dofs_per_cell = fe_iv.n_current_interface_dofs();
127 cell_jacobian.reinit(dofs_per_cell, dofs_per_cell);
128 if (n_extractors > 0) extractor_cell_jacobian.reinit(dofs_per_cell, n_extractors);
129 joint_dof_indices = fe_iv.get_interface_dof_indices();
130 }
131 };
132
133 FullMatrix<NumberType> cell_jacobian;
134 FullMatrix<NumberType> extractor_cell_jacobian;
135 FullMatrix<NumberType> cell_mass_jacobian;
136 std::vector<types::global_dof_index> local_dof_indices;
137 std::vector<CopyDataFace_J> face_data;
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 cell_mass_jacobian.reinit(dofs_per_cell, dofs_per_cell);
144 local_dof_indices.resize(dofs_per_cell);
145 cell->get_dof_indices(local_dof_indices);
146 face_data.reserve(cell->n_faces());
147 }
148 };
149
150 template <typename NumberType> struct CopyData_I {
152 std::array<uint, 2> cell_indices;
153 std::array<double, 2> values;
154 };
155 std::vector<CopyFaceData_I> face_data;
156 double value;
158 };
159 } // namespace internal
160
166 template <typename Discretization_, typename Model_> class Assembler : public FEMAssembler<Discretization_, Model_>
167 {
169
170 public:
171 using Discretization = Discretization_;
172 using Model = Model_;
175
177 static constexpr uint n_components = Components::count_fe_functions(0);
178 static constexpr uint dim = Discretization::dim;
179
181 : Base(discretization, model, json),
182 quadrature(fe.degree + 1 + json.get_uint("/discretization/overintegration")),
183 quadrature_face(fe.degree + 1 + json.get_uint("/discretization/overintegration"))
184 {
185 static_assert(Components::count_fe_subsystems() == 1, "A DG model cannot have multiple submodels!");
186 reinit();
187 }
188
189 virtual void reinit_vector(VectorType &vec) const override { vec.reinit(dof_handler.n_dofs()); }
190
191 virtual void reinit() override
192 {
193 Timer timer;
194
195 Base::reinit();
196
197 // Mass sparsity pattern
198 {
199 DynamicSparsityPattern dsp(dof_handler.n_dofs());
200 DoFTools::make_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
201 /*keep_constrained_dofs = */ true);
202 sparsity_pattern_mass.copy_from(dsp);
204 MatrixCreator::create_mass_matrix(dof_handler, quadrature, mass_matrix, (Function<dim, NumberType> *)nullptr,
205 discretization.get_constraints());
206 }
207 // Jacobian sparsity pattern
208 {
209 DynamicSparsityPattern dsp(dof_handler.n_dofs());
210 DoFTools::make_flux_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
211 /*keep_constrained_dofs = */ true);
212 sparsity_pattern_jacobian.copy_from(dsp);
213 }
214
215 timings_reinit.push_back(timer.wall_time());
216 }
217
218 virtual void rebuild_jacobian_sparsity() override
219 {
220 // Jacobian sparsity pattern
221 DynamicSparsityPattern dsp(dof_handler.n_dofs());
222 DoFTools::make_flux_sparsity_pattern(dof_handler, dsp, discretization.get_constraints(),
223 /*keep_constrained_dofs = */ true);
224 for (uint row = 0; row < dsp.n_rows(); ++row)
225 for (const auto &col : extractor_dof_indices)
226 dsp.add(row, col);
227 sparsity_pattern_jacobian.copy_from(dsp);
228 }
229
230 virtual const SparsityPattern &get_sparsity_pattern_jacobian() const override
231 {
233 }
234 virtual const SparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
235
242 virtual void refinement_indicator(Vector<double> &indicator, const VectorType &solution_global) override
243 {
244 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
246 using CopyData = internal::CopyData_I<NumberType>;
247
248 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData &copy_data) {
249 scratch_data.fe_values.reinit(t_cell);
250 const auto &fe_v = scratch_data.fe_values;
251 copy_data.cell_index = t_cell->active_cell_index();
252 copy_data.value = 0;
253
254 const auto &JxW = fe_v.get_JxW_values();
255 const auto &q_points = fe_v.get_quadrature_points();
256 const auto &q_indices = fe_v.quadrature_point_indices();
257
258 auto &solution = scratch_data.solution;
259 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad;
260 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess;
261 solution_grad.resize(q_points.size(), std::vector<Tensor<1, dim, NumberType>>(n_components));
262 solution_hess.resize(q_points.size(), std::vector<Tensor<2, dim, NumberType>>(n_components));
263 fe_v.get_function_values(solution_global, solution);
264 fe_v.get_function_gradients(solution_global, solution_grad);
265 fe_v.get_function_hessians(solution_global, solution_hess);
266
267 double local_indicator = 0.;
268 for (const auto &q_index : q_indices) {
269 const auto &x_q = q_points[q_index];
270 model.cell_indicator(local_indicator, x_q,
271 i_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index]));
272 copy_data.value += JxW[q_index] * local_indicator;
273 }
274 };
275 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
276 const uint &nf, const unsigned int &nsf, Scratch &scratch_data,
277 CopyData &copy_data) {
278 scratch_data.fe_interface_values.reinit(t_cell, f, sf, t_ncell, nf, nsf);
279 const auto &fe_iv = scratch_data.fe_interface_values;
280 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
281 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
282
283 auto &copy_data_face = copy_data.face_data.emplace_back();
284 copy_data_face.cell_indices[0] = t_cell->active_cell_index();
285 copy_data_face.cell_indices[1] = t_ncell->active_cell_index();
286 copy_data_face.values[0] = 0;
287 copy_data_face.values[1] = 0;
288
289 const auto &JxW = fe_iv.get_JxW_values();
290 const auto &q_points = fe_iv.get_quadrature_points();
291 const auto &q_indices = fe_iv.quadrature_point_indices();
292 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
293
294 auto &solution_s = scratch_data.solution_interface[0];
295 auto &solution_n = scratch_data.solution_interface[1];
296 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad_s;
297 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess_s;
298 solution_grad_s.resize(q_points.size(), std::vector<Tensor<1, dim, NumberType>>(n_components));
299 solution_hess_s.resize(q_points.size(), std::vector<Tensor<2, dim, NumberType>>(n_components));
300 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad_n;
301 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess_n;
302 solution_grad_n.resize(q_points.size(), std::vector<Tensor<1, dim, NumberType>>(n_components));
303 solution_hess_n.resize(q_points.size(), std::vector<Tensor<2, dim, NumberType>>(n_components));
304 fe_iv_s.get_function_values(solution_global, solution_s);
305 fe_iv_n.get_function_values(solution_global, solution_n);
306 fe_iv_s.get_function_gradients(solution_global, solution_grad_s);
307 fe_iv_n.get_function_gradients(solution_global, solution_grad_n);
308 fe_iv_s.get_function_hessians(solution_global, solution_hess_s);
309 fe_iv_n.get_function_hessians(solution_global, solution_hess_n);
310
311 array<double, 2> local_indicator{};
312 for (const auto &q_index : q_indices) {
313 const auto &x_q = q_points[q_index];
314 model.face_indicator(local_indicator, normals[q_index], x_q,
315 i_tie(solution_s[q_index], solution_grad_s[q_index], solution_hess_s[q_index]),
316 i_tie(solution_n[q_index], solution_grad_n[q_index], solution_hess_n[q_index]));
317
318 copy_data_face.values[0] += JxW[q_index] * local_indicator[0] * (1. + t_cell->at_boundary());
319 copy_data_face.values[1] += JxW[q_index] * local_indicator[1] * (1. + t_ncell->at_boundary());
320 }
321 };
322 const auto copier = [&](const CopyData &c) {
323 for (auto &cdf : c.face_data)
324 for (uint j = 0; j < 2; ++j)
325 indicator[cdf.cell_indices[j]] += cdf.values[j];
326 indicator[c.cell_index] += c.value;
327 };
328 const UpdateFlags update_flags =
329 update_values | update_gradients | update_quadrature_points | update_JxW_values | update_hessians;
330 const UpdateFlags interface_update_flags = update_values | update_gradients | update_quadrature_points |
331 update_JxW_values | update_normal_vectors | update_hessians;
332
333 Scratch scratch_data(mapping, fe, quadrature, quadrature_face, update_flags, interface_update_flags);
334 CopyData copy_data;
335 MeshWorker::AssembleFlags assemble_flags = MeshWorker::assemble_own_cells;
336
337 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
338 copy_data, assemble_flags, nullptr, face_worker, threads, batch_size);
339 }
340
341 virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot,
342 NumberType weight) override
343 {
344 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
346 using CopyData = internal::CopyData_R<NumberType>;
347 const auto &constraints = discretization.get_constraints();
348
349 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
350 scratch_data.fe_values.reinit(cell);
351 const auto &fe_v = scratch_data.fe_values;
352 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
353
354 copy_data.reinit(cell, n_dofs);
355 const auto &JxW = fe_v.get_JxW_values();
356 const auto &q_points = fe_v.get_quadrature_points();
357 const auto &q_indices = fe_v.quadrature_point_indices();
358
359 auto &solution = scratch_data.solution;
360 auto &solution_dot = scratch_data.solution_dot;
361 fe_v.get_function_values(solution_global, solution);
362 fe_v.get_function_values(solution_global_dot, solution_dot);
363
364 const auto &comp = scratch_data.comp;
365
366 array<NumberType, n_components> mass{};
367 for (const auto &q_index : q_indices) {
368 const auto &x_q = q_points[q_index];
369 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
370
371 for (uint i = 0; i < n_dofs; ++i) {
372 const auto component_i = comp[i];
373 copy_data.cell_residual(i) += weight * JxW[q_index] *
374 fe_v.shape_value_component(i, q_index, component_i) *
375 mass[component_i]; // +phi_i(x_q) * mass(x_q, u_q)
376 }
377 }
378 };
379 const auto copier = [&](const CopyData &c) {
380 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, mass);
381 };
382
383 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
384 CopyData copy_data;
385 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
386
387 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
388 copy_data, flags, nullptr, nullptr, threads, batch_size);
389 }
390
391 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
392 const VectorType &solution_global_dot, NumberType weight_mass,
393 const VectorType &variables = VectorType()) override
394 {
395 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
397 using CopyData = internal::CopyData_R<NumberType>;
398 const auto &constraints = discretization.get_constraints();
399
400 // Find the EoM and extract whatever data is needed for the model.
401 std::array<NumberType, Components::count_extractors()> __extracted_data{{}};
402 if constexpr (Components::count_extractors() > 0)
403 this->extract(__extracted_data, solution_global, variables, true, false, true);
404 const auto &extracted_data = __extracted_data;
405
406 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
407 scratch_data.fe_values.reinit(cell);
408 const auto &fe_v = scratch_data.fe_values;
409 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
410
411 copy_data.reinit(cell, n_dofs);
412 const auto &JxW = fe_v.get_JxW_values();
413 const auto &q_points = fe_v.get_quadrature_points();
414 const auto &q_indices = fe_v.quadrature_point_indices();
415
416 auto &solution = scratch_data.solution;
417 auto &solution_dot = scratch_data.solution_dot;
418 fe_v.get_function_values(solution_global, solution);
419 fe_v.get_function_values(solution_global_dot, solution_dot);
420
421 const auto &comp = scratch_data.comp;
422
423 array<NumberType, n_components> mass{};
424 array<Tensor<1, dim, NumberType>, n_components> flux{};
425 array<NumberType, n_components> source{};
426 for (const auto &q_index : q_indices) {
427 const auto &x_q = q_points[q_index];
428 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
429 model.flux(flux, x_q, fe_tie(solution[q_index], extracted_data, variables));
430 model.source(source, x_q, fe_tie(solution[q_index], extracted_data, variables));
431
432 for (uint i = 0; i < n_dofs; ++i) {
433 const auto component_i = comp[i];
434 copy_data.cell_mass(i) += weight_mass * JxW[q_index] *
435 fe_v.shape_value_component(i, q_index, component_i) *
436 mass[component_i]; // +phi_i(x_q) * mass(x_q, u_q)
437 copy_data.cell_residual(i) += JxW[q_index] * weight * // dx *
438 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
439 flux[component_i]) // -dphi_i(x_q) * flux(x_q, u_q)
440 + fe_v.shape_value_component(i, q_index, component_i) *
441 source[component_i]); // -phi_i(x_q) * source(x_q, u_q)
442 }
443 }
444 };
445 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
446 CopyData &copy_data) {
447 scratch_data.fe_interface_values.reinit(cell, face_no);
448 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
449 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
450
451 const auto &JxW = fe_fv.get_JxW_values();
452 const auto &q_points = fe_fv.get_quadrature_points();
453 const auto &q_indices = fe_fv.quadrature_point_indices();
454 auto &solution = scratch_data.solution_interface[0];
455
456 array<Tensor<1, dim, NumberType>, n_components> numflux{};
457
458 fe_fv.get_function_values(solution_global, solution);
459 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
460
461 const auto &comp = scratch_data.comp;
462
463 for (const auto &q_index : q_indices) {
464 const auto &x_q = q_points[q_index];
465 model.boundary_numflux(numflux, normals[q_index], x_q,
466 fe_tie(solution[q_index], extracted_data, variables));
467
468 for (uint i = 0; i < n_dofs; ++i) {
469 const auto component_i = comp[i];
470 copy_data.cell_residual(i) +=
471 weight * JxW[q_index] * // dx
472 (fe_fv.shape_value_component(i, q_index, component_i) *
473 scalar_product(numflux[component_i], normals[q_index])); // phi_i(x_q) * numflux(x_q, u_q) * n(x_q)
474 }
475 }
476 };
477 const auto face_worker = [&](const Iterator &cell, const uint &f, const uint &sf, const Iterator &ncell,
478 const uint &nf, const uint &nsf, Scratch &scratch_data, CopyData &copy_data) {
479 scratch_data.fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
480 const auto &fe_iv = scratch_data.fe_interface_values;
481 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
482 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
483 const uint n_dofs = fe_iv.n_current_interface_dofs();
484
485 copy_data.face_data.emplace_back();
486 auto &copy_data_face = copy_data.face_data.back();
487 copy_data_face.reinit(fe_iv);
488
489 const auto &JxW = fe_iv.get_JxW_values();
490 const auto &q_points = fe_iv.get_quadrature_points();
491 const auto &q_indices = fe_iv.quadrature_point_indices();
492 auto &solution_s = scratch_data.solution_interface[0];
493 auto &solution_n = scratch_data.solution_interface[1];
494
495 array<Tensor<1, dim, NumberType>, n_components> numflux{};
496
497 fe_iv_s.get_function_values(solution_global, solution_s);
498 fe_iv_n.get_function_values(solution_global, solution_n);
499 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
500
501 // Pre-compute component indices for interface DoFs
502 const auto &comp = scratch_data.comp;
503 std::vector<uint> iface_comp(n_dofs);
504 for (uint i = 0; i < n_dofs; ++i) {
505 const auto &cd_i = fe_iv.interface_dof_to_dof_indices(i);
506 iface_comp[i] = cd_i[0] == numbers::invalid_unsigned_int ? comp[cd_i[1]] : comp[cd_i[0]];
507 }
508
509 for (const auto &q_index : q_indices) {
510 const auto &x_q = q_points[q_index];
511 model.numflux(numflux, normals[q_index], x_q, fe_tie(solution_s[q_index], extracted_data, variables),
512 fe_tie(solution_n[q_index], extracted_data, variables));
513
514 for (uint i = 0; i < n_dofs; ++i) {
515 const auto component_i = iface_comp[i];
516 copy_data_face.cell_residual(i) +=
517 weight * JxW[q_index] * // dx
518 (fe_iv.jump_in_shape_values(i, q_index, component_i) *
519 scalar_product(numflux[component_i],
520 normals[q_index])); // [[phi_i(x_q)]] * numflux(x_q, u_q) * n(x_q)
521 }
522 }
523 };
524 const auto copier = [&](const CopyData &c) {
525 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, residual);
526 constraints.distribute_local_to_global(c.cell_mass, c.local_dof_indices, residual);
527 for (auto &cdf : c.face_data)
528 constraints.distribute_local_to_global(cdf.cell_residual, cdf.joint_dof_indices, residual);
529 };
530
531 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
532 CopyData copy_data;
533 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
534 MeshWorker::assemble_own_interior_faces_once;
535
536 Timer timer;
537 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
538 copy_data, flags, boundary_worker, face_worker, threads, batch_size);
539 timings_residual.push_back(timer.wall_time());
540 }
541
542 virtual void jacobian_mass(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
543 const VectorType &solution_global_dot, NumberType alpha = 1.,
544 NumberType beta = 1.) override
545 {
546 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
548 using CopyData = internal::CopyData_J<NumberType>;
549 const auto &constraints = discretization.get_constraints();
550
551 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
552 scratch_data.fe_values.reinit(cell);
553 const auto &fe_v = scratch_data.fe_values;
554 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
555
556 copy_data.reinit(cell, n_dofs, Components::count_extractors());
557 const auto &JxW = fe_v.get_JxW_values();
558 const auto &q_points = fe_v.get_quadrature_points();
559 const auto &q_indices = fe_v.quadrature_point_indices();
560
561 auto &solution = scratch_data.solution;
562 auto &solution_dot = scratch_data.solution_dot;
563 fe_v.get_function_values(solution_global, solution);
564 fe_v.get_function_values(solution_global_dot, solution_dot);
565
566 const auto &comp = scratch_data.comp;
567
570 for (const auto &q_index : q_indices) {
571 const auto &x_q = q_points[q_index];
572 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
573 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
574
575 for (uint i = 0; i < n_dofs; ++i) {
576 const auto component_i = comp[i];
577 for (uint j = 0; j < n_dofs; ++j) {
578 const auto component_j = comp[j];
579 copy_data.cell_jacobian(i, j) +=
580 JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) *
581 fe_v.shape_value_component(i, q_index, component_i) *
582 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
583 }
584 }
585 }
586 };
587 const auto copier = [&](const CopyData &c) {
588 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
589 };
590
591 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
592 CopyData copy_data;
593 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
594
595 Timer timer;
596 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
597 copy_data, flags, nullptr, nullptr, threads, batch_size);
598 timings_jacobian.push_back(timer.wall_time());
599 }
600
601 virtual void jacobian(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global, NumberType weight,
602 const VectorType &solution_global_dot, NumberType alpha, NumberType beta,
603 const VectorType &variables = VectorType()) override
604 {
605 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
607 using CopyData = internal::CopyData_J<NumberType>;
608 const auto &constraints = discretization.get_constraints();
609
610 // Find the EoM and extract whatever data is needed for the model.
611 std::array<NumberType, Components::count_extractors()> extracted_data{{}};
612 if constexpr (Components::count_extractors() > 0) {
613 this->extract(extracted_data, solution_global, variables, true, true, true);
614 if (this->jacobian_extractors(this->extractor_jacobian, solution_global, variables))
616 }
617
618 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
619 scratch_data.fe_values.reinit(cell);
620 const auto &fe_v = scratch_data.fe_values;
621 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
622
623 copy_data.reinit(cell, n_dofs, Components::count_extractors());
624 const auto &JxW = fe_v.get_JxW_values();
625 const auto &q_points = fe_v.get_quadrature_points();
626 const auto &q_indices = fe_v.quadrature_point_indices();
627
628 auto &solution = scratch_data.solution;
629 auto &solution_dot = scratch_data.solution_dot;
630 fe_v.get_function_values(solution_global, solution);
631 fe_v.get_function_values(solution_global_dot, solution_dot);
632
633 const auto &comp = scratch_data.comp;
634
638 SimpleMatrix<Tensor<1, dim>, n_components, Components::count_extractors()> j_extr_flux;
640 SimpleMatrix<NumberType, n_components, Components::count_extractors()> j_extr_source;
641 for (const auto &q_index : q_indices) {
642 const auto &x_q = q_points[q_index];
643 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
644 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
645 model.template jacobian_flux_source<0, 0>(j_flux, j_source, x_q,
646 fe_tie(solution[q_index], extracted_data, variables));
647 if constexpr (Components::count_extractors() > 0) {
648 model.template jacobian_flux_source_extr<1>(j_extr_flux, j_extr_source, x_q,
649 fe_tie(solution[q_index], extracted_data, variables));
650 }
651
652 for (uint i = 0; i < n_dofs; ++i) {
653 const auto component_i = comp[i];
654 for (uint j = 0; j < n_dofs; ++j) {
655 const auto component_j = comp[j];
656 copy_data.cell_jacobian(i, j) += weight * JxW[q_index] *
657 fe_v.shape_value_component(j, q_index, component_j) * // dx * phi_j * (
658 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
659 j_flux(component_i, component_j)) // -dphi_i * jflux
660 + fe_v.shape_value_component(i, q_index, component_i) *
661 j_source(component_i, component_j)); // -phi_i * jsource)
662 copy_data.cell_mass_jacobian(i, j) +=
663 JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) *
664 fe_v.shape_value_component(i, q_index, component_i) *
665 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
666 }
667 // extractor contribution
668 if constexpr (Components::count_extractors() > 0)
669 for (uint e = 0; e < Components::count_extractors(); ++e)
670 copy_data.extractor_cell_jacobian(i, e) +=
671 weight * JxW[q_index] * // dx * phi_j * (
672 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
673 j_extr_flux(component_i, e)) // -dphi_i * jflux
674 + fe_v.shape_value_component(i, q_index, component_i) *
675 j_extr_source(component_i, e)); // -phi_i * jsource)
676 }
677 }
678 };
679 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
680 CopyData &copy_data) {
681 scratch_data.fe_interface_values.reinit(cell, face_no);
682 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
683 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
684
685 const auto &JxW = fe_fv.get_JxW_values();
686 const auto &q_points = fe_fv.get_quadrature_points();
687 const auto &q_indices = fe_fv.quadrature_point_indices();
688 auto &solution = scratch_data.solution_interface[0];
689
690 SimpleMatrix<Tensor<1, dim>, n_components> j_boundary_numflux;
691 SimpleMatrix<Tensor<1, dim>, n_components, Components::count_extractors()> j_extr_boundary_numflux;
692
693 fe_fv.get_function_values(solution_global, solution);
694 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
695
696 const auto &comp = scratch_data.comp;
697
698 for (const auto &q_index : q_indices) {
699 const auto &x_q = q_points[q_index];
700 model.template jacobian_boundary_numflux<0, 0>(j_boundary_numflux, normals[q_index], x_q,
701 fe_tie(solution[q_index], extracted_data, variables));
702 if constexpr (Components::count_extractors() > 0)
703 model.template jacobian_boundary_numflux_extr<1>(j_extr_boundary_numflux, normals[q_index], x_q,
704 fe_tie(solution[q_index], extracted_data, variables));
705
706 for (uint i = 0; i < n_dofs; ++i) {
707 const auto component_i = comp[i];
708 for (uint j = 0; j < n_dofs; ++j) {
709 const auto component_j = comp[j];
710 copy_data.cell_jacobian(i, j) +=
711 weight * JxW[q_index] * fe_fv.shape_value_component(j, q_index, component_j) * // dx * phi_j(x_q)
712 (fe_fv.shape_value_component(i, q_index, component_i) *
713 scalar_product(j_boundary_numflux(component_i, component_j),
714 normals[q_index])); // phi_i(x_q) * j_numflux(x_q, u_q) * n(x_q)
715 }
716 // extractor contribution
717 if constexpr (Components::count_extractors() > 0)
718 for (uint e = 0; e < Components::count_extractors(); ++e)
719 copy_data.extractor_cell_jacobian(i, e) +=
720 weight * JxW[q_index] * // dx * phi_j(x_q)
721 (fe_fv.shape_value_component(i, q_index, component_i) *
722 scalar_product(j_extr_boundary_numflux(component_i, e),
723 normals[q_index])); // phi_i(x_q) * j_numflux(x_q, u_q) * n(x_q)
724 }
725 }
726 };
727 const auto face_worker = [&](const Iterator &cell, const uint &f, const uint &sf, const Iterator &ncell,
728 const uint &nf, const uint &nsf, Scratch &scratch_data, CopyData &copy_data) {
729 scratch_data.fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
730 const auto &fe_iv = scratch_data.fe_interface_values;
731 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
732 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
733 const uint n_dofs = fe_iv.n_current_interface_dofs();
734
735 copy_data.face_data.emplace_back();
736 auto &copy_data_face = copy_data.face_data.back();
737 copy_data_face.reinit(fe_iv, Components::count_extractors());
738
739 const auto &JxW = fe_iv.get_JxW_values();
740 const auto &q_points = fe_iv.get_quadrature_points();
741 const auto &q_indices = fe_iv.quadrature_point_indices();
742 auto &solution_s = scratch_data.solution_interface[0];
743 auto &solution_n = scratch_data.solution_interface[1];
744
745 fe_iv_s.get_function_values(solution_global, solution_s);
746 fe_iv_n.get_function_values(solution_global, solution_n);
747 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
748
749 // Pre-compute component indices and face numbers for interface DoFs
750 const auto &comp = scratch_data.comp;
751 std::vector<uint> iface_comp(n_dofs);
752 std::vector<uint> iface_face_no(n_dofs);
753 std::vector<uint> iface_local_dof(n_dofs);
754 for (uint i = 0; i < n_dofs; ++i) {
755 const auto &cd_i = fe_iv.interface_dof_to_dof_indices(i);
756 const uint fno = cd_i[0] == numbers::invalid_unsigned_int ? 1 : 0;
757 iface_face_no[i] = fno;
758 iface_local_dof[i] = cd_i[fno];
759 iface_comp[i] = comp[cd_i[fno]];
760 }
761
762 array<SimpleMatrix<Tensor<1, dim>, n_components>, 2> j_numflux;
763 array<SimpleMatrix<Tensor<1, dim>, n_components, Components::count_extractors()>, 2> j_extr_numflux;
764 for (const auto &q_index : q_indices) {
765 const auto &x_q = q_points[q_index];
766 model.template jacobian_numflux<0, 0>(j_numflux, normals[q_index], x_q,
767 fe_tie(solution_s[q_index], extracted_data, variables),
768 fe_tie(solution_n[q_index], extracted_data, variables));
769 if constexpr (Components::count_extractors() > 0)
770 model.template jacobian_numflux_extr<1>(j_extr_numflux, normals[q_index], x_q,
771 fe_tie(solution_s[q_index], extracted_data, variables),
772 fe_tie(solution_n[q_index], extracted_data, variables));
773
774 for (uint i = 0; i < n_dofs; ++i) {
775 const auto component_i = iface_comp[i];
776 const uint face_no_i = iface_face_no[i];
777 for (uint j = 0; j < n_dofs; ++j) {
778 const auto component_j = iface_comp[j];
779 const uint face_no_j = iface_face_no[j];
780
781 copy_data_face.cell_jacobian(i, j) +=
782 weight * JxW[q_index] *
783 fe_iv.get_fe_face_values(face_no_j).shape_value_component(iface_local_dof[j], q_index,
784 component_j) * // dx * phi_j(x_q)
785 (fe_iv.jump_in_shape_values(i, q_index, component_i) *
786 scalar_product(j_numflux[face_no_j](component_i, component_j),
787 normals[q_index])); // [[phi_i(x_q)]] * j_numflux(x_q, u_q)
788 }
789 // extractor contribution
790 if constexpr (Components::count_extractors() > 0)
791 for (uint e = 0; e < Components::count_extractors(); ++e)
792 copy_data_face.extractor_cell_jacobian(i, e) +=
793 weight * JxW[q_index] * // dx * phi_j(x_q)
794 (fe_iv.jump_in_shape_values(i, q_index, component_i) *
795 scalar_product(j_extr_numflux[face_no_i](component_i, e),
796 normals[q_index])); // [[phi_i(x_q)]] * j_numflux(x_q, u_q)
797 }
798 }
799 };
800 const auto copier = [&](const CopyData &c) {
801 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
802 constraints.distribute_local_to_global(c.cell_mass_jacobian, c.local_dof_indices, jacobian);
803 for (auto &cdf : c.face_data) {
804 constraints.distribute_local_to_global(cdf.cell_jacobian, cdf.joint_dof_indices, jacobian);
805 if constexpr (Components::count_extractors() > 0) {
806 FullMatrix<NumberType> extractor_dependence(cdf.joint_dof_indices.size(), extractor_dof_indices.size());
807 cdf.extractor_cell_jacobian.mmult(extractor_dependence, this->extractor_jacobian);
808 constraints.distribute_local_to_global(extractor_dependence, cdf.joint_dof_indices, extractor_dof_indices,
809 jacobian);
810 }
811 }
812 if constexpr (Components::count_extractors() > 0) {
813 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices.size(), extractor_dof_indices.size());
814 c.extractor_cell_jacobian.mmult(extractor_dependence, this->extractor_jacobian);
815 constraints.distribute_local_to_global(extractor_dependence, c.local_dof_indices, extractor_dof_indices,
816 jacobian);
817 }
818 };
819
820 Scratch scratch_data(mapping, fe, quadrature, quadrature_face);
821 CopyData copy_data;
822 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
823 MeshWorker::assemble_own_interior_faces_once;
824
825 Timer timer;
826 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
827 copy_data, flags, boundary_worker, face_worker, threads, batch_size);
828 timings_jacobian.push_back(timer.wall_time());
829 }
830
831 void log(const std::string logger)
832 {
833 std::stringstream ss;
834 ss << "DG Assembler: " << std::endl;
835 ss << " Reinit: " << average_time_reinit() * 1000 << "ms (" << num_reinits() << ")" << std::endl;
836 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
837 << std::endl;
838 ss << " Jacobian: " << average_time_jacobian_assembly() * 1000 << "ms (" << num_jacobians() << ")"
839 << std::endl;
840 spdlog::get(logger)->info(ss.str());
841 }
842
843 double average_time_reinit() const
844 {
845 double t = 0.;
846 double n = timings_reinit.size();
847 for (const auto &t_ : timings_reinit)
848 t += t_ / n;
849 return t;
850 }
851 uint num_reinits() const { return timings_reinit.size(); }
852
854 {
855 double t = 0.;
856 double n = timings_residual.size();
857 for (const auto &t_ : timings_residual)
858 t += t_ / n;
859 return t;
860 }
861 uint num_residuals() const { return timings_residual.size(); }
862
864 {
865 double t = 0.;
866 double n = timings_jacobian.size();
867 for (const auto &t_ : timings_jacobian)
868 t += t_ / n;
869 return t;
870 }
871 uint num_jacobians() const { return timings_jacobian.size(); }
872
873 protected:
875 using Base::dof_handler;
876 using Base::fe;
877 using Base::mapping;
878 using Base::model;
879
880 QGauss<dim> quadrature;
881 QGauss<dim - 1> quadrature_face;
882 using Base::batch_size;
883 using Base::threads;
884
885 SparsityPattern sparsity_pattern_mass;
887 SparseMatrix<NumberType> mass_matrix;
888
889 std::vector<double> timings_reinit;
890 std::vector<double> timings_residual;
891 std::vector<double> timings_jacobian;
892
894 };
895 } // namespace DG
896} // namespace DiFfRG
The basic assembler that can be used for any standard DG scheme with flux and source.
Definition dg.hh:167
static constexpr uint n_components
Definition dg.hh:177
Model & model
Definition common.hh:338
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 dg.hh:601
SparsityPattern sparsity_pattern_mass
Definition dg.hh:885
const Mapping< dim > & mapping
Definition common.hh:341
virtual const SparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition dg.hh:230
void log(const std::string logger)
Definition dg.hh:831
SparsityPattern sparsity_pattern_jacobian
Definition dg.hh:886
const FiniteElement< dim > & fe
Definition common.hh:339
typename Discretization::NumberType NumberType
Definition dg.hh:173
double average_time_residual_assembly()
Definition dg.hh:853
QGauss< dim - 1 > quadrature_face
Definition dg.hh:881
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition dg.hh:191
std::vector< double > timings_reinit
Definition dg.hh:889
typename Discretization::Components Components
Definition dg.hh:176
SparseMatrix< NumberType > mass_matrix
Definition dg.hh:887
uint num_residuals() const
Definition dg.hh:861
virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight) override
Definition dg.hh:341
typename Discretization::VectorType VectorType
Definition dg.hh:174
uint threads
Definition common.hh:343
uint num_jacobians() const
Definition dg.hh:871
static constexpr uint dim
Definition dg.hh:178
Assembler(Discretization &discretization, Model &model, const JSONValue &json)
Definition dg.hh:180
Discretization_ Discretization
Definition dg.hh:171
uint num_reinits() const
Definition dg.hh:851
virtual void jacobian_mass(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha=1., NumberType beta=1.) override
Definition dg.hh:542
Model_ Model
Definition dg.hh:172
const DoFHandler< dim > & dof_handler
Definition common.hh:340
double average_time_reinit() const
Definition dg.hh:843
Discretization & discretization
Definition common.hh:337
virtual void rebuild_jacobian_sparsity() override
Definition dg.hh:218
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 dg.hh:242
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 dg.hh:391
uint batch_size
Definition common.hh:344
std::vector< types::global_dof_index > extractor_dof_indices
Definition common.hh:355
virtual const SparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition dg.hh:234
virtual void reinit_vector(VectorType &vec) const override
Definition dg.hh:189
std::vector< double > timings_jacobian
Definition dg.hh:891
double average_time_jacobian_assembly()
Definition dg.hh:863
std::vector< double > timings_residual
Definition dg.hh:890
QGauss< dim > quadrature
Definition dg.hh:880
static constexpr uint dim
Definition dg.hh:39
Vector< NumberType > VectorType
Definition dg.hh:36
NumberType_ NumberType
Definition dg.hh:35
Components_ Components
Definition dg.hh:34
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
auto fe_tie(T &&...t)
Definition dg.hh:16
auto i_tie(T &&...t)
Definition dg.hh:21
Definition complex_math.hh:10
unsigned int uint
Definition utils.hh:24
std::array< uint, 2 > cell_indices
Definition dg.hh:152
std::array< double, 2 > values
Definition dg.hh:153
std::vector< CopyFaceData_I > face_data
Definition dg.hh:155
uint cell_index
Definition dg.hh:157
double value
Definition dg.hh:156
void reinit(const FEInterfaceValues< dim > &fe_iv, uint n_extractors)
Definition dg.hh:124
FullMatrix< NumberType > extractor_cell_jacobian
Definition dg.hh:121
FullMatrix< NumberType > cell_jacobian
Definition dg.hh:120
std::vector< types::global_dof_index > joint_dof_indices
Definition dg.hh:122
FullMatrix< NumberType > cell_mass_jacobian
Definition dg.hh:135
std::vector< types::global_dof_index > local_dof_indices
Definition dg.hh:136
FullMatrix< NumberType > extractor_cell_jacobian
Definition dg.hh:134
FullMatrix< NumberType > cell_jacobian
Definition dg.hh:133
void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
Definition dg.hh:139
std::vector< CopyDataFace_J > face_data
Definition dg.hh:137
Vector< NumberType > cell_residual
Definition dg.hh:92
std::vector< types::global_dof_index > joint_dof_indices
Definition dg.hh:93
void reinit(const FEInterfaceValues< dim > &fe_iv)
Definition dg.hh:95
std::vector< CopyDataFace_R > face_data
Definition dg.hh:105
Vector< NumberType > cell_mass
Definition dg.hh:103
void reinit(const Iterator &cell, uint dofs_per_cell)
Definition dg.hh:107
std::vector< types::global_dof_index > local_dof_indices
Definition dg.hh:104
Vector< NumberType > cell_residual
Definition dg.hh:102
Class to hold data for each assembly thread, i.e. FEValues for cells, interfaces, as well as pre-allo...
Definition dg.hh:33
array< std::vector< VectorType >, 2 > solution_interface
Definition dg.hh:85
typename Discretization::NumberType NumberType
Definition dg.hh:35
FEInterfaceValues< dim > fe_interface_values
Definition dg.hh:81
Vector< NumberType > VectorType
Definition dg.hh:36
std::vector< VectorType > solution
Definition dg.hh:83
ScratchData(const ScratchData< Discretization > &scratch_data)
Definition dg.hh:57
std::vector< VectorType > solution_dot
Definition dg.hh:84
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, const UpdateFlags interface_update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values|update_normal_vectors)
Definition dg.hh:38
FEValues< dim > fe_values
Definition dg.hh:80
static constexpr int dim
Definition dg.hh:34
const uint n_components
Definition dg.hh:78
std::vector< uint > comp
Definition dg.hh:86
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