Local solver for bilaplacian equation on uniform hypergraph.
More...
#include <bilaplacian_uniform_ldgh.hxx>
|
| BilaplacianUniform (const constructor_value_type &tau=1.) |
| Constructor for local solver. More...
|
|
template<typename SmallMatInT , typename SmallMatOutT > |
SmallMatOutT & | trace_to_flux (const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, const lSol_float_t UNUSED(time)=0.) const |
| Evaluate local contribution to matrix–vector multiplication. More...
|
|
template<typename SmallMatInT , typename SmallMatOutT > |
SmallMatOutT & | residual_flux (const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, const lSol_float_t UNUSED(time)=0.) const |
| Evaluate local contribution to matrix–vector residual. More...
|
|
template<typename SmallMatT > |
std::array< lSol_float_t, 1U > | errors (const SmallMatT &UNUSED(lambda_values), const lSol_float_t UNUSED(time)=0.) const |
| Evaluate local squared L2 error. More...
|
|
template<typename abscissa_float_t , std::size_t sizeT, class input_array_t > |
std::array< std::array< lSol_float_t, Hypercube< hyEdge_dimT >::pow(sizeT)>, system_dim > | bulk_values (const std::array< abscissa_float_t, sizeT > &abscissas, const input_array_t &lambda_values, const lSol_float_t UNUSED(time)=0.) const |
| Evaluate local reconstruction at tensorial products of abscissas. More...
|
|
|
static constexpr unsigned int | hyEdge_dim () |
| Dimension of hyper edge type that this object solves on. More...
|
|
static constexpr unsigned int | n_glob_dofs_per_node () |
| Evaluate amount of global degrees of freedom per hypernode. More...
|
|
static constexpr unsigned int | system_dimension () |
| Dimension of of the solution evaluated with respect to a hyperedge. More...
|
|
static constexpr unsigned int | node_system_dimension () |
| Dimension of of the solution evaluated with respect to a hypernode. More...
|
|
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
class LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >
Local solver for bilaplacian equation on uniform hypergraph.
- Template Parameters
-
hyEdge_dimT | Dimension of a hyperedge, i.e., 1 is for PDEs defined on graphs, 2 is for PDEs defined on surfaces, and 3 is for PDEs defined on volumes. |
poly_deg | The polynomial degree of test and trial functions. |
quad_deg | The order of the quadrature rule. |
lSol_float_t | The floating point type calculations are executed in. Defaults to double. |
- Authors
- Guido Kanschat, Heidelberg University, 2019–2020.
-
Andreas Rupp, Heidelberg University, 2019–2020.
◆ constructor_value_type
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Class is constructed using a single double indicating the penalty parameter.
◆ integrator
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
An integrator helps to easily evaluate integrals (e.g. via quadrature).
◆ BilaplacianUniform()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Constructor for local solver.
- Parameters
-
tau | Penalty parameter of HDG scheme. |
◆ assemble_loc_matrix()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Assemble local matrix for the local solver.
The local solver neither depends on the geometry, nor on global functions. Thus, its local matrix is the same for all hyperedges and can be assembled once in the constructor. This is done in this function.
The function is static, since it is used in the constructor's initializer list.
- Parameters
-
tau | Penalty parameter for HDG. |
- Return values
-
loc_mat | Matrix of the local solver. |
◆ assemble_rhs()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
template<typename SmallMatT >
Assemble local right hand for the local solver.
The right hand side needs the values of the global degrees of freedom. Thus, it needs to be constructed individually for every hyperedge.
- Template Parameters
-
SmallMatT | The cata type of the lambda_values . |
- Parameters
-
lambda_values | Global degrees of freedom associated to the hyperedge. |
- Return values
-
loc_rhs | Local right hand side of the locasl solver. |
◆ bulk_values()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
template<typename abscissa_float_t , std::size_t sizeT, class input_array_t >
std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(sizeT)>, system_dim> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::bulk_values |
( |
const std::array< abscissa_float_t, sizeT > & |
abscissas, |
|
|
const input_array_t & |
lambda_values, |
|
|
const lSol_float_t |
UNUSEDtime = 0. |
|
) |
| const |
Evaluate local reconstruction at tensorial products of abscissas.
- Template Parameters
-
absc_float_t | Floating type for the abscissa values. |
abscissas_sizeT | Size of the array of array of abscissas. |
input_array_t | Type of input array. |
- Parameters
-
abscissas | Abscissas of the supporting points. |
lambda_values | The values of the skeletal variable's coefficients. |
time | Time — this parameter is redundant for this local solver. |
- Return values
-
function_values | Function values at quadrature points. |
◆ dual_at_boundary()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Evaluate dual variable at boundary.
Function to evaluate dual variable of the solution. This function is needed to calculate the local numerical fluxes.
- Parameters
-
coeffs | Coefficients of the local solution. |
- Return values
-
bdr_coeffs | Coefficients of respective (dim-1) dimensional function at boundaries. |
◆ errors()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
template<typename SmallMatT >
std::array<lSol_float_t, 1U> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::errors |
( |
const SmallMatT & |
UNUSEDlambda_values, |
|
|
const lSol_float_t |
UNUSEDtime = 0. |
|
) |
| const |
|
inline |
Evaluate local squared L2 error.
- Parameters
-
lambda_values | The values of the skeletal variable's coefficients. |
time | Time — this parameter is redundant for this local solver. |
- Return values
-
vec_b | Local part of vector b. |
◆ hyEdge_dim()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Dimension of hyper edge type that this object solves on.
◆ n_glob_dofs_per_node()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
◆ node_system_dimension()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Dimension of of the solution evaluated with respect to a hypernode.
◆ primal_at_boundary()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Evaluate primal variable at boundary.
Function to evaluate primal variable of the solution. This function is needed to calculate the local numerical fluxes.
- Parameters
-
coeffs | Coefficients of the local solution. |
- Return values
-
bdr_coeffs | Coefficients of respective (dim-1) dimensional function at boundaries. |
◆ residual_flux()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
template<typename SmallMatInT , typename SmallMatOutT >
SmallMatOutT& LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::residual_flux |
( |
const SmallMatInT & |
lambda_values_in, |
|
|
SmallMatOutT & |
lambda_values_out, |
|
|
const lSol_float_t |
UNUSEDtime = 0. |
|
) |
| const |
|
inline |
Evaluate local contribution to matrix–vector residual.
- Template Parameters
-
SmallMatInT | Data type of lambda_values_in . |
SmallMatOutT | Data type of lambda_values_out . |
- Parameters
-
lambda_values_in | Local part of vector x. |
lambda_values_out | Local part that will be added to A * x. |
time | Time — this parameter is redundant for this local solver. |
- Return values
-
vecAx | Local part of vector A * x. |
◆ solve_local_problem()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
template<typename SmallMatT >
Solve local problem.
- Template Parameters
-
SmallMatT | The cata type of the lambda_values . |
- Parameters
-
lambda_values | Global degrees of freedom associated to the hyperedge. |
- Return values
-
loc_sol | Solution of the local problem. |
◆ system_dimension()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Dimension of of the solution evaluated with respect to a hyperedge.
◆ trace_to_flux()
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
template<typename SmallMatInT , typename SmallMatOutT >
SmallMatOutT& LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::trace_to_flux |
( |
const SmallMatInT & |
lambda_values_in, |
|
|
SmallMatOutT & |
lambda_values_out, |
|
|
const lSol_float_t |
UNUSEDtime = 0. |
|
) |
| const |
|
inline |
Evaluate local contribution to matrix–vector multiplication.
Execute matrix–vector multiplication y = A * x, where x represents the vector containing the skeletal variable (adjacent to one hyperedge), and A is the condensed matrix arising from the HDG discretization. This function does this multiplication (locally) for one hyperedge. The hyperedge is no parameter, since all hyperedges are assumed to have the same properties.
- Template Parameters
-
SmallMatInT | Data type of lambda_values_in . |
SmallMatOutT | Data type of lambda_values_out . |
- Parameters
-
lambda_values_in | Local part of vector x. |
lambda_values_out | Local part that will be added to A * x. |
time | Time — this parameter is redundant for this local solver. |
- Return values
-
vecAx | Local part of vector A * x. |
◆ loc_mat_
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Local matrix for the local solver.
◆ n_loc_dofs_
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Number of (local) degrees of freedom per hyperedge.
◆ n_shape_bdr_
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Number oflocal shape functions (with respect to a face / hypernode).
◆ n_shape_fct_
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Number of local shape functions (with respect to all spatial dimensions).
◆ node_system_dim
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Dimension of of the solution evaluated with respect to a hypernode.
This allows to the use of this quantity as template parameter in member functions.
◆ system_dim
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
Dimension of of the solution evaluated with respect to a hypernode.
This allows to the use of this quantity as template parameter in member functions.
◆ tau_
template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
(Globally constant) penalty parameter for HDG scheme.
The documentation for this class was generated from the following file: