HyperHDG
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes | Static Private Attributes | List of all members
LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t > Class Template Reference

Local solver for bilaplacian equation on uniform hypergraph. More...

#include <bilaplacian_uniform_ldgh.hxx>

Collaboration diagram for LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >:
Collaboration graph
[legend]

Classes

struct  data_type
 Define type of (hyperedge related) data that is stored in HyDataContainer. More...
 
struct  error_def
 Define how errors are evaluated. More...
 
struct  node_element
 Define type of node elements, especially with respect to nodal shape functions. More...
 

Public Types

typedef lSol_float_t constructor_value_type
 Class is constructed using a single double indicating the penalty parameter. More...
 

Public Member Functions

 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_dimbulk_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 Public Member Functions

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...
 

Private Types

typedef TPP::Quadrature::Tensorial< TPP::Quadrature::GaussLegendre< quad_deg >, TPP::ShapeFunction< TPP::ShapeType::Tensorial< TPP::ShapeType::Legendre< poly_deg >, hyEdge_dimT > >, lSol_float_t > integrator
 An integrator helps to easily evaluate integrals (e.g. via quadrature). More...
 

Private Member Functions

template<typename SmallMatT >
SmallVec< n_loc_dofs_, lSol_float_t > assemble_rhs (const SmallMatT &lambda_values) const
 Assemble local right hand for the local solver. More...
 
template<typename SmallMatT >
std::array< lSol_float_t, n_loc_dofs_solve_local_problem (const SmallMatT &lambda_values) const
 Solve local problem. More...
 
std::array< std::array< lSol_float_t, 2 *n_shape_bdr_ >, 2 *hyEdge_dimT > primal_at_boundary (const std::array< lSol_float_t, n_loc_dofs_ > &coeffs) const
 Evaluate primal variable at boundary. More...
 
std::array< std::array< lSol_float_t, 2 *n_shape_bdr_ >, 2 *hyEdge_dimT > dual_at_boundary (const std::array< lSol_float_t, 2 *(hyEdge_dimT+1) *n_shape_fct_ > &coeffs) const
 Evaluate dual variable at boundary. More...
 

Static Private Member Functions

static SmallSquareMat< n_loc_dofs_, lSol_float_t > assemble_loc_matrix (const lSol_float_t tau)
 Assemble local matrix for the local solver. More...
 

Private Attributes

const lSol_float_t tau_
 (Globally constant) penalty parameter for HDG scheme. More...
 
const SmallSquareMat< n_loc_dofs_, lSol_float_t > loc_mat_
 Local matrix for the local solver. More...
 

Static Private Attributes

static constexpr unsigned int n_shape_fct_ = Hypercube<hyEdge_dimT>::pow(poly_deg + 1)
 Number of local shape functions (with respect to all spatial dimensions). More...
 
static constexpr unsigned int n_shape_bdr_ = Hypercube<hyEdge_dimT - 1>::pow(poly_deg + 1)
 Number oflocal shape functions (with respect to a face / hypernode). More...
 
static constexpr unsigned int n_loc_dofs_ = 2 * (hyEdge_dimT + 1) * n_shape_fct_
 Number of (local) degrees of freedom per hyperedge. More...
 
static constexpr unsigned int system_dim = system_dimension()
 Dimension of of the solution evaluated with respect to a hypernode. More...
 
static constexpr unsigned int node_system_dim = node_system_dimension()
 Dimension of of the solution evaluated with respect to a hypernode. More...
 

Detailed Description

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_dimTDimension 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_degThe polynomial degree of test and trial functions.
quad_degThe order of the quadrature rule.
lSol_float_tThe floating point type calculations are executed in. Defaults to double.
Authors
Guido Kanschat, Heidelberg University, 2019–2020.
Andreas Rupp, Heidelberg University, 2019–2020.

Member Typedef Documentation

◆ constructor_value_type

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
typedef lSol_float_t LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::constructor_value_type

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>
typedef TPP::Quadrature::Tensorial< TPP::Quadrature::GaussLegendre<quad_deg>, TPP::ShapeFunction<TPP::ShapeType::Tensorial<TPP::ShapeType::Legendre<poly_deg>, hyEdge_dimT> >, lSol_float_t> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::integrator
private

An integrator helps to easily evaluate integrals (e.g. via quadrature).


Constructor & Destructor Documentation

◆ BilaplacianUniform()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::BilaplacianUniform ( const constructor_value_type tau = 1.)
inline

Constructor for local solver.


Parameters
tauPenalty parameter of HDG scheme.

Member Function Documentation

