DiFfRG
Loading...
Searching...
No Matches
DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim > Class Template Referenceabstract

This is the general assembler interface for any kind of discretization. An assembler is responsible for calculating residuals and their jacobians for any given discretization, including both the spatial part and any further variables. Any assembler for a specific spatial discretization must fully implement this interface. More...

#include <abstract_assembler.hh>

Public Types

using NumberType = typename get_type::NumberType<VectorType>
 

Public Member Functions

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 from the solution and write it to the output file. This includes both derivatives and other spatial functions, as well as single values that can be appended to the .csv file.
 
virtual void reinit ()=0
 Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
 
virtual void set_time (double t)=0
 Set the current time. The assembler should usually just forward this to the numerical model.
 
virtual const get_type::SparsityPattern< SparseMatrixType > & get_sparsity_pattern_jacobian () const =0
 Obtain the sparsity pattern of the jacobian matrix.
 
virtual IndexSet get_differential_indices () const =0
 Obtain the dofs which contain time derivatives.
 
virtual void reinit_vector (VectorType &vector) const =0
 Reinitialize an arbitrary vector so that it has the correct size and structure.
 
virtual const SparseMatrixType & get_mass_matrix () const =0
 Obtain the mass matrix.
 
Residual and jacobian functions
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 residual for the additional variables.
 
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 jacobian for the additional variables.
 
void mass (VectorType &mass, const VectorType &solution_global, NumberType weight)
 Calculates the mass \(m_i(u)\) for an ODE.
 
virtual void mass (VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight)=0
 Calculates the mass \(m(u, \partial_t u)\) for a DAE.
 
void residual (VectorType &residual, const VectorType &solution_global, NumberType weight, NumberType weight_mass, const VectorType &variables=VectorType())
 Calculates the residual for an ODE.
 
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.
 
void jacobian_mass (SparseMatrixType &jacobian, const VectorType &solution_global, NumberType mass_weight=1.)
 Calculates the jacobian of the mass function for an ODE.
 
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.
 
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.
 
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.
 

Detailed Description

template<typename VectorType, typename SparseMatrixType, uint dim>
class DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >

This is the general assembler interface for any kind of discretization. An assembler is responsible for calculating residuals and their jacobians for any given discretization, including both the spatial part and any further variables. Any assembler for a specific spatial discretization must fully implement this interface.

In general, we have either a system of differential and algebraic equations which can be written as

\[\large m_i(u, \partial_t u) + \partial_x F_i(u, \dots) + s_i(u, \dots) = 0\,,i\in \{0, \dots, N_f-1\}\,, \]

or a system of differential equations which can be written as

\[\large m_{ij}\partial_t u_j + \partial_x F_i(u, \dots) + s_i(u, \dots) = 0\,,i\in \{0, \dots, N_f-1\}\,. \]

Here, \(u_i\) is the solution vector, \(\partial_t u_i\) is the time derivative of the solution vector, \(m_i\) is the mass function and \(m_{ij}\) the mass matrix, \(F_i\) is called the flux, and \(s_i\) is the source term. We will henceforth call the first equation a DAE (differential algebraic equation) and the second equation an ODE (ordinary differential equation).

For more information, please also refer to the guide and the documentation of DiFfRG::def::AbstractModel.

The second equation is used in all residual and jacobian functions which do not take the time derivative as an explicit argument, while the first equation is used in all other methods.

Template Parameters
VectorTypeThe vector type used in the spatial discretization
dimThe dimension of the spatial discretization

Member Typedef Documentation

◆ NumberType

template<typename VectorType , typename SparseMatrixType , uint dim>
using DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::NumberType = typename get_type::NumberType<VectorType>

Member Function Documentation

◆ attach_data_output()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::attach_data_output ( DataOutput< dim, VectorType > & data_out,
const VectorType & solution,
const VectorType & variables = VectorType(),
const VectorType & dt_solution = VectorType(),
const VectorType & residual = VectorType() )
pure virtual

Attach any data output to the DataOutput object provided. This can be used to extract additional data from the solution and write it to the output file. This includes both derivatives and other spatial functions, as well as single values that can be appended to the .csv file.

