Metric terms and transformations¶

Coordinate transformations¶

grudge.geometry.forward_metric_nth_derivative(actx: ArrayContext, dcoll: DiscretizationCollection, xyz_axis: int, ref_axes: int | Tuple[Tuple[int, int], ...], dd: DOFDesc | None = None, *, _use_geoderiv_connection=False) DOFArray[source]¶

Pointwise metric derivatives representing repeated derivatives of the physical coordinate enumerated by xyz_axis: \(x_{\mathrm{xyz\_axis}}\) with respect to the coordiantes on the reference element \(\xi_i\):

\[D^\alpha x_{\mathrm{xyz\_axis}} = \frac{\partial^{|\alpha|} x_{\mathrm{xyz\_axis}} }{ \partial \xi_1^{\alpha_1}\cdots \partial \xi_m^{\alpha_m}}\]

where \(\alpha\) is a multi-index described by ref_axes.

Parameters:
  • xyz_axis – an integer denoting which physical coordinate to differentiate.

  • ref_axes – a tuple of tuples indicating indices of coordinate axes of the reference element to the number of derivatives which will be taken. For example, the value ((0, 2), (1, 1)) indicates taking the second derivative with respect to the first axis and the first derivative with respect to the second axis. Each axis must occur only once and the tuple must be sorted by the axis index.

  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – If True, process returned DOFArrays through _base_to_geoderiv_connection(). This should be set based on whether the code using the result of this function is able to make use of these arrays. (This is an internal argument and not intended for use outside grudge.)

Returns:

a DOFArray containing the pointwise metric derivative at each nodal coordinate.

grudge.geometry.forward_metric_derivative_mat(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None, *, _use_geoderiv_connection=False) ndarray[source]¶

Computes the forward metric derivative matrix, also commonly called the Jacobian matrix, with entries defined as the forward metric derivatives:

\[J = \left\lbrack \frac{\partial x_i}{\partial \xi_j} \right\rbrack_{(0, 0) \leq (i, j) \leq (n, m)}\]

where \(x_1, \dots, x_n\) denote the physical coordinates and \(\xi_1, \dots, \xi_m\) denote coordinates on the reference element. Note that, in the case of immersed manifolds, J is not necessarily a square matrix.

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

a matrix containing the evaluated forward metric derivatives of each physical coordinate, with respect to each reference coordinate.

grudge.geometry.inverse_metric_derivative_mat(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None, *, _use_geoderiv_connection=False) ndarray[source]¶

Computes the inverse metric derivative matrix, which is the inverse of the Jacobian (forward metric derivative) matrix.

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

a matrix containing the evaluated inverse metric derivatives.

grudge.geometry.first_fundamental_form(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None, *, _use_geoderiv_connection=False) ndarray[source]¶

Computes the first fundamental form using the Jacobian matrix:

\[\begin{split}\begin{bmatrix} E & F \\ F & G \end{bmatrix} := \begin{bmatrix} (\partial_u x)^2 & \partial_u x \partial_v x \\ \partial_u x \partial_v x & (\partial_v x)^2 \end{bmatrix} = J^T \cdot J\end{split}\]

where \(u, v\) are coordinates on the parameterized surface and \(x(u, v)\) defines a parameterized region. Here, \(J\) is the corresponding Jacobian matrix.

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

a matrix containing coefficients of the first fundamental form.

grudge.geometry.inverse_first_fundamental_form(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None, *, _use_geoderiv_connection=False) ndarray[source]¶

Computes the inverse of the first fundamental form:

\[\begin{split}\begin{bmatrix} E & F \\ F & G \end{bmatrix}^{-1} = \frac{1}{E G - F^2} \begin{bmatrix} G & -F \\ -F & E \end{bmatrix}\end{split}\]

where \(E, F, G\) are coefficients of the first fundamental form.

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

a matrix containing coefficients of the inverse of the first fundamental form.

Geometry terms¶

grudge.geometry.inverse_surface_metric_derivative(actx: ArrayContext, dcoll: DiscretizationCollection, rst_axis, xyz_axis, dd: DOFDesc | None = None, *, _use_geoderiv_connection=False)[source]¶

Computes the inverse surface metric derivative of the physical coordinate enumerated by xyz_axis with respect to the reference axis rst_axis. These geometric terms are used in the transformation of physical gradients.

This function does not cache its results.

Parameters:
  • rst_axis – an integer denoting the reference coordinate axis.

  • xyz_axis – an integer denoting the physical coordinate axis.

  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

a DOFArray containing the inverse metric derivative at each nodal coordinate.