◆ assemble_loc_matrix()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
static SmallSquareMat<n_loc_dofs_, lSol_float_t> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::assemble_loc_matrix ( const lSol_float_t  tau)
staticprivate

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
tauPenalty parameter for HDG.
Return values
loc_matMatrix 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 >
SmallVec<n_loc_dofs_, lSol_float_t> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::assemble_rhs ( const SmallMatT &  lambda_values) const
inlineprivate

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
SmallMatTThe cata type of the lambda_values.
Parameters
lambda_valuesGlobal degrees of freedom associated to the hyperedge.
Return values
loc_rhsLocal 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_tFloating type for the abscissa values.
abscissas_sizeTSize of the array of array of abscissas.
input_array_tType of input array.
Parameters
abscissasAbscissas of the supporting points.
lambda_valuesThe values of the skeletal variable's coefficients.
timeTime — this parameter is redundant for this local solver.
Return values
function_valuesFunction 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>
std::array<std::array<lSol_float_t, 2 * n_shape_bdr_>, 2 * hyEdge_dimT> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::dual_at_boundary ( const std::array< lSol_float_t, 2 *(hyEdge_dimT+1) *n_shape_fct_ > &  coeffs) const
inlineprivate

Evaluate dual variable at boundary.


Function to evaluate dual variable of the solution. This function is needed to calculate the local numerical fluxes.

Parameters
coeffsCoefficients of the local solution.
Return values
bdr_coeffsCoefficients 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_valuesThe values of the skeletal variable's coefficients.
timeTime — this parameter is redundant for this local solver.
Return values
vec_bLocal part of vector b.

◆ hyEdge_dim()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
static constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::hyEdge_dim ( )
inlinestaticconstexpr

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>
static constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::n_glob_dofs_per_node ( )
inlinestaticconstexpr

Evaluate amount of global degrees of freedom per hypernode.


This number must be equal to HyperNodeFactory::n_dofs_per_node() of the HyperNodeFactory cooperating with this object.

Return values
n_dofsNumber of global degrees of freedom per hypernode.

◆ node_system_dimension()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
static constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::node_system_dimension ( )
inlinestaticconstexpr

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>
std::array<std::array<lSol_float_t, 2 * n_shape_bdr_>, 2 * hyEdge_dimT> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::primal_at_boundary ( const std::array< lSol_float_t, n_loc_dofs_ > &  coeffs) const
inlineprivate

Evaluate primal variable at boundary.


Function to evaluate primal variable of the solution. This function is needed to calculate the local numerical fluxes.

Parameters
coeffsCoefficients of the local solution.
Return values
bdr_coeffsCoefficients 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
SmallMatInTData type of lambda_values_in.
SmallMatOutTData type of lambda_values_out.
Parameters
lambda_values_inLocal part of vector x.
lambda_values_outLocal part that will be added to A * x.
timeTime — this parameter is redundant for this local solver.
Return values
vecAxLocal 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 >
std::array<lSol_float_t, n_loc_dofs_> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::solve_local_problem ( const SmallMatT &  lambda_values) const
inlineprivate

Solve local problem.


Template Parameters
SmallMatTThe cata type of the lambda_values.
Parameters
lambda_valuesGlobal degrees of freedom associated to the hyperedge.
Return values
loc_solSolution of the local problem.

◆ system_dimension()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
static constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::system_dimension ( )
inlinestaticconstexpr

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
SmallMatInTData type of lambda_values_in.
SmallMatOutTData type of lambda_values_out.
Parameters
lambda_values_inLocal part of vector x.
lambda_values_outLocal part that will be added to A * x.
timeTime — this parameter is redundant for this local solver.
Return values
vecAxLocal part of vector A * x.

Member Data Documentation

◆ loc_mat_

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, typename lSol_float_t = double>
const SmallSquareMat<n_loc_dofs_, lSol_float_t> LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::loc_mat_
private

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>
constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::n_loc_dofs_ = 2 * (hyEdge_dimT + 1) * n_shape_fct_
staticconstexprprivate

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>
constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::n_shape_bdr_ = Hypercube<hyEdge_dimT - 1>::pow(poly_deg + 1)
staticconstexprprivate

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>
constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::n_shape_fct_ = Hypercube<hyEdge_dimT>::pow(poly_deg + 1)
staticconstexprprivate

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>
constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::node_system_dim = node_system_dimension()
staticconstexprprivate

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>
constexpr unsigned int LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::system_dim = system_dimension()
staticconstexprprivate

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>
const lSol_float_t LocalSolver::BilaplacianUniform< hyEdge_dimT, poly_deg, quad_deg, lSol_float_t >::tau_
private

(Globally constant) penalty parameter for HDG scheme.



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