Parameters
data_outThe DataOutput object to attach the data to
solutionThe spatial solution vector
variablesThe additional variables vector

Implemented in DiFfRG::Variables::Assembler< Model_ >.

◆ get_differential_indices()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual IndexSet DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::get_differential_indices ( ) const
pure virtual

Obtain the dofs which contain time derivatives.

Returns
IndexSet The indices of the dofs which contain time derivatives

Implemented in DiFfRG::FEMAssembler< Discretization_, Model_ >, DiFfRG::FV::KurganovTadmor::Assembler< Discretization_, Model_ >, DiFfRG::LDG::LDGAssemblerBase< Discretization_, Model_ >, and DiFfRG::Variables::Assembler< Model_ >.

◆ get_mass_matrix()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual const SparseMatrixType & DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::get_mass_matrix ( ) const
pure virtual

◆ get_sparsity_pattern_jacobian()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual const get_type::SparsityPattern< SparseMatrixType > & DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::get_sparsity_pattern_jacobian ( ) const
pure virtual

◆ jacobian() [1/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::jacobian ( SparseMatrixType & jacobian,
const VectorType & solution_global,
NumberType weight,
const VectorType & solution_global_dot,
NumberType alpha,
NumberType beta,
const VectorType & variables = VectorType() )
pure virtual

Calculates the jacobian of the residual function for a DAE.

Parameters
jacobian
solution_global
weight
solution_global_dot
alpha
beta

◆ jacobian() [2/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::jacobian ( SparseMatrixType & jacobian,
const VectorType & solution_global,
NumberType weight,
NumberType mass_weight,
const VectorType & variables = VectorType() )
inline

Calculates the jacobian of the residual function for an ODE.

Parameters
jacobian
solution_global
weight
mass_weight

◆ jacobian_mass() [1/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::jacobian_mass ( SparseMatrixType & jacobian,
const VectorType & solution_global,
const VectorType & solution_global_dot,
NumberType alpha = 1.,
NumberType beta = 1. )
pure virtual

Calculates the jacobian of the mass function for a DAE.

This function calculates the jacobian of the mass, i.e.

\[\large \alpha\, \int_\Omega\, \frac{\partial m_i(u_k, \partial_t u_k)}{\partial (\partial_t u_j)} \phi_i + \beta\, \int_\Omega\, \frac{\partial m_i(u_k, \partial_t u_k)}{\partial u_j} \phi_i \]

.

Parameters
jacobianThe jacobian matrix to be filled. This function adds to the matrix.
solution_globalThe spatial solution vector
solution_global_dotThe spatial solution vector time derivative
alphaThe weight \( \alpha \) for the derivative with respect to \( \partial_t \)
betaThe weight \( \beta \) for the derivative with respect to \( u \)

◆ jacobian_mass() [2/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::jacobian_mass ( SparseMatrixType & jacobian,
const VectorType & solution_global,
NumberType mass_weight = 1. )
inline

Calculates the jacobian of the mass function for an ODE.

This function calculates the jacobian of the mass, i.e.

\[\large w_m\, \int_\Omega\, m_{ij} \phi_i \]

.

Parameters
jacobianThe jacobian matrix to be filled. This function adds to the matrix.
solution_globalThe spatial solution vector
mass_weightThe weight for the mass \( w_m \)

◆ jacobian_variables()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::jacobian_variables ( FullMatrix< NumberType > & jacobian,
const VectorType & variables,
const VectorType & spatial_solution )
inlinevirtual

When coupling the spatial discretization to additional variables, this function should calculate the jacobian for the additional variables.

Parameters
jacobianThe jacobian matrix of the additional variables
variablesThe current additional variables vector
spatial_solutionThe spatial solution vector, which is needed for the calculation of the jacobian

◆ mass() [1/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::mass ( VectorType & mass,
const VectorType & solution_global,
const VectorType & solution_global_dot,
NumberType weight )
pure virtual

Calculates the mass \(m(u, \partial_t u)\) for a DAE.

This function calculates the mass as \( w_m \, \int_\Omega\, m_i(u_j, \partial_t u_j) \phi_i\).

Parameters
massThe mass vector to be filled
solution_globalThe spatial solution vector
solution_global_dotThe spatial solution vector time derivative
weightThe weight for the mass \( w_m \)

◆ mass() [2/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::mass ( VectorType & mass,
const VectorType & solution_global,
NumberType weight )
inline

Calculates the mass \(m_i(u)\) for an ODE.

This function calculates the mass as \( w_m \, \int_\Omega\, m_ij \partial_t u_j \phi_i\).

Parameters
massThe mass vector to be filled
solution_globalThe spatial solution vector
weightThe weight for the mass \( w_m \)

◆ reinit()

◆ reinit_vector()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::reinit_vector ( VectorType & vector) const
pure virtual

Reinitialize an arbitrary vector so that it has the correct size and structure.

Parameters
vectorThe vector to be reinitialized

Implemented in DiFfRG::Variables::Assembler< Model_ >.

◆ residual() [1/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::residual ( VectorType & residual,
const VectorType & solution_global,
NumberType weight,
const VectorType & solution_global_dot,
NumberType weight_mass,
const VectorType & variables = VectorType() )
pure virtual

Calculates the residual for a DAE.

This function calculates the sum

\[\large \int_\Omega\, \bigg(w_m m_i(u_j, \partial_t u_j) \phi_i - w_r F_i(u_j, \dots) \partial_x \phi_i + w_r s_i(u_j, \dots) \phi_i \bigg) + \int_{\partial\Omega}\, w_r \widehat{F}_i(u_j, \dots) \phi_i \]

.

Here, \(w_r\) is the residual weight, \(w_m\,\,\) is the weight for the mass, and \(\widehat{F}_i\) are the boundary terms picked up by the partial integration of the flux.

Parameters
residualThe residual vector to be filled
solution_globalThe spatial solution vector
weightThe weight for the residual \( w_r \)
solution_global_dotThe spatial solution vector time derivative
weight_massThe weight for the mass \( w_m \)

◆ residual() [2/2]

template<typename VectorType , typename SparseMatrixType , uint dim>
void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::residual ( VectorType & residual,
const VectorType & solution_global,
NumberType weight,
NumberType weight_mass,
const VectorType & variables = VectorType() )
inline

Calculates the residual for an ODE.

This function calculates the sum

\[\large \int_\Omega\, \bigg(w_m m_ij \partial_t u_j \phi_i - w_r F_i(u_j, \dots) \partial_x \phi_i + w_r s_i(u_j, \dots) \phi_i \bigg) + \int_{\partial\Omega}\, w_r \widehat{F}_i(u_j, \dots) \phi_i \]

.

Here, \(w_r\) is the residual weight, \(w_m\,\,\) is the weight for the mass, and \(\widehat{F}_i\) are any boundary terms picked up by the partial integration of the flux.

Parameters
residualThe residual vector to be filled
solution_globalThe spatial solution vector
weightThe weight for the residual \( w_r \)
weight_massThe weight for the mass \( w_m \)

◆ residual_variables()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::residual_variables ( VectorType & residual,
const VectorType & variables,
const VectorType & spatial_solution )
inlinevirtual

When coupling the spatial discretization to additional variables, this function should calculate the residual for the additional variables.

Parameters
residualThe residual vector of the additional variables
variablesThe current additional variables vector
spatial_solutionThe spatial solution vector, which is needed for the calculation of the residual

Reimplemented in DiFfRG::Variables::Assembler< Model_ >.

◆ set_time()

template<typename VectorType , typename SparseMatrixType , uint dim>
virtual void DiFfRG::AbstractAssembler< VectorType, SparseMatrixType, dim >::set_time ( double t)
pure virtual

Set the current time. The assembler should usually just forward this to the numerical model.

Parameters
tThe current time

Implemented in DiFfRG::FEMAssembler< Discretization_, Model_ >, DiFfRG::FV::KurganovTadmor::Assembler< Discretization_, Model_ >, DiFfRG::LDG::LDGAssemblerBase< Discretization_, Model_ >, and DiFfRG::Variables::Assembler< Model_ >.


The documentation for this class was generated from the following file: