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

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/variables/assembler/variables.hh Source File
DiFfRG
variables.hh
Go to the documentation of this file.
1#pragma once
2
3// standard library
4#include <sstream>
5
6// external libraries
7#include <deal.II/base/multithread_info.h>
8#include <deal.II/base/timer.h>
9#include <deal.II/dofs/dof_handler.h>
10#include <deal.II/lac/full_matrix.h>
11#include <deal.II/lac/sparse_matrix.h>
12#include <deal.II/lac/vector.h>
13
14// DiFfRG
18
19namespace DiFfRG
20{
21 namespace Variables
22 {
23 using namespace dealii;
24
25 template <typename... T> auto fe_tie(T &&...t)
26 {
27 return named_tuple<std::tuple<T &...>, StringSet<"variables">>(std::tie(t...));
28 }
29
35 template <typename Model_> class Assembler : public AbstractAssembler<Vector<double>, SparseMatrix<double>, 0>
36 {
37 constexpr static int nothing = 0;
38
39 public:
40 using Model = Model_;
41 using NumberType = double;
42 using VectorType = Vector<double>;
43
44 using Components = typename Model_::Components;
45 static constexpr uint dim = 0;
46
47 Assembler(Model &model, const JSONValue &json) : model(model), threads(json.get_uint("/discretization/threads"))
48 {
49 if (threads == 0) threads = dealii::MultithreadInfo::n_threads() / 2;
50 static_assert(Components::count_fe_functions() == 0, "The pure variable assembler cannot handle FE functions!");
51 reinit();
52 }
53
54 virtual void reinit_vector(VectorType &vec) const override { vec.reinit(0); }
55
56 virtual IndexSet get_differential_indices() const override { return IndexSet(); }
57
58 virtual void attach_data_output(DataOutput<dim, VectorType> &data_out, const VectorType &solution,
59 const VectorType &variables, const VectorType &dt_solution = VectorType(),
61 {
62 (void)dt_solution;
63 (void)residual;
64 readouts(data_out, solution, variables);
65 }
66
67 virtual void reinit() override {}
68
69 virtual void set_time(double t) override { model.set_time(t); }
70
71 virtual const SparsityPattern &get_sparsity_pattern_jacobian() const override
72 {
74 }
75 virtual const SparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
76
77 virtual void residual_variables(VectorType &residual, const VectorType &variables, const VectorType &) override
78 {
79 Timer timer;
80 model.dt_variables(residual, fe_tie(variables));
81 Kokkos::fence();
82 timings_residual.push_back(timer.wall_time());
83 };
84
85 virtual void jacobian_variables(FullMatrix<NumberType> &jacobian, const VectorType &variables,
86 const VectorType &) override
87 {
88 Timer timer;
89 model.template jacobian_variables<0>(jacobian, fe_tie(variables));
90 Kokkos::fence();
91 timings_jacobian.push_back(timer.wall_time());
92 };
93
94 void readouts(DataOutput<dim, VectorType> &data_out, const VectorType &, const VectorType &variables) const
95 {
96 auto helper = [&](auto EoMfun, auto outputter) {
97 (void)EoMfun;
98 outputter(data_out, Point<0>(), fe_tie(variables));
99 };
100 model.readouts_multiple(helper, data_out);
101 }
102
103 virtual void mass(VectorType &, const VectorType &, const VectorType &, NumberType) override {}
104
105 virtual void residual(VectorType &, const VectorType &, NumberType, const VectorType &, NumberType,
106 const VectorType &variables = VectorType()) override
107 {
108 (void)variables;
109 }
110
111 virtual void jacobian_mass(SparseMatrix<NumberType> &, const VectorType &, const VectorType &, NumberType,
112 NumberType) override
113 {
114 }
115
116 virtual void jacobian(SparseMatrix<NumberType> &, const VectorType &, NumberType, const VectorType &, NumberType,
117 NumberType, const VectorType &variables = VectorType()) override
118 {
119 (void)variables;
120 }
121
122 void log(const std::string logger) const
123 {
124 std::stringstream ss;
125 ss << "Variable Assembler: " << std::endl;
126 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
127 << std::endl;
128 ss << " Jacobian: " << average_time_jacobian_assembly() * 1000 << "ms (" << num_jacobians() << ")"
129 << std::endl;
130 spdlog::get(logger)->info(ss.str());
131 }
132
134 {
135 double t = 0.;
136 double n = timings_residual.size();
137 for (const auto &t_ : timings_residual)
138 t += t_ / n;
139 return t;
140 }
141 uint num_residuals() const { return timings_residual.size(); }
142
144 {
145 double t = 0.;
146 double n = timings_jacobian.size();
147 for (const auto &t_ : timings_jacobian)
148 t += t_ / n;
149 return t;
150 }
151 uint num_jacobians() const { return timings_jacobian.size(); }
152
153 private:
155
157
158 SparsityPattern sparsity_pattern_mass;
160 SparseMatrix<NumberType> mass_matrix;
161
162 std::vector<double> timings_residual;
163 std::vector<double> timings_jacobian;
164 };
165 } // namespace Variables
166} // 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
A wrapper around the boost json value class.
Definition json.hh:19
The basic assembler that can be used for any standard CG scheme with flux and source.
Definition variables.hh:36
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition variables.hh:67
virtual void reinit_vector(VectorType &vec) const override
Reinitialize an arbitrary vector so that it has the correct size and structure.
Definition variables.hh:54
virtual void residual_variables(VectorType &residual, const VectorType &variables, const VectorType &) override
When coupling the spatial discretization to additional variables, this function should calculate the ...
Definition variables.hh:77
std::vector< double > timings_jacobian
Definition variables.hh:163
virtual void set_time(double t) override
Set the current time. The assembler should usually just forward this to the numerical model.
Definition variables.hh:69
virtual void residual(VectorType &, const VectorType &, NumberType, const VectorType &, NumberType, const VectorType &variables=VectorType()) override
Definition variables.hh:105
virtual IndexSet get_differential_indices() const override
Obtain the dofs which contain time derivatives.
Definition variables.hh:56
void log(const std::string logger) const
Definition variables.hh:122
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())
Attach any data output to the DataOutput object provided. This can be used to extract additional data...
Definition variables.hh:58
uint threads
Definition variables.hh:156
Model & model
Definition variables.hh:154
SparsityPattern sparsity_pattern_mass
Definition variables.hh:158
static constexpr uint dim
Definition variables.hh:45
uint num_jacobians() const
Definition variables.hh:151
uint num_residuals() const
Definition variables.hh:141
void readouts(DataOutput< dim, VectorType > &data_out, const VectorType &, const VectorType &variables) const
Definition variables.hh:94
virtual void jacobian_mass(SparseMatrix< NumberType > &, const VectorType &, const VectorType &, NumberType, NumberType) override
Definition variables.hh:111
virtual const SparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition variables.hh:75
typename Model_::Components Components
Definition variables.hh:44
SparsityPattern sparsity_pattern_jacobian
Definition variables.hh:159
Model_ Model
Definition variables.hh:40
double average_time_jacobian_assembly()
Definition variables.hh:143
std::vector< double > timings_residual
Definition variables.hh:162
Assembler(Model &model, const JSONValue &json)
Definition variables.hh:47
SparseMatrix< NumberType > mass_matrix
Definition variables.hh:160
virtual void jacobian(SparseMatrix< NumberType > &, const VectorType &, NumberType, const VectorType &, NumberType, NumberType, const VectorType &variables=VectorType()) override
Definition variables.hh:116
virtual void mass(VectorType &, const VectorType &, const VectorType &, NumberType) override
Definition variables.hh:103
virtual const SparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition variables.hh:71
Vector< double > VectorType
Definition variables.hh:42
static constexpr int nothing
Definition variables.hh:37
double NumberType
Definition variables.hh:41
double average_time_residual_assembly()
Definition variables.hh:133
virtual void jacobian_variables(FullMatrix< NumberType > &jacobian, const VectorType &variables, const VectorType &) override
Definition variables.hh:85
auto fe_tie(T &&...t)
Definition variables.hh:25
Definition complex_math.hh:10
unsigned int uint
Definition utils.hh:24
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