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::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t > Class Template Reference

Local solver for the nonlinear eigenvalue problem induced by the diffusion equation. More...

#include <diffusion_eigs_ldgh.hxx>

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

 DiffusionEigs (const constructor_value_type &tau=1.)
 Constructor for local solver. More...
 
template<typename hyEdgeT , typename SmallMatInT , typename SmallMatOutT >
SmallMatOutT & trace_to_flux (const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, hyEdgeT &hyper_edge, const lSol_float_t eig=0.) const
 Evaluate local contribution to matrix–vector multiplication. More...
 
template<class hyEdgeT >
std::array< unsigned int, 2 *hyEdge_dimT > node_types (hyEdgeT &hyper_edge) const
 Fill array with 1 if the node is of Dirichlet type and 0 elsewhen. More...
 
template<class hyEdgeT , typename SmallMatT >
SmallMatT & make_initial (SmallMatT &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
 L2 project given analytical functio to skeletal space. More...
 
template<class hyEdgeT , typename SmallMatInT , typename SmallMatOutT >
SmallMatOutT & jacobian_of_trace_to_flux (const SmallMatInT &lambda_values_in, SmallMatOutT &lambda_values_out, const lSol_float_t eig, const SmallMatInT &lambda_vals, const lSol_float_t eig_val, hyEdgeT &hyper_edge) const
 Evaluate local contribution to matrix–vector multiplication. More...
 
template<class hyEdgeT >
std::array< lSol_float_t, 1U > errors (const std::array< std::array< lSol_float_t, n_shape_bdr_ >, 2 *hyEdge_dimT > &lambda_values, hyEdgeT &hy_edge, const lSol_float_t eig=0.) const
 Local squared contribution to the L2 error. More...
 
template<typename abscissa_float_t , std::size_t abscissas_sizeT, class input_array_t , class hyEdgeT >
std::array< std::array< lSol_float_t, Hypercube< hyEdge_dimT >::pow(abscissas_sizeT)>, system_dimbulk_values (const std::array< abscissa_float_t, abscissas_sizeT > &abscissas, const input_array_t &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t time=0.) const
 Evaluate local 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 hyEdgeT >
SmallSquareMat< n_loc_dofs_, lSol_float_t > assemble_loc_matrix (const lSol_float_t tau, hyEdgeT &hyper_edge, const lSol_float_t eig) const
 Assemble local matrix for the local solver. More...
 
template<typename hyEdgeT , typename SmallMatT >
SmallVec< n_loc_dofs_, lSol_float_t > assemble_rhs_from_lambda (const SmallMatT &lambda_values, hyEdgeT &hyper_edge) const
 Assemble local right-hand for the local solver (from skeletal). More...
 
template<typename hyEdgeT , typename SmallMatT >
std::array< lSol_float_t, n_loc_dofs_solve_local_problem (const SmallMatT &lambda_values, hyEdgeT &hyper_edge, const lSol_float_t eig) const
 Solve local problem (with right-hand side from skeletal). More...
 
template<typename hyEdgeT >
std::array< std::array< lSol_float_t, n_shape_bdr_ >, 2 *hyEdge_dimT > primal_at_boundary (const std::array< lSol_float_t, n_loc_dofs_ > &coeffs, hyEdgeT &hyper_edge) const
 Evaluate primal variable at boundary. More...
 
template<typename hyEdgeT >
std::array< std::array< lSol_float_t, n_shape_bdr_ >, 2 *hyEdge_dimT > dual_at_boundary (const std::array< lSol_float_t,(hyEdge_dimT+1) *n_shape_fct_ > &coeffs, hyEdgeT &hyper_edge) const
 Evaluate dual variable at boundary. More...
 

Static Private Member Functions

template<typename parameters >
static constexpr bool is_dirichlet (const unsigned int node_type)
 Find out whether a node is of Dirichlet type. More...
 

Private Attributes

const lSol_float_t tau_
 (Globally constant) penalty parameter for HDG scheme. More...
 

Static Private Attributes

static constexpr unsigned int n_shape_fct_ = n_glob_dofs_per_node() * (poly_deg + 1)
 Number of local shape functions (with respect to all spatial dimensions). More...
 
static constexpr unsigned int n_shape_bdr_ = n_glob_dofs_per_node()
 Number oflocal shape functions (with respect to a face / hypernode). More...
 
static constexpr unsigned int n_loc_dofs_ = (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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
class LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >

Local solver for the nonlinear eigenvalue problem induced by the diffusion equation.


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.
space_dimTThe dimension of the surrounding space.
poly_degThe polynomial degree of test and trial functions.
quad_degThe order of the quadrature rule.
parametersTStruct depending on templates space_dimTP and lSol_float_TP that contains static parameter functions. Defaults to above functions included in DiffusionParametersDefault.
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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
typedef lSol_float_t LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, 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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, 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::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::integrator
private

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


Constructor & Destructor Documentation

◆ DiffusionEigs()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::DiffusionEigs ( 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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename hyEdgeT >
SmallSquareMat<n_loc_dofs_, lSol_float_t> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::assemble_loc_matrix ( const lSol_float_t  tau,
hyEdgeT &  hyper_edge,
const lSol_float_t  eig 
) const
inlineprivate

Assemble local matrix for the local solver.


Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
Parameters
tauPenalty parameter.
hyper_edgeThe geometry of the considered hyperedge (of typename GeomT).
eigEigenvalue for which the matrix is assembled.
Return values
loc_matMatrix of the local solver.

◆ assemble_rhs_from_lambda()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename hyEdgeT , typename SmallMatT >
SmallVec<n_loc_dofs_, lSol_float_t> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::assemble_rhs_from_lambda ( const SmallMatT &  lambda_values,
hyEdgeT &  hyper_edge 
) const
inlineprivate

Assemble local right-hand for the local solver (from skeletal).


The right hand side needs the values of the global degrees of freedom. Note that we basically follow the lines of

B. Cockburn, J. Gopalakrishnan, and R. Lazarov. Unified hybridization of discontinuous Galerkin, mixed, and continuous Galerkin methods for second order elliptic problems. SIAM Journal on Numerical Analysis, 47(2):1319–1365, 2009

and discriminate between local solvers with respect to the skeletal variable and with respect to the global right-hand side. This assembles the local right-hand side with respect to the skeletal variable.

Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
SmallMatTThe data type of the lambda_values.
Parameters
lambda_valuesGlobal degrees of freedom associated to the hyperedge.
hyper_edgeThe geometry of the considered hyperedge (of typename GeomT).
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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename abscissa_float_t , std::size_t abscissas_sizeT, class input_array_t , class hyEdgeT >
std::array<std::array<lSol_float_t, Hypercube<hyEdge_dimT>::pow(abscissas_sizeT)>, system_dim> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::bulk_values ( const std::array< abscissa_float_t, abscissas_sizeT > &  abscissas,
const input_array_t &  lambda_values,
hyEdgeT &  hyper_edge,
const lSol_float_t  time = 0. 
) const

Evaluate local 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.
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
Parameters
abscissasAbscissas of the supporting points.
lambda_valuesThe values of the skeletal variable's coefficients.
hyper_edgeThe geometry of the considered hyperedge (of typename GeomT).
timeTime.
Return values
func_valuesFunction values at tensorial points.

◆ dual_at_boundary()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename hyEdgeT >
std::array<std::array<lSol_float_t, n_shape_bdr_>, 2 * hyEdge_dimT> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::dual_at_boundary ( const std::array< lSol_float_t,(hyEdge_dimT+1) *n_shape_fct_ > &  coeffs,
hyEdgeT &  hyper_edge 
) 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.

Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
Parameters
coeffsCoefficients of the local solution.
hyper_edgeThe geometry of the considered hyperedge (of typename GeomT).
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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<class hyEdgeT >
std::array<lSol_float_t, 1U> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::errors ( const std::array< std::array< lSol_float_t, n_shape_bdr_ >, 2 *hyEdge_dimT > &  lambda_values,
hyEdgeT &  hy_edge,
const lSol_float_t  eig = 0. 
) const
inline

Local squared contribution to the L2 error.


Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
Parameters
lambda_valuesThe values of the skeletal variable's coefficients.
hy_edgeThe geometry of the considered hyperedge (of typename GeomT).
eigApproximated eigenvalue.
Return values
vec_bLocal part of vector b.
Here is the call graph for this function:

◆ hyEdge_dim()

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

Dimension of hyper edge type that this object solves on.


◆ is_dirichlet()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename parameters >
static constexpr bool LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::is_dirichlet ( const unsigned int  node_type)
inlinestaticconstexprprivate

Find out whether a node is of Dirichlet type.


◆ jacobian_of_trace_to_flux()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<class hyEdgeT , typename SmallMatInT , typename SmallMatOutT >
SmallMatOutT& LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::jacobian_of_trace_to_flux ( const SmallMatInT &  lambda_values_in,
SmallMatOutT &  lambda_values_out,
const lSol_float_t  eig,
const SmallMatInT &  lambda_vals,
const lSol_float_t  eig_val,
hyEdgeT &  hyper_edge 
) const
inline

Evaluate local contribution to matrix–vector multiplication.


This corresponds to the evaluation of the Jacobian evaluated at lambda_vals, eig_val in the direction lambda_values, eig.

Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
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.
eigEigenvalue in whose direction Jacobian is evaluated.
lambda_valsLambda at which Jacobian is evaluated.
eig_valEigenvalue at which Jacobian is evaluated.
hyper_edgeThe geometry of the considered hyperedge (of typename hyEdgeT).
Return values
vecAxLocal part of vector A * x.
Here is the call graph for this function:

◆ make_initial()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<class hyEdgeT , typename SmallMatT >
SmallMatT& LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::make_initial ( SmallMatT &  lambda_values,
hyEdgeT &  hyper_edge,
const lSol_float_t  time = 0. 
) const
inline

L2 project given analytical functio to skeletal space.


Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
SmallMatTThe data tyepe of lambda_values.
Parameters
lambda_valuesLocal part of vector x. (Redundant)
hyper_edgeThe geometry of the considered hyperedge (of typename hyEdgeT).
timeTime. (Redundant)
Return values
bdr_valuesCoefficients of L2 projection.

◆ n_glob_dofs_per_node()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
static constexpr unsigned int LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, 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.
Here is the call graph for this function:

◆ node_system_dimension()

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

Dimension of of the solution evaluated with respect to a hypernode.


◆ node_types()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<class hyEdgeT >
std::array<unsigned int, 2 * hyEdge_dimT> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::node_types ( hyEdgeT &  hyper_edge) const
inline

Fill array with 1 if the node is of Dirichlet type and 0 elsewhen.


Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
Parameters
hyper_edgeThe geometry of the considered hyperedge (of typename hyEdgeT).
Return values
resultArray encoding edge types.

◆ primal_at_boundary()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename hyEdgeT >
std::array<std::array<lSol_float_t, n_shape_bdr_>, 2 * hyEdge_dimT> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::primal_at_boundary ( const std::array< lSol_float_t, n_loc_dofs_ > &  coeffs,
hyEdgeT &  hyper_edge 
) 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.

Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
Parameters
coeffsCoefficients of the local solution.
hyper_edgeThe geometry of the considered hyperedge (of typename GeomT).
Return values
bdr_coeffsCoefficients of respective (dim-1) dimensional function at boundaries.

◆ solve_local_problem()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename hyEdgeT , typename SmallMatT >
std::array<lSol_float_t, n_loc_dofs_> LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::solve_local_problem ( const SmallMatT &  lambda_values,
hyEdgeT &  hyper_edge,
const lSol_float_t  eig 
) const
inlineprivate

Solve local problem (with right-hand side from skeletal).


Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
SmallMatTThe data type of the lambda_values.
Parameters
lambda_valuesGlobal degrees of freedom associated to the hyperedge.
hyper_edgeThe geometry of the considered hyperedge (of typename GeomT).
eigEigenvalue for which the local problem is solved.
Return values
loc_solSolution of the local problem.
Here is the call graph for this function:

◆ system_dimension()

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
static constexpr unsigned int LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, 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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
template<typename hyEdgeT , typename SmallMatInT , typename SmallMatOutT >
SmallMatOutT& LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::trace_to_flux ( const SmallMatInT &  lambda_values_in,
SmallMatOutT &  lambda_values_out,
hyEdgeT &  hyper_edge,
const lSol_float_t  eig = 0. 
) const
inline

Evaluate local contribution to matrix–vector multiplication.


This corresponds to an evaluation of the residual in the nonlinear context

Template Parameters
hyEdgeTThe geometry type / typename of the considered hyEdge's geometry.
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.
hyper_edgeThe geometry of the considered hyperedge (of typename hyEdgeT).
eigEigenvalue.
Return values
vecAxLocal part of vector A * x.
Here is the call graph for this function:

Member Data Documentation

◆ n_loc_dofs_

template<unsigned int hyEdge_dimT, unsigned int poly_deg, unsigned int quad_deg, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
constexpr unsigned int LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::n_loc_dofs_ = (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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
constexpr unsigned int LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::n_shape_bdr_ = n_glob_dofs_per_node()
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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
constexpr unsigned int LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::n_shape_fct_ = n_glob_dofs_per_node() * (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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
constexpr unsigned int LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, 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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
constexpr unsigned int LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, 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, template< unsigned int, typename > typename parametersT = DiffusionParametersDefault, typename lSol_float_t = double>
const lSol_float_t LocalSolver::DiffusionEigs< hyEdge_dimT, poly_deg, quad_deg, parametersT, lSol_float_t >::tau_
private

(Globally constant) penalty parameter for HDG scheme.



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