DiFfRG
Loading...
Searching...
No Matches
abstract_assembler.hh
Go to the documentation of this file.
1#pragma once
2
3// DiFfRG
6
7namespace DiFfRG
8{
9 using namespace dealii;
10
38 template <typename VectorType, typename SparseMatrixType, uint dim> class AbstractAssembler
39 {
40 public:
42
52 virtual void attach_data_output(DataOutput<dim, VectorType> &data_out, const VectorType &solution,
53 const VectorType &variables = VectorType(),
54 const VectorType &dt_solution = VectorType(),
55 const VectorType &residual = VectorType()) = 0;
56
61 virtual void reinit() = 0;
62
68 virtual void set_time(double t) = 0;
69
76
82 virtual IndexSet get_differential_indices() const = 0;
83
89 virtual void reinit_vector(VectorType &vector) const = 0;
90
96 virtual const SparseMatrixType &get_mass_matrix() const = 0;
97
102
111 virtual void residual_variables([[maybe_unused]] VectorType &residual, [[maybe_unused]] const VectorType &variables,
112 [[maybe_unused]] const VectorType &spatial_solution)
113 {
114 throw std::runtime_error("Not implemented!");
115 };
116
125 virtual void jacobian_variables([[maybe_unused]] FullMatrix<NumberType> &jacobian,
126 [[maybe_unused]] const VectorType &variables,
127 [[maybe_unused]] const VectorType &spatial_solution)
128 {
129 throw std::runtime_error("Not implemented!");
130 };
131
141 void mass(VectorType &mass, const VectorType &solution_global, NumberType weight)
142 {
143 this->mass(mass, solution_global, solution_global, weight);
144 }
145
156 virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot,
157 NumberType weight) = 0;
158
177 void residual(VectorType &residual, const VectorType &solution_global, NumberType weight, NumberType weight_mass,
178 const VectorType &variables = VectorType())
179 {
180 this->residual(residual, solution_global, weight, solution_global, weight_mass, variables);
181 }
182
202 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
203 const VectorType &solution_global_dot, NumberType weight_mass,
204 const VectorType &variables = VectorType()) = 0;
205
218 void jacobian_mass(SparseMatrixType &jacobian, const VectorType &solution_global, NumberType mass_weight = 1.)
219 {
220 this->jacobian_mass(jacobian, solution_global, solution_global, mass_weight, 0.);
221 }
222
238 virtual void jacobian_mass(SparseMatrixType &jacobian, const VectorType &solution_global,
239 const VectorType &solution_global_dot, NumberType alpha = 1., NumberType beta = 1.) = 0;
240
251 void jacobian(SparseMatrixType &jacobian, const VectorType &solution_global, NumberType weight,
252 NumberType mass_weight, const VectorType &variables = VectorType())
253 {
254 this->jacobian(jacobian, solution_global, weight, solution_global, mass_weight, mass_weight, variables);
255 }
256
267 virtual void jacobian(SparseMatrixType &jacobian, const VectorType &solution_global, NumberType weight,
268 const VectorType &solution_global_dot, NumberType alpha, NumberType beta,
269 const VectorType &variables = VectorType()) = 0;
271 };
272} // namespace DiFfRG
This is the general assembler interface for any kind of discretization. An assembler is responsible f...
Definition abstract_assembler.hh:39
void mass(VectorType &mass, const VectorType &solution_global, NumberType weight)
Calculates the mass for an ODE.
Definition abstract_assembler.hh:141
virtual void jacobian_variables(FullMatrix< NumberType > &jacobian, const VectorType &variables, const VectorType &spatial_solution)
When coupling the spatial discretization to additional variables, this function should calculate the ...
Definition abstract_assembler.hh:125
void jacobian_mass(SparseMatrixType &jacobian, const VectorType &solution_global, NumberType mass_weight=1.)
Calculates the jacobian of the mass function for an ODE.
Definition abstract_assembler.hh:218
virtual void reinit_vector(VectorType &vector) const =0
Reinitialize an arbitrary vector so that it has the correct size and structure.
void jacobian(SparseMatrixType &jacobian, const VectorType &solution_global, NumberType weight, NumberType mass_weight, const VectorType &variables=VectorType())
Calculates the jacobian of the residual function for an ODE.
Definition abstract_assembler.hh:251
virtual void jacobian_mass(SparseMatrixType &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha=1., NumberType beta=1.)=0
Calculates the jacobian of the mass function for a DAE.
virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight)=0
Calculates the mass for a DAE.
typename get_type::NumberType< VectorType > NumberType
Definition abstract_assembler.hh:41
void residual(VectorType &residual, const VectorType &solution_global, NumberType weight, NumberType weight_mass, const VectorType &variables=VectorType())
Calculates the residual for an ODE.
Definition abstract_assembler.hh:177
virtual void residual_variables(VectorType &residual, const VectorType &variables, const VectorType &spatial_solution)
When coupling the spatial discretization to additional variables, this function should calculate the ...
Definition abstract_assembler.hh:111
virtual void attach_data_output(DataOutput< dim, VectorType > &data_out, const VectorType &solution, const VectorType &variables=VectorType(), const VectorType &dt_solution=VectorType(), const VectorType &residual=VectorType())=0
Attach any data output to the DataOutput object provided. This can be used to extract additional data...
virtual IndexSet get_differential_indices() const =0
Obtain the dofs which contain time derivatives.
virtual void jacobian(SparseMatrixType &jacobian, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType alpha, NumberType beta, const VectorType &variables=VectorType())=0
Calculates the jacobian of the residual function for a DAE.
virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType weight_mass, const VectorType &variables=VectorType())=0
Calculates the residual for a DAE.
virtual void set_time(double t)=0
Set the current time. The assembler should usually just forward this to the numerical model.
virtual const SparseMatrixType & get_mass_matrix() const =0
Obtain the mass matrix.
virtual const get_type::SparsityPattern< SparseMatrixType > & get_sparsity_pattern_jacobian() const =0
Obtain the sparsity pattern of the jacobian matrix.
virtual void reinit()=0
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Class to manage writing to files. FEM functions are written to vtk files and other data is written to...
Definition data_output.hh:20
typename internal::_NumberType< VectorType >::value NumberType
Definition types.hh:61
typename internal::_SparsityPattern< SparseMatrixType >::value SparsityPattern
Definition types.hh:64
Definition complex_math.hh:14