grudge.geometry.inverse_surface_metric_derivative_mat(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None, *, times_area_element=False, _use_geoderiv_connection=False)[source]¶

Computes the matrix of inverse surface metric derivatives, indexed by (xyz_axis, rst_axis). It returns all values of inverse_surface_metric_derivative_mat() in cached matrix form.

This function caches its results.

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • times_area_element – If True, each entry of the matrix is premultiplied with the value of area_element(), reflecting the typical use of the matrix in integrals evaluating weak derivatives.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

a DOFArray containing the inverse metric derivatives in per-group arrays of shape (xyz_dimension, rst_dimension, nelements, ndof).

grudge.geometry.pseudoscalar(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None, *, _use_geoderiv_connection=False) MultiVector[source]¶

Computes the field of pseudoscalars for the domain/discretization identified by dd.

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

A MultiVector of DOFArrays.

grudge.geometry.area_element(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None, *, _use_geoderiv_connection=False) DOFArray[source]¶

Computes the scale factor used to transform integrals from reference to global space.

This function caches its results.

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization.

  • _use_geoderiv_connection – For internal use. See forward_metric_nth_derivative() for an explanation.

Returns:

a DOFArray containing the transformed volumes for each element.

Normal vectors¶

grudge.geometry.mv_normal(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc, *, _use_geoderiv_connection=False) MultiVector[source]¶

Exterior unit normal as a MultiVector. This supports both volume discretizations (where ambient == topological dimension) and surface discretizations (where ambient == topological dimension + 1). In the latter case, extra processing ensures that the returned normal is in the local tangent space of the element at the point where the normal is being evaluated.

This function caches its results.

Parameters:
  • dd – a DOFDesc as the surface discretization.

  • _use_geoderiv_connection – If True, process returned DOFArrays through _base_to_geoderiv_connection(). This should be set based on whether the code using the result of this function is able to make use of these arrays. The default value agrees with actx.supports_nonscalar_broadcasting. (This is an internal argument and not intended for use outside grudge.)

Returns:

a MultiVector containing the unit normals.

grudge.geometry.normal(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc, *, _use_geoderiv_connection=None)[source]¶

Get the unit normal to the specified surface discretization, dd. This supports both volume discretizations (where ambient == topological dimension) and surface discretizations (where ambient == topological dimension + 1). In the latter case, extra processing ensures that the returned normal is in the local tangent space of the element at the point where the normal is being evaluated.

This function may be treated as if it caches its results.

Parameters:
  • dd – a DOFDesc as the surface discretization.

  • _use_geoderiv_connection – See mv_normal() for a full description. (This is an internal argument and not intended for use outside grudge.)

Returns:

an object array of DOFArray containing the unit normals at each nodal location.

Curvature tensors¶

grudge.geometry.second_fundamental_form(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None) ndarray[source]¶

Computes the second fundamental form:

\[\begin{split}S(x) = \begin{bmatrix} \partial_{uu} x\cdot n & \partial_{uv} x\cdot n \\ \partial_{uv} x\cdot n & \partial_{vv} x\cdot n \end{bmatrix}\end{split}\]

where \(n\) is the surface normal, \(x(u, v)\) defines a parameterized surface, and \(u,v\) are coordinates on the parameterized surface.

Parameters:

dd – a DOFDesc, or a value convertible to one.

Returns:

a rank-2 object array describing second fundamental form.

grudge.geometry.shape_operator(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None) ndarray[source]¶

Computes the shape operator (also called the curvature tensor) containing second order derivatives:

\[\begin{split}C(x) = \begin{bmatrix} \partial_{uu} x & \partial_{uv} x \\ \partial_{uv} x & \partial_{vv} x \end{bmatrix}\end{split}\]

where \(x(u, v)\) defines a parameterized surface, and \(u,v\) are coordinates on the parameterized surface.

Parameters:

dd – a DOFDesc, or a value convertible to one.

Returns:

a rank-2 object array describing the shape operator.

grudge.geometry.summed_curvature(actx: ArrayContext, dcoll: DiscretizationCollection, dd: DOFDesc | None = None) DOFArray[source]¶

Computes the sum of the principal curvatures:

\[\kappa = \operatorname{Trace}(C(x))\]

where \(x(u, v)\) defines a parameterized surface, \(u,v\) are coordinates on the parameterized surface, and \(C(x)\) is the shape operator.

Parameters:

dd – a DOFDesc, or a value convertible to one.

Returns:

a DOFArray containing the summed curvature at each nodal coordinate.