metrics

Additional metrics for evaluating a set of nodes

recursivenodes.metrics.mass_matrix_condition(d, n, nodes, domain='biunit')[source]

Compute the condition number of the mass matrix \(M_{ij} = \int_{\Delta^d} \varphi_i \varphi_j\ dx\), where \(\{\varphi_i\}\) is the set of shape functions associated with a set of nodes.

Parameters:
  • d (int) – The spatial dimension

  • n (int) – The polynomial degree

  • nodes (ndarray) – The set of polynomial nodes, with shape \((\binom{n+d}{d}, d)\) (unless domain=barycentric, then with shape \((\binom{n+d}{d}, d+1))\)

  • domain (str, optional) – The domain of the nodes. See “The domain optional argument”.

Returns:

the condition number \(\kappa_2(M)\).

Return type:

float

recursivenodes.metrics.nodal_gradient_condition(d, n, nodes, domain='biunit')[source]

Compute the condition number of the nodal gradient matrix \(G_{(d i + k)j} = \nabla_k \varphi_j(x_i)\), where \(\{\varphi_i\}\) is the set of shape functions associated with a set of nodes.

Parameters:
  • d (int) – The spatial dimension

  • n (int) – The polynomial degree

  • nodes (ndarray) – The set of polynomial nodes, with shape \((\binom{n+d}{d}, d)\) (unless domain=barycentric, then with shape \((\binom{n+d}{d}, d+1))\)

  • domain (str, optional) – The domain of the nodes. See “The domain optional argument”. The condition number will be computed for the nodes mapped to the biunit simplex.

Returns:

the condition number \(\kappa_2(G)\) (with respect to the pseudoinverse).

Return type:

float

recursivenodes.metrics.nodal_laplacian_condition(d, n, nodes, domain='biunit')[source]

Compute the condition number of the nodal strong-form Laplacian matrix \(L_{ij} = \Delta \varphi_j(x_i)\), where \(\{\varphi_i\}\) is the set of shape functions associated with a set of nodes.

Parameters:
  • d (int) – The spatial dimension

  • n (int) – The polynomial degree

  • nodes (ndarray) – The set of polynomial nodes, with shape \((\binom{n+d}{d}, d)\) (unless domain=barycentric, then with shape \((\binom{n+d}{d}, d+1))\)

  • domain (str, optional) – The domain of the nodes. See “The domain optional argument”. The condition number will be computed for the nodes mapped to the biunit simplex.

Returns:

the condition number \(\kappa_2(L)\) (with respect to the pseudoinverse).

Return type:

float

recursivenodes.metrics.strong_laplacian_condition(d, n, nodes, domain='biunit')[source]

Compute the condition number of the strong-form Laplacian matrix \(L_{ij} = \int_{\Delta^d} \varphi_i \cdot \Delta \varphi_j\ dx\), where \(\{\varphi_i\}\) is the set of shape functions associated with a set of nodes.

Parameters:
  • d (int) – The spatial dimension

  • n (int) – The polynomial degree

  • nodes (ndarray) – The set of polynomial nodes, with shape \((\binom{n+d}{d}, d)\) (unless domain=barycentric, then with shape \((\binom{n+d}{d}, d+1))\)

  • domain (str, optional) – The domain of the nodes. See “The domain optional argument”. The condition number will be computed for the nodes mapped to the biunit simplex.

Returns:

the condition number \(\kappa_2(L)\) (with respect to the pseudoinverse).

Return type:

float

recursivenodes.metrics.weak_laplacian_condition(d, n, nodes, domain='biunit')[source]

Compute the condition number of the weak-form Laplacian matrix \(K_{ij} = \int_{\Delta^d} \nabla \varphi_i \cdot \nabla \varphi_j\ dx\), where \(\{\varphi_i\}\) is the set of shape functions associated with a set of nodes.

Parameters:
  • d (int) – The spatial dimension

  • n (int) – The polynomial degree

  • nodes (ndarray) – The set of polynomial nodes, with shape \((\binom{n+d}{d}, d)\) (unless domain=barycentric, then with shape \((\binom{n+d}{d}, d+1)\))

  • domain (str, optional) – The domain of the nodes. See “The domain optional argument”. The condition number will be computed for the nodes mapped to the biunit simplex.

Returns:

the condition number \(\kappa_2(K)\) (with respect to the pseudoinverse).

Return type:

float