Discontinuous Galerkin operators#

Core DG routines#

Elementwise differentiation#

grudge.op.local_grad(dcoll: DiscretizationCollection, *args, nested=False) ArrayOrContainer[source]#

Return the element-local gradient of a function \(f\) represented by vec:

\[\nabla|_E f = \left( \partial_x|_E f, \partial_y|_E f, \partial_z|_E f \right)\]

May be called with (vec) or (dd_in, vec).

Parameters:
  • vec – a DOFArray or an ArrayContainer of them.

  • dd_in – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization if not provided.

  • nested – return nested object arrays instead of a single multidimensional array if vec is non-scalar.

Returns:

an object array (possibly nested) of DOFArrays or ArrayContainer of object arrays.

grudge.op.local_d_dx(dcoll: DiscretizationCollection, xyz_axis, *args) ArrayOrContainer[source]#

Return the element-local derivative along axis xyz_axis of a function \(f\) represented by vec:

\[\frac{\partial f}{\partial \lbrace x,y,z\rbrace}\Big|_E\]

May be called with (vec) or (dd, vec).

Parameters:
  • xyz_axis – an integer indicating the axis along which the derivative is taken.

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

  • vec – a DOFArray or an ArrayContainer of them.

Returns:

a DOFArray or an ArrayContainer of them.

grudge.op.local_div(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Return the element-local divergence of the vector function \(\mathbf{f}\) represented by vecs:

\[\nabla|_E \cdot \mathbf{f} = \sum_{i=1}^d \partial_{x_i}|_E \mathbf{f}_i\]

May be called with (vecs) or (dd, vecs).

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

  • vecs – an object array of DOFArrays or an ArrayContainer object with object array entries. The last axis of the array must have length matching the volume dimension.

Returns:

a DOFArray or an ArrayContainer of them.

Weak derivative operators#

grudge.op.weak_local_grad(dcoll: DiscretizationCollection, *args, nested=False) ArrayOrContainer[source]#

Return the element-local weak gradient of the volume function represented by vec.

May be called with (vec) or (dd_in, vec).

Specifically, the function returns an object array where the \(i\)-th component is the weak derivative with respect to the \(i\)-th coordinate of a scalar function \(f\). See weak_local_d_dx() for further information. For non-scalar \(f\), the function will return a nested object array containing the component-wise weak derivatives.

Parameters:
  • dd_in – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization if not provided.

  • vec – a DOFArray or an ArrayContainer of them.

  • nested – return nested object arrays instead of a single multidimensional array if vec is non-scalar

Returns:

an object array (possibly nested) of DOFArrays or ArrayContainer of object arrays.

grudge.op.weak_local_d_dx(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Return the element-local weak derivative along axis xyz_axis of the volume function represented by vec.

May be called with (xyz_axis, vec) or (dd_in, xyz_axis, vec).

Specifically, this function computes the volume contribution of the weak derivative in the \(i\)-th component (specified by xyz_axis) of a function \(f\), in each element \(E\), with respect to polynomial test functions \(\phi\):

\[\int_E \partial_i\phi\,f\,\mathrm{d}x \sim \mathbf{D}_{E,i}^T \mathbf{M}_{E}^T\mathbf{f}|_E,\]

where \(\mathbf{D}_{E,i}\) is the polynomial differentiation matrix on an \(E\) for the \(i\)-th spatial coordinate, \(\mathbf{M}_E\) is the elemental mass matrix (see mass() for more information), and \(\mathbf{f}|_E\) is a vector of coefficients for \(f\) on \(E\).

Parameters:
  • dd_in – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization if not provided.

  • xyz_axis – an integer indicating the axis along which the derivative is taken.

  • vec – a DOFArray or an ArrayContainer of them.

Returns:

a DOFArray or an ArrayContainer of them.

grudge.op.weak_local_div(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Return the element-local weak divergence of the vector volume function represented by vecs.

May be called with (vecs) or (dd, vecs).

Specifically, this function computes the volume contribution of the weak divergence of a vector function \(\mathbf{f}\), in each element \(E\), with respect to polynomial test functions \(\phi\):

\[\int_E \nabla \phi \cdot \mathbf{f}\,\mathrm{d}x \sim \sum_{i=1}^d \mathbf{D}_{E,i}^T \mathbf{M}_{E}^T\mathbf{f}_i|_E,\]

where \(\mathbf{D}_{E,i}\) is the polynomial differentiation matrix on an \(E\) for the \(i\)-th spatial coordinate, and \(\mathbf{M}_E\) is the elemental mass matrix (see mass() for more information).

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

  • vecs – an object array of DOFArrays or an ArrayContainer object with object array entries. The last axis of the array must have length matching the volume dimension.

Returns:

a DOFArray or an ArrayContainer like vec.

Mass, inverse mass, and face mass operators#

grudge.op.mass(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Return the action of the DG mass matrix on a vector (or vectors) of DOFArrays, vec. In the case of vec being an ArrayContainer, the mass operator is applied component-wise.

May be called with (vec) or (dd_in, vec).

Specifically, this function applies the mass matrix elementwise on a vector of coefficients \(\mathbf{f}\) via: \(\mathbf{M}_{E}\mathbf{f}|_E\), where

\[\left(\mathbf{M}_{E}\right)_{ij} = \int_E \phi_i \cdot \phi_j\,\mathrm{d}x,\]

where \(\phi_i\) are local polynomial basis functions on \(E\).

Parameters:
  • dd_in – a DOFDesc, or a value convertible to one. Defaults to the base volume discretization if not provided.

  • vec – a DOFArray or an ArrayContainer of them.

Returns:

a DOFArray or an ArrayContainer like vec.

grudge.op.inverse_mass(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Return the action of the DG mass matrix inverse on a vector (or vectors) of DOFArrays, vec. In the case of vec being an ArrayContainer, the inverse mass operator is applied component-wise.

For affine elements \(E\), the element-wise mass inverse is computed directly as the inverse of the (physical) mass matrix:

\[\left(\mathbf{M}_{J^e}\right)_{ij} = \int_{\widehat{E}} \widehat{\phi}_i\cdot\widehat{\phi}_j J^e \mathrm{d}\widehat{x},\]

where \(\widehat{\phi}_i\) are basis functions over the reference element \(\widehat{E}\), and \(J^e\) is the (constant) Jacobian scaling factor (see grudge.geometry.area_element()).

For non-affine \(E\), \(J^e\) is not constant. In this case, a weight-adjusted approximation is used instead following [Chan_2016]:

\[\mathbf{M}_{J^e}^{-1} \approx \widehat{\mathbf{M}}^{-1}\mathbf{M}_{1/J^e}\widehat{\mathbf{M}}^{-1},\]

where \(\widehat{\mathbf{M}}\) is the reference mass matrix on \(\widehat{E}\).

May be called with (vec) or (dd, vec).

Parameters:
  • vec – a DOFArray or an ArrayContainer of them.

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

Returns:

a DOFArray or an ArrayContainer like vec.

grudge.op.face_mass(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Return the action of the DG face mass matrix on a vector (or vectors) of DOFArrays, vec. In the case of vec being an arbitrary ArrayContainer, the face mass operator is applied component-wise.

May be called with (vec) or (dd_in, vec).

Specifically, this function applies the face mass matrix elementwise on a vector of coefficients \(\mathbf{f}\) as the sum of contributions for each face \(f \subset \partial E\):

\[\sum_{f=1}^{N_{\text{faces}} } \mathbf{M}_{f, E}\mathbf{f}|_f,\]

where

\[\left(\mathbf{M}_{f, E}\right)_{ij} = \int_{f \subset \partial E} \phi_i(s)\psi_j(s)\,\mathrm{d}s,\]

where \(\phi_i\) are (volume) polynomial basis functions on \(E\) evaluated on the face \(f\), and \(\psi_j\) are basis functions for a polynomial space defined on \(f\).

Parameters:
  • dd – a DOFDesc, or a value convertible to one. Defaults to the base "all_faces" discretization if not provided.

  • vec – a DOFArray or an ArrayContainer of them.

Returns:

a DOFArray or an ArrayContainer like vec.

Working around documentation tool awkwardness#

class grudge.op.TracePair#

See grudge.trace_pair.TracePair.

Trace Pairs#

Container class and auxiliary functionality#

class grudge.trace_pair.TracePair(dd: DOFDesc, *, interior: Array | ArrayContainer, exterior: Array | ArrayContainer)[source]#

A container class for data (both interior and exterior restrictions) on the boundaries of mesh elements.

dd#

an instance of grudge.dof_desc.DOFDesc describing the discretization on which int and ext live.

int#

A DOFArray or ArrayContainer of them representing the interior value to be used for the flux.

ext#

A DOFArray or ArrayContainer of them representing the exterior value to be used for the flux.

avg#

A DOFArray or ArrayContainer of them representing the average of the interior and exterior values.

diff#

A DOFArray or ArrayContainer of them representing the difference (exterior - interior) of the pair values.

__getattr__(name)[source]#

Return a new TracePair resulting from executing attribute lookup with name on int and ext.

__getitem__(index)[source]#

Return a new TracePair resulting from executing subscripting with index on int and ext.

__len__()[source]#

Return the total number of arrays associated with the int and ext restrictions of the TracePair. Note that both must be the same.

class grudge.op.project_tracepair(dcoll: DiscretizationCollection, new_dd: DOFDesc, tpair: TracePair)[source]#

Return a new TracePair project()‘ed onto new_dd.

class grudge.op.tracepair_with_discr_tag(dcoll, discr_tag, tpair: TracePair)[source]#

Return a new TracePair project()‘ed onto a DOFDesc with discretization tag discr_tag.

Boundary trace functions#

grudge.op.bdry_trace_pair(dcoll: DiscretizationCollection, dd: Any, interior, exterior) TracePair[source]#

Returns a trace pair defined on the exterior boundary. Input arguments are assumed to already be defined on the boundary denoted by dd. If the input arguments interior and exterior are ArrayContainer objects, they must both have the same internal structure.

Parameters:
  • dd – a DOFDesc, or a value convertible to one, which describes the boundary discretization.

  • interior – a DOFArray or an ArrayContainer of them that contains data already on the boundary representing the interior value to be used for the flux.

  • exterior – a DOFArray or an ArrayContainer of them that contains data that already lives on the boundary representing the exterior value to be used for the flux.

Returns:

a TracePair on the boundary.

grudge.op.bv_trace_pair(dcoll: DiscretizationCollection, dd: Any, interior, exterior) TracePair[source]#

Returns a trace pair defined on the exterior boundary. The interior argument is assumed to be defined on the volume discretization, and will therefore be restricted to the boundary dd prior to creating a TracePair. If the input arguments interior and exterior are ArrayContainer objects, they must both have the same internal structure.

Parameters:
  • dd – a DOFDesc, or a value convertible to one, which describes the boundary discretization.

  • interior – a DOFArray or an ArrayContainer that contains data defined in the volume, which will be restricted to the boundary denoted by dd. The result will be used as the interior value for the flux.

  • exterior – a DOFArray or an ArrayContainer that contains data that already lives on the boundary representing the exterior value to be used for the flux.

Returns:

a TracePair on the boundary.

Interior and cross-rank trace functions#

grudge.op.interior_trace_pairs(dcoll: DiscretizationCollection, vec, *, comm_tag: Hashable | None = None, volume_dd: DOFDesc | None = None) List[TracePair][source]#

Return a list of TracePair objects defined on the interior faces of dcoll and any faces connected to a parallel boundary.

Note that local_interior_trace_pair() provides the rank-local contributions if those are needed in isolation. Similarly, cross_rank_trace_pairs() provides only the trace pairs defined on cross-rank boundaries.

Parameters:
  • vec – a DOFArray or an ArrayContainer of them.

  • comm_tag – a hashable object used to match sent and received data across ranks. Communication will only match if both endpoints specify objects that compare equal. A generalization of MPI communication tags to arbitrary, potentially composite objects.

Returns:

a list of TracePair objects.

grudge.op.local_interior_trace_pair(dcoll: DiscretizationCollection, vec, *, volume_dd: DOFDesc | None = None) TracePair[source]#

Return a TracePair for the interior faces of dcoll with a discretization tag specified by discr_tag. This does not include interior faces on different MPI ranks.

Parameters:

vec – a DOFArray or an ArrayContainer of them.

For certain applications, it may be useful to distinguish between rank-local and cross-rank trace pairs. For example, avoiding unnecessary communication of derived quantities (i.e. temperature) on partition boundaries by computing them directly. Having the ability for user applications to distinguish between rank-local and cross-rank contributions can also help enable overlapping communication with computation. :returns: a TracePair object.

grudge.op.cross_rank_trace_pairs(dcoll: DiscretizationCollection, ary: Array | ArrayContainer, tag: Hashable = None, *, comm_tag: Hashable = None, volume_dd: DOFDesc | None = None) List[TracePair][source]#

Get a list of ary trace pairs for each partition boundary.

For each partition boundary, the field data values in ary are communicated to/from the neighboring partition. Presumably, this communication is MPI (but strictly speaking, may not be, and this routine is agnostic to the underlying communication).

For each face on each partition boundary, a TracePair is created with the locally, and remotely owned partition boundary face data as the internal, and external components, respectively. Each of the TracePair components are structured like ary.

If ary is a number, rather than a DOFArray or an ArrayContainer of them, it is assumed that the same number is being communicated on every rank.

Parameters:
  • ary – a DOFArray or an ArrayContainer of them.

  • comm_tag – a hashable object used to match sent and received data across ranks. Communication will only match if both endpoints specify objects that compare equal. A generalization of MPI communication tags to arbitrary, potentially composite objects.

Returns:

a list of TracePair objects.

Transfering data between discretizations#

Projections#

grudge.op.project(dcoll: DiscretizationCollection, src: ConvertibleToDOFDesc, tgt: ConvertibleToDOFDesc, vec) ArrayOrContainer[source]#

Project from one discretization to another, e.g. from the volume to the boundary, or from the base to the an overintegrated quadrature discretization.

Parameters:
Returns:

a DOFArray or an ArrayContainer like vec.

Reductions#

Nodal Reductions#

Note

In a distributed-memory setting, these reductions automatically reduce over all ranks involved, and return the same value on all ranks, in the manner of an MPI allreduce.

grudge.op.norm(dcoll: DiscretizationCollection, vec, p, dd=None) int | float | complex | generic[source]#

Return the vector p-norm of a function represented by its vector of degrees of freedom vec.

Parameters:
  • vec – a DOFArray or an ArrayContainer of them.

  • p – an integer denoting the order of the integral norm. Currently, only values of 2 or numpy.inf are supported.

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

Returns:

a nonegative scalar denoting the norm.

grudge.op.nodal_sum(dcoll: DiscretizationCollection, dd, vec) int | float | complex | generic[source]#

Return the nodal sum of a vector of degrees of freedom vec.

Parameters:
Returns:

a device scalar denoting the nodal sum.

grudge.op.nodal_min(dcoll: DiscretizationCollection, dd, vec, *, initial=None) int | float | complex | generic[source]#

Return the nodal minimum of a vector of degrees of freedom vec.

Parameters:
  • dd – a DOFDesc, or a value convertible to one.

  • vec – a DOFArray or an ArrayContainer of them.

  • initial – an optional initial value. Defaults to numpy.inf.

Returns:

a device scalar denoting the nodal minimum.

grudge.op.nodal_max(dcoll: DiscretizationCollection, dd, vec, *, initial=None) int | float | complex | generic[source]#

Return the nodal maximum of a vector of degrees of freedom vec.

Parameters:
  • dd – a DOFDesc, or a value convertible to one.

  • vec – a DOFArray or an ArrayContainer of them.

  • initial – an optional initial value. Defaults to -numpy.inf.

Returns:

a device scalar denoting the nodal maximum.

grudge.op.integral(dcoll: DiscretizationCollection, dd, vec) int | float | complex | generic[source]#

Numerically integrates a function represented by a DOFArray of degrees of freedom.

Parameters:
Returns:

a device scalar denoting the evaluated integral.

Rank-local reductions#

grudge.op.nodal_sum_loc(dcoll: DiscretizationCollection, dd, vec) int | float | complex | generic[source]#

Return the rank-local nodal sum of a vector of degrees of freedom vec.

Parameters:
Returns:

a scalar denoting the rank-local nodal sum.

grudge.op.nodal_min_loc(dcoll: DiscretizationCollection, dd, vec, *, initial=None) int | float | complex | generic[source]#

Return the rank-local nodal minimum of a vector of degrees of freedom vec.

Parameters:
  • dd – a DOFDesc, or a value convertible to one.

  • vec – a DOFArray or an ArrayContainer of them.

  • initial – an optional initial value. Defaults to numpy.inf.

Returns:

a scalar denoting the rank-local nodal minimum.

grudge.op.nodal_max_loc(dcoll: DiscretizationCollection, dd, vec, *, initial=None) int | float | complex | generic[source]#

Return the rank-local nodal maximum of a vector of degrees of freedom vec.

Parameters:
  • dd – a DOFDesc, or a value convertible to one.

  • vec – a DOFArray or an ArrayContainer.

  • initial – an optional initial value. Defaults to -numpy.inf.

Returns:

a scalar denoting the rank-local nodal maximum.

Elementwise reductions#

grudge.op.elementwise_sum(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Returns a vector of DOFs with all entries on each element set to the sum of DOFs on that element.

May be called with (vec) or (dd, vec).

The input vec can either be a DOFArray or an ArrayContainer with DOFArray entries. If the underlying array context (see arraycontext.ArrayContext) for vec supports nonscalar broadcasting, all DOFArray entries will contain a single value for each element. Otherwise, the entries will have the same number of degrees of freedom as vec, but set to the same value.

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

  • vec – a DOFArray or an ArrayContainer of them

Returns:

a DOFArray or an ArrayContainer like vec whose entries denote the element-wise sum of vec.

grudge.op.elementwise_max(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Returns a vector of DOFs with all entries on each element set to the maximum over all DOFs on that element.

May be called with (vec) or (dd, vec).

The input vec can either be a DOFArray or an ArrayContainer with DOFArray entries. If the underlying array context (see arraycontext.ArrayContext) for vec supports nonscalar broadcasting, all DOFArray entries will contain a single value for each element. Otherwise, the entries will have the same number of degrees of freedom as vec, but set to the same value.

Parameters:
Returns:

a DOFArray or an ArrayContainer like vec whose entries denote the element-wise max of vec.

grudge.op.elementwise_min(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Returns a vector of DOFs with all entries on each element set to the minimum over all DOFs on that element.

May be called with (vec) or (dd, vec).

The input vec can either be a DOFArray or an ArrayContainer with DOFArray entries. If the underlying array context (see arraycontext.ArrayContext) for vec supports nonscalar broadcasting, all DOFArray entries will contain a single value for each element. Otherwise, the entries will have the same number of degrees of freedom as vec, but set to the same value.

Parameters:
Returns:

a DOFArray or an ArrayContainer like vec whose entries denote the element-wise min of vec.

grudge.op.elementwise_integral(dcoll: DiscretizationCollection, *args) ArrayOrContainer[source]#

Numerically integrates a function represented by a DOFArray of degrees of freedom in each element of a discretization, given by dd.

May be called with (vec) or (dd, vec).

The input vec can either be a DOFArray or an ArrayContainer with DOFArray entries. If the underlying array context (see arraycontext.ArrayContext) for vec supports nonscalar broadcasting, all DOFArray entries will contain a single value for each element. Otherwise, the entries will have the same number of degrees of freedom as vec, but set to the same value.

Parameters:
Returns:

a DOFArray or an ArrayContainer like vec containing the elementwise integral if vec.