/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FV/assembler/KurganovTadmor.hh Source File#

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FV/assembler/KurganovTadmor.hh Source File
DiFfRG
KurganovTadmor.hh
Go to the documentation of this file.
1#pragma once
2
3// external libraries
4
5// DiFfRG
6#include <deal.II/base/multithread_info.h>
7#include <deal.II/base/quadrature_lib.h>
8#include <deal.II/base/timer.h>
9#include <deal.II/dofs/dof_handler.h>
10#include <deal.II/dofs/dof_tools.h>
11#include <deal.II/fe/fe_interface_values.h>
12#include <deal.II/fe/fe_values.h>
13#include <deal.II/lac/full_matrix.h>
14#include <deal.II/lac/sparse_matrix.h>
15#include <deal.II/lac/vector.h>
16#include <deal.II/meshworker/mesh_loop.h>
17#include <deal.II/numerics/fe_field_function.h>
18#include <deal.II/numerics/matrix_tools.h>
19#include <deal.II/numerics/vector_tools.h>
20#include <spdlog/spdlog.h>
21#include <tbb/tbb.h>
22
25
27
28namespace DiFfRG
29{
30 namespace FV
31 {
32 namespace KurganovTadmor
33 {
34 using namespace dealii;
35
36 namespace internal
37 {
42 template <typename Discretization> struct ScratchData {
43 static constexpr int dim = Discretization::dim;
45 using VectorType = Vector<NumberType>;
46
47 ScratchData(const Mapping<dim> &mapping, const FiniteElement<dim> &fe, const Quadrature<dim> &quadrature,
48 const UpdateFlags update_flags = update_values | update_gradients | update_quadrature_points |
49 update_JxW_values)
50 : n_components(fe.n_components()), fe_values(mapping, fe, quadrature, update_flags)
51 {
52 solution.resize(quadrature.size(), VectorType(n_components));
53 solution_dot.resize(quadrature.size(), VectorType(n_components));
54 }
55
57 : n_components(scratch_data.fe_values.get_fe().n_components()),
58 fe_values(scratch_data.fe_values.get_mapping(), scratch_data.fe_values.get_fe(),
59 scratch_data.fe_values.get_quadrature(), scratch_data.fe_values.get_update_flags())
60 {
61 solution.resize(scratch_data.fe_values.get_quadrature().size(), VectorType(n_components));
62 solution_dot.resize(scratch_data.fe_values.get_quadrature().size(), VectorType(n_components));
63 }
64
66
67 FEValues<dim> fe_values;
68
69 std::vector<VectorType> solution;
70 std::vector<VectorType> solution_dot;
71 };
72
73 // TODO fewer memory allocations
74 template <typename NumberType> struct CopyData_R {
76 Vector<NumberType> cell_residual;
77 std::vector<types::global_dof_index> joint_dof_indices;
78
79 template <int dim> void reinit(const FEInterfaceValues<dim> &fe_iv)
80 {
81 cell_residual.reinit(fe_iv.n_current_interface_dofs());
82 joint_dof_indices = fe_iv.get_interface_dof_indices();
83 }
84 };
85
86 Vector<NumberType> cell_residual;
87 Vector<NumberType> cell_mass;
88 std::vector<types::global_dof_index> local_dof_indices;
89 std::vector<CopyDataFace_R> face_data;
90
91 template <class Iterator> void reinit(const Iterator &cell, uint dofs_per_cell)
92 {
93 cell_residual.reinit(dofs_per_cell);
94 cell_mass.reinit(dofs_per_cell);
95 local_dof_indices.resize(dofs_per_cell);
96 cell->get_dof_indices(local_dof_indices);
97 }
98 };
99
100 // TODO fewer memory allocations
101 template <typename NumberType> struct CopyData_J {
103 FullMatrix<NumberType> cell_jacobian;
104 FullMatrix<NumberType> extractor_cell_jacobian;
105 std::vector<types::global_dof_index> joint_dof_indices;
106
107 template <int dim> void reinit(const FEInterfaceValues<dim> &fe_iv, uint n_extractors)
108 {
109 uint dofs_per_cell = fe_iv.n_current_interface_dofs();
110 cell_jacobian.reinit(dofs_per_cell, dofs_per_cell);
111 if (n_extractors > 0) extractor_cell_jacobian.reinit(dofs_per_cell, n_extractors);
112 joint_dof_indices = fe_iv.get_interface_dof_indices();
113 }
114 };
115
116 FullMatrix<NumberType> cell_jacobian;
117 FullMatrix<NumberType> extractor_cell_jacobian;
118 FullMatrix<NumberType> cell_mass_jacobian;
119 std::vector<types::global_dof_index> local_dof_indices;
120 std::vector<CopyDataFace_J> face_data;
121
122 template <class Iterator> void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
123 {
124 cell_jacobian.reinit(dofs_per_cell, dofs_per_cell);
125 if (n_extractors > 0) extractor_cell_jacobian.reinit(dofs_per_cell, n_extractors);
126 cell_mass_jacobian.reinit(dofs_per_cell, dofs_per_cell);
127 local_dof_indices.resize(dofs_per_cell);
128 cell->get_dof_indices(local_dof_indices);
129 }
130 };
131
132 template <typename NumberType> struct CopyData_I {
134 std::array<uint, 2> cell_indices;
135 std::array<double, 2> values;
136 };
137 std::vector<CopyFaceData_I> face_data;
138 double value;
140 };
141 } // namespace internal
142
143 template <typename Discretization_, typename Model_>
145 class Assembler : public AbstractAssembler<typename Discretization_::VectorType,
146 typename Discretization_::SparseMatrixType, Discretization_::dim>
147 {
148 protected:
149 template <typename... T> auto fv_tie(T &&...t)
150 {
151 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions">>(std::tie(t...));
152 }
153
154 template <typename... T> static constexpr auto v_tie(T &&...t)
155 {
156 return named_tuple<std::tuple<T &...>, StringSet<"variables", "extractors">>(std::tie(t...));
157 }
158
159 template <typename... T> static constexpr auto e_tie(T &&...t)
160 {
161 return named_tuple<std::tuple<T &...>,
162 StringSet<"fe_functions", "fe_derivatives", "fe_hessians", "extractors", "variables">>(
163 std::tie(t...));
164 }
165
166 public:
167 using Discretization = Discretization_;
168 using Model = Model_;
171
173 static constexpr uint dim = Discretization::dim;
174 static constexpr uint n_components = Components::count_fe_functions(0);
175
178 mapping(discretization.get_mapping()), triangulation(discretization.get_triangulation()), json(json),
179 threads(json.get_uint("/discretization/threads")),
180 batch_size(json.get_uint("/discretization/batch_size")),
181 EoM_abs_tol(json.get_double("/discretization/EoM_abs_tol")),
182 EoM_max_iter(json.get_uint("/discretization/EoM_max_iter")),
183 quadrature(1 + json.get_uint("/discretization/overintegration"))
184 {
185 if (this->threads == 0) this->threads = dealii::MultithreadInfo::n_threads() / 2;
186 spdlog::get("log")->info("FV: Using {} threads for assembly.", threads);
187
188 reinit();
189 }
190
191 virtual void reinit_vector(VectorType &vec) const override
192 {
193 const auto block_structure = discretization.get_block_structure();
194 vec.reinit(block_structure[0]);
195 }
196
197 virtual IndexSet get_differential_indices() const override { return IndexSet(); }
198
199 virtual void attach_data_output(DataOutput<dim, VectorType> &data_out, const VectorType &solution,
200 const VectorType &variables, const VectorType &dt_solution = VectorType(),
201 const VectorType &residual = VectorType())
202 {
203 const auto fe_function_names = Components::FEFunction_Descriptor::get_names_vector();
204 std::vector<std::string> fe_function_names_residual;
205 for (const auto &name : fe_function_names)
206 fe_function_names_residual.push_back(name + "_residual");
207 std::vector<std::string> fe_function_names_dot;
208 for (const auto &name : fe_function_names)
209 fe_function_names_dot.push_back(name + "_dot");
210
211 auto &fe_out = data_out.fe_output();
212 fe_out.attach(dof_handler, solution, fe_function_names);
213 if (dt_solution.size() > 0) fe_out.attach(dof_handler, dt_solution, fe_function_names_dot);
214 if (residual.size() > 0) fe_out.attach(dof_handler, residual, fe_function_names_residual);
215
216 // readouts(data_out, solution, variables);
217 }
218
219 virtual void reinit() override
220 {
221 Timer timer;
222
223 // Mass sparsity pattern is diagonal
224 auto dynamic_sparsity_pattern_mass = DynamicSparsityPattern(dof_handler.n_dofs());
225 for (uint i = 0; i < dof_handler.n_dofs(); ++i)
226 dynamic_sparsity_pattern_mass.add(i, i);
227 sparsity_pattern_mass.copy_from(dynamic_sparsity_pattern_mass);
228
229 mass_matrix = SparseMatrix<NumberType>(sparsity_pattern_mass, IdentityMatrix(dof_handler.n_dofs()));
230
231 constexpr uint stencil = 1;
233
234 timings_reinit.push_back(timer.wall_time());
235 }
236
237 virtual void set_time(double t) override { model.set_time(t); }
238
239 virtual const SparsityPattern &get_sparsity_pattern_jacobian() const override
240 {
242 }
243 virtual const SparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
244
245 virtual void residual_variables(VectorType &residual, const VectorType &variables, const VectorType &) override
246 {
247 model.dt_variables(residual, fv_tie(variables));
248 };
249
250 virtual void jacobian_variables(FullMatrix<NumberType> &jacobian, const VectorType &variables,
251 const VectorType &) override
252 {
253 model.template jacobian_variables<0>(jacobian, fv_tie(variables));
254 };
255
256 void readouts(DataOutput<dim, VectorType> &data_out, const VectorType &, const VectorType &variables) const
257 {
258 auto helper = [&](auto EoMfun, auto outputter) {
259 (void)EoMfun;
260 outputter(data_out, Point<0>(), fv_tie(variables));
261 };
262 model.template readouts_multiple(helper, data_out);
263 }
264
265 virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot,
266 NumberType weight) override
267 {
268 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
270 using CopyData = internal::CopyData_R<NumberType>;
271 const auto &constraints = discretization.get_constraints();
272
273 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
274 scratch_data.fe_values.reinit(cell);
275 const auto &fe_v = scratch_data.fe_values;
276 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
277
278 copy_data.reinit(cell, n_dofs);
279 const auto &JxW = fe_v.get_JxW_values();
280 const auto &q_points = fe_v.get_quadrature_points();
281 const auto &q_indices = fe_v.quadrature_point_indices();
282
283 auto &solution = scratch_data.solution;
284 auto &solution_dot = scratch_data.solution_dot;
285 fe_v.get_function_values(solution_global, solution);
286 fe_v.get_function_values(solution_global_dot, solution_dot);
287
288 std::array<NumberType, n_components> mass{};
289 for (const auto &q_index : q_indices) {
290 const auto &x_q = q_points[q_index];
291 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
292
293 for (uint i = 0; i < n_dofs; ++i) {
294 const auto component_i = fe_v.get_fe().system_to_component_index(i).first;
295 copy_data.cell_residual(i) += weight * JxW[q_index] *
296 fe_v.shape_value_component(i, q_index, component_i) *
297 mass[component_i]; // +phi_i(x_q) * mass(x_q, u_q)
298 }
299 }
300 };
301 const auto copier = [&](const CopyData &c) {
302 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, mass);
303 };
304
305 Scratch scratch_data(mapping, discretization.get_fe(), quadrature);
306 CopyData copy_data;
307 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
308
309 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
310 copy_data, flags, nullptr, nullptr, threads, batch_size);
311 }
312
313 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
314 const VectorType &solution_global_dot, NumberType weight_mass,
315 const VectorType &variables = VectorType()) override
316 {
317 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
319 using CopyData = internal::CopyData_R<NumberType>;
320 const auto &constraints = discretization.get_constraints();
321
322 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
323 scratch_data.fe_values.reinit(cell);
324 const auto &fe_v = scratch_data.fe_values;
325 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
326
327 copy_data.reinit(cell, n_dofs);
328 const auto &JxW = fe_v.get_JxW_values();
329 const auto &q_points = fe_v.get_quadrature_points();
330 const auto &q_indices = fe_v.quadrature_point_indices();
331
332 auto &solution = scratch_data.solution;
333 auto &solution_dot = scratch_data.solution_dot;
334 fe_v.get_function_values(solution_global, solution);
335 fe_v.get_function_values(solution_global_dot, solution_dot);
336
337 std::array<NumberType, n_components> mass{};
338 std::array<NumberType, n_components> source{};
339 for (const auto &q_index : q_indices) {
340 const auto &x_q = q_points[q_index];
341 model.mass(mass, x_q, solution[q_index], solution_dot[q_index]);
342 model.source(source, x_q, fv_tie(solution[q_index]));
343
344 for (uint i = 0; i < n_dofs; ++i) {
345 const auto component_i = fe_v.get_fe().system_to_component_index(i).first;
346 copy_data.cell_mass(i) += weight_mass * JxW[q_index] *
347 fe_v.shape_value_component(i, q_index, component_i) *
348 mass[component_i]; // +phi_i(x_q) * mass(x_q, u_q)
349 copy_data.cell_residual(i) += JxW[q_index] * weight * // dx *
350 (fe_v.shape_value_component(i, q_index, component_i) *
351 source[component_i]); // -phi_i(x_q) * source(x_q, u_q)
352 }
353 }
354 };
355 const auto copier = [&](const CopyData &c) {
356 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices, residual);
357 constraints.distribute_local_to_global(c.cell_mass, c.local_dof_indices, residual);
358 };
359
360 Scratch scratch_data(mapping, discretization.get_fe(), quadrature);
361 CopyData copy_data;
362 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
363
364 Timer timer;
365 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
366 copy_data, flags, nullptr, nullptr, threads, batch_size);
367 timings_residual.push_back(timer.wall_time());
368 }
369
370 virtual void jacobian_mass(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
371 const VectorType &solution_global_dot, NumberType alpha = 1.,
372 NumberType beta = 1.) override
373 {
374 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
376 using CopyData = internal::CopyData_J<NumberType>;
377 const auto &constraints = discretization.get_constraints();
378
379 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
380 scratch_data.fe_values.reinit(cell);
381 const auto &fe_v = scratch_data.fe_values;
382 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
383
384 copy_data.reinit(cell, n_dofs, Components::count_extractors());
385 const auto &JxW = fe_v.get_JxW_values();
386 const auto &q_points = fe_v.get_quadrature_points();
387 const auto &q_indices = fe_v.quadrature_point_indices();
388
389 auto &solution = scratch_data.solution;
390 auto &solution_dot = scratch_data.solution_dot;
391 fe_v.get_function_values(solution_global, solution);
392 fe_v.get_function_values(solution_global_dot, solution_dot);
393
396 for (const auto &q_index : q_indices) {
397 const auto &x_q = q_points[q_index];
398 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
399 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
400
401 for (uint i = 0; i < n_dofs; ++i) {
402 const auto component_i = fe_v.get_fe().system_to_component_index(i).first;
403 for (uint j = 0; j < n_dofs; ++j) {
404 const auto component_j = fe_v.get_fe().system_to_component_index(j).first;
405 copy_data.cell_jacobian(i, j) +=
406 JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) *
407 fe_v.shape_value_component(i, q_index, component_i) *
408 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
409 }
410 }
411 }
412 };
413 const auto copier = [&](const CopyData &c) {
414 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
415 };
416
417 Scratch scratch_data(mapping, discretization.get_fe(), quadrature);
418 CopyData copy_data;
419 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
420
421 Timer timer;
422 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
423 copy_data, flags, nullptr, nullptr, threads, batch_size);
424 timings_jacobian.push_back(timer.wall_time());
425 }
426
427 virtual void jacobian(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global, NumberType weight,
428 const VectorType &solution_global_dot, NumberType alpha, NumberType beta,
429 const VectorType &variables = VectorType()) override
430 {
431 using Iterator = typename DoFHandler<dim>::active_cell_iterator;
433 using CopyData = internal::CopyData_J<NumberType>;
434 const auto &constraints = discretization.get_constraints();
435
436 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData &copy_data) {
437 scratch_data.fe_values.reinit(cell);
438 const auto &fe_v = scratch_data.fe_values;
439 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
440
441 copy_data.reinit(cell, n_dofs, Components::count_extractors());
442 const auto &JxW = fe_v.get_JxW_values();
443 const auto &q_points = fe_v.get_quadrature_points();
444 const auto &q_indices = fe_v.quadrature_point_indices();
445
446 auto &solution = scratch_data.solution;
447 auto &solution_dot = scratch_data.solution_dot;
448 fe_v.get_function_values(solution_global, solution);
449 fe_v.get_function_values(solution_global_dot, solution_dot);
450
454 for (const auto &q_index : q_indices) {
455 const auto &x_q = q_points[q_index];
456 model.template jacobian_mass<0>(j_mass, x_q, solution[q_index], solution_dot[q_index]);
457 model.template jacobian_mass<1>(j_mass_dot, x_q, solution[q_index], solution_dot[q_index]);
458 model.template jacobian_source<0, 0>(j_source, x_q, fv_tie(solution[q_index]));
459
460 for (uint i = 0; i < n_dofs; ++i) {
461 const auto component_i = fe_v.get_fe().system_to_component_index(i).first;
462 for (uint j = 0; j < n_dofs; ++j) {
463 const auto component_j = fe_v.get_fe().system_to_component_index(j).first;
464 copy_data.cell_jacobian(i, j) +=
465 weight * JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) * // dx * phi_j * (
466 (fe_v.shape_value_component(i, q_index, component_i) *
467 j_source(component_i, component_j)); // -phi_i * jsource)
468 copy_data.cell_mass_jacobian(i, j) +=
469 JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) *
470 fe_v.shape_value_component(i, q_index, component_i) *
471 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
472 }
473 }
474 }
475 };
476 const auto copier = [&](const CopyData &c) {
477 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices, jacobian);
478 constraints.distribute_local_to_global(c.cell_mass_jacobian, c.local_dof_indices, jacobian);
479 };
480
481 Scratch scratch_data(mapping, discretization.get_fe(), quadrature);
482 CopyData copy_data;
483 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells;
484
485 Timer timer;
486 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
487 copy_data, flags, nullptr, nullptr, threads, batch_size);
488 timings_jacobian.push_back(timer.wall_time());
489 }
490
491 void build_sparsity(SparsityPattern &sparsity_pattern, const DoFHandler<dim> &to_dofh,
492 const DoFHandler<dim> &from_dofh, const int stencil = 1,
493 bool add_extractor_dofs = false) const
494 {
495 const auto &triangulation = discretization.get_triangulation();
496
497 DynamicSparsityPattern dsp(dof_handler.n_dofs(), dof_handler.n_dofs());
498
499 const auto to_dofs_per_cell = to_dofh.get_fe().dofs_per_cell;
500 const auto from_dofs_per_cell = from_dofh.get_fe().dofs_per_cell;
501
502 for (const auto &t_cell : triangulation.active_cell_iterators()) {
503 std::vector<types::global_dof_index> to_dofs(to_dofs_per_cell);
504 std::vector<types::global_dof_index> from_dofs(from_dofs_per_cell);
505 const auto to_cell = t_cell->as_dof_handler_iterator(to_dofh);
506 const auto from_cell = t_cell->as_dof_handler_iterator(from_dofh);
507 to_cell->get_dof_indices(to_dofs);
508 from_cell->get_dof_indices(from_dofs);
509
510 std::function<void(decltype(from_cell) &, const int)> add_all_neighbor_dofs = [&](const auto &from_cell,
511 const int stencil_level) {
512 for (const auto face_no : from_cell->face_indices()) {
513 const auto face = from_cell->face(face_no);
514 if (!face->at_boundary()) {
515 auto neighbor_cell = from_cell->neighbor(face_no);
516
517 if (dim == 1)
518 while (neighbor_cell->has_children())
519 neighbor_cell = neighbor_cell->child(face_no == 0 ? 1 : 0);
520
521 // add all children
522 else if (neighbor_cell->has_children()) {
523 throw std::runtime_error("not yet implemented lol");
524 }
525
526 if (!neighbor_cell->is_active()) continue;
527
528 std::vector<types::global_dof_index> tmp(from_dofs_per_cell);
529 neighbor_cell->get_dof_indices(tmp);
530
531 from_dofs.insert(std::end(from_dofs), std::begin(tmp), std::end(tmp));
532
533 if (stencil_level < stencil) add_all_neighbor_dofs(neighbor_cell, stencil_level + 1);
534 }
535 }
536 };
537
538 add_all_neighbor_dofs(from_cell, 1);
539
540 for (const auto i : to_dofs)
541 for (const auto j : from_dofs)
542 dsp.add(i, j);
543 }
544
545 // if (add_extractor_dofs)
546 // for (uint row = 0; row < dsp.n_rows(); ++row)
547 // dsp.add_row_entries(row, extractor_dof_indices);
548 sparsity_pattern.copy_from(dsp);
549 }
550
551 void log(const std::string logger)
552 {
553 std::stringstream ss;
554 ss << "FV Assembler: " << std::endl;
555 ss << " Reinit: " << average_time_reinit() * 1000 << "ms (" << num_reinits() << ")" << std::endl;
556 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
557 << std::endl;
558 ss << " Jacobian: " << average_time_jacobian_assembly() * 1000 << "ms (" << num_jacobians() << ")"
559 << std::endl;
560 spdlog::get(logger)->info(ss.str());
561 }
562
563 double average_time_reinit() const
564 {
565 double t = 0.;
566 double n = timings_reinit.size();
567 for (const auto &t_ : timings_reinit)
568 t += t_ / n;
569 return t;
570 }
571 uint num_reinits() const { return timings_reinit.size(); }
572
574 {
575 double t = 0.;
576 double n = timings_residual.size();
577 for (const auto &t_ : timings_residual)
578 t += t_ / n;
579 return t;
580 }
581 uint num_residuals() const { return timings_residual.size(); }
582
584 {
585 double t = 0.;
586 double n = timings_jacobian.size();
587 for (const auto &t_ : timings_jacobian)
588 t += t_ / n;
589 return t;
590 }
591 uint num_jacobians() const { return timings_jacobian.size(); }
592
593 protected:
596 const DoFHandler<dim> &dof_handler;
597 const Mapping<dim> &mapping;
598 const Triangulation<dim> &triangulation;
599
601
604
605 mutable typename Triangulation<dim>::cell_iterator EoM_cell;
606 typename Triangulation<dim>::cell_iterator old_EoM_cell;
607 const double EoM_abs_tol;
609
610 const QGauss<dim> quadrature;
611
612 SparsityPattern sparsity_pattern_mass;
614 SparseMatrix<NumberType> mass_matrix;
615
616 std::vector<double> timings_reinit;
617 std::vector<double> timings_residual;
618 std::vector<double> timings_jacobian;
619 };
620 } // namespace KurganovTadmor
621 } // namespace FV
622} // namespace DiFfRG
This is the general assembler interface for any kind of discretization. An assembler is responsible f...
Definition abstract_assembler.hh:39
Class to manage writing to files. FEM functions are written to vtk files and other data is written to...
Definition data_output.hh:20
FEOutput< dim, VectorType > & fe_output()
Returns a reference to the FEOutput object used to write FEM functions to .vtu files and ....
static constexpr uint dim
Definition discretization.hh:38
Vector< NumberType > VectorType
Definition discretization.hh:35
NumberType_ NumberType
Definition discretization.hh:34
Components_ Components
Definition discretization.hh:33
Definition KurganovTadmor.hh:147
Assembler(Discretization &discretization, Model &model, const JSONValue &json)
Definition KurganovTadmor.hh:176
double average_time_jacobian_assembly() const
Definition KurganovTadmor.hh:583
uint num_residuals() const
Definition KurganovTadmor.hh:581
Triangulation< dim >::cell_iterator old_EoM_cell
Definition KurganovTadmor.hh:606
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition KurganovTadmor.hh:219
Model & model
Definition KurganovTadmor.hh:595
const DoFHandler< dim > & dof_handler
Definition KurganovTadmor.hh:596
virtual void set_time(double t) override
Set the current time. The assembler should usually just forward this to the numerical model.
Definition KurganovTadmor.hh:237
const QGauss< dim > quadrature
Definition KurganovTadmor.hh:610
const Mapping< dim > & mapping
Definition KurganovTadmor.hh:597
virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight) override
Definition KurganovTadmor.hh:265
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 KurganovTadmor.hh:427
const uint EoM_max_iter
Definition KurganovTadmor.hh:608
Model_ Model
Definition KurganovTadmor.hh:168
typename Discretization::VectorType VectorType
Definition KurganovTadmor.hh:170
double average_time_residual_assembly() const
Definition KurganovTadmor.hh:573
SparsityPattern sparsity_pattern_jacobian
Definition KurganovTadmor.hh:613
SparsityPattern sparsity_pattern_mass
Definition KurganovTadmor.hh:612
static constexpr auto e_tie(T &&...t)
Definition KurganovTadmor.hh:159
virtual void jacobian_mass(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha=1., NumberType beta=1.) override
Definition KurganovTadmor.hh:370
static constexpr auto v_tie(T &&...t)
Definition KurganovTadmor.hh:154
const uint batch_size
Definition KurganovTadmor.hh:603
double average_time_reinit() const
Definition KurganovTadmor.hh:563
static constexpr uint n_components
Definition KurganovTadmor.hh:174
uint num_reinits() const
Definition KurganovTadmor.hh:571
uint num_jacobians() const
Definition KurganovTadmor.hh:591
typename Discretization::NumberType NumberType
Definition KurganovTadmor.hh:169
virtual void attach_data_output(DataOutput< dim, VectorType > &data_out, const VectorType &solution, const VectorType &variables, const VectorType &dt_solution=VectorType(), const VectorType &residual=VectorType())
Definition KurganovTadmor.hh:199
Triangulation< dim >::cell_iterator EoM_cell
Definition KurganovTadmor.hh:605
SparseMatrix< NumberType > mass_matrix
Definition KurganovTadmor.hh:614
virtual const SparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition KurganovTadmor.hh:243
std::vector< double > timings_jacobian
Definition KurganovTadmor.hh:618
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 KurganovTadmor.hh:313
Discretization & discretization
Definition KurganovTadmor.hh:594
uint threads
Definition KurganovTadmor.hh:602
void log(const std::string logger)
Definition KurganovTadmor.hh:551
typename Discretization::Components Components
Definition KurganovTadmor.hh:172
static constexpr uint dim
Definition KurganovTadmor.hh:173
virtual IndexSet get_differential_indices() const override
Obtain the dofs which contain time derivatives.
Definition KurganovTadmor.hh:197
virtual void reinit_vector(VectorType &vec) const override
Definition KurganovTadmor.hh:191
const double EoM_abs_tol
Definition KurganovTadmor.hh:607
virtual const SparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition KurganovTadmor.hh:239
std::vector< double > timings_residual
Definition KurganovTadmor.hh:617
virtual void residual_variables(VectorType &residual, const VectorType &variables, const VectorType &) override
Definition KurganovTadmor.hh:245
void build_sparsity(SparsityPattern &sparsity_pattern, const DoFHandler< dim > &to_dofh, const DoFHandler< dim > &from_dofh, const int stencil=1, bool add_extractor_dofs=false) const
Definition KurganovTadmor.hh:491
const JSONValue & json
Definition KurganovTadmor.hh:600
Discretization_ Discretization
Definition KurganovTadmor.hh:167
std::vector< double > timings_reinit
Definition KurganovTadmor.hh:616
virtual void jacobian_variables(FullMatrix< NumberType > &jacobian, const VectorType &variables, const VectorType &) override
Definition KurganovTadmor.hh:250
auto fv_tie(T &&...t)
Definition KurganovTadmor.hh:149
const Triangulation< dim > & triangulation
Definition KurganovTadmor.hh:598
void readouts(DataOutput< dim, VectorType > &data_out, const VectorType &, const VectorType &variables) const
Definition KurganovTadmor.hh:256
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
Definition types.hh:6
Definition complex_math.hh:10
unsigned int uint
Definition utils.hh:24
std::array< double, 2 > values
Definition KurganovTadmor.hh:135
std::array< uint, 2 > cell_indices
Definition KurganovTadmor.hh:134
Definition KurganovTadmor.hh:132
std::vector< CopyFaceData_I > face_data
Definition KurganovTadmor.hh:137
uint cell_index
Definition KurganovTadmor.hh:139
double value
Definition KurganovTadmor.hh:138
FullMatrix< NumberType > cell_jacobian
Definition KurganovTadmor.hh:103
void reinit(const FEInterfaceValues< dim > &fe_iv, uint n_extractors)
Definition KurganovTadmor.hh:107
std::vector< types::global_dof_index > joint_dof_indices
Definition KurganovTadmor.hh:105
FullMatrix< NumberType > extractor_cell_jacobian
Definition KurganovTadmor.hh:104
Definition KurganovTadmor.hh:101
FullMatrix< NumberType > extractor_cell_jacobian
Definition KurganovTadmor.hh:117
std::vector< CopyDataFace_J > face_data
Definition KurganovTadmor.hh:120
FullMatrix< NumberType > cell_jacobian
Definition KurganovTadmor.hh:116
FullMatrix< NumberType > cell_mass_jacobian
Definition KurganovTadmor.hh:118
std::vector< types::global_dof_index > local_dof_indices
Definition KurganovTadmor.hh:119
void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
Definition KurganovTadmor.hh:122
Vector< NumberType > cell_residual
Definition KurganovTadmor.hh:76
std::vector< types::global_dof_index > joint_dof_indices
Definition KurganovTadmor.hh:77
void reinit(const FEInterfaceValues< dim > &fe_iv)
Definition KurganovTadmor.hh:79
void reinit(const Iterator &cell, uint dofs_per_cell)
Definition KurganovTadmor.hh:91
Vector< NumberType > cell_mass
Definition KurganovTadmor.hh:87
std::vector< CopyDataFace_R > face_data
Definition KurganovTadmor.hh:89
std::vector< types::global_dof_index > local_dof_indices
Definition KurganovTadmor.hh:88
Vector< NumberType > cell_residual
Definition KurganovTadmor.hh:86
Class to hold data for each assembly thread, i.e. FEValues for cells, interfaces, as well as pre-allo...
Definition KurganovTadmor.hh:42
typename Discretization::NumberType NumberType
Definition KurganovTadmor.hh:44
Vector< NumberType > VectorType
Definition KurganovTadmor.hh:45
const uint n_components
Definition KurganovTadmor.hh:65
std::vector< VectorType > solution
Definition KurganovTadmor.hh:69
std::vector< VectorType > solution_dot
Definition KurganovTadmor.hh:70
ScratchData(const ScratchData< Discretization > &scratch_data)
Definition KurganovTadmor.hh:56
ScratchData(const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values)
Definition KurganovTadmor.hh:47
FEValues< dim > fe_values
Definition KurganovTadmor.hh:67
static constexpr int dim
Definition KurganovTadmor.hh:43
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