pymor.discretizers.builtin.grids.interfaces

Module Contents

Classes

ReferenceElement

Defines a reference element.

Grid

Topological grid with geometry where each codim-0 entity is affinely mapped to the same

GridWithOrthogonalCenters

Grid with an additional orthogonal_centers method.

BoundaryInfo

Provides boundary types for the boundaries of a given Grid.

class pymor.discretizers.builtin.grids.interfaces.ReferenceElement[source]

Bases: pymor.core.cache.CacheableObject

Defines a reference element.

All reference elements have the property that all subentities of a given codimension are of the same type. I.e. a three-dimensional reference element cannot have triangles and rectangles as faces at the same time.

dim[source]

The dimension of the reference element

volume[source]

The volume of the reference element

cache_region = memory[source]
abstract size(self, codim)[source]

Number of subentities of codimension codim.

abstract subentities(self, codim, subentity_codim)[source]

subentities(c,sc)[i,j] is, with respect to the indexing inside the reference element, the index of the j-th codim-subentity_codim subentity of the i-th codim-codim subentity of the reference element.

abstract subentity_embedding(self, subentity_codim)[source]

Returns a tuple (A, B) which defines the embedding of the codim-subentity_codim subentities into the reference element.

For subentity_codim > 1', the embedding is by default given recursively via `subentity_embedding(subentity_codim - 1) and sub_reference_element(subentity_codim - 1).subentity_embedding(1) choosing always the superentity with smallest index.

_subentity_embedding(self, subentity_codim)[source]
abstract sub_reference_element(self, codim)[source]

Returns the reference element of the codim-codim subentities.

_sub_reference_element(self, codim)[source]
__call__(self, codim)[source]

Returns the reference element of the codim-codim subentities.

abstract unit_outer_normals(self)[source]

retval[e] is the unit outer-normal vector to the codim-1 subentity with index e.

abstract center(self)[source]

Coordinates of the barycenter.

abstract mapped_diameter(self, A)[source]

The diameter of the reference element after transforming it with the matrix A (vectorized).

abstract quadrature(self, order=None, npoints=None, quadrature_type='default')[source]

Returns tuple (P, W) where P is an array of quadrature points with corresponding weights W.

The quadrature is of order order or has npoints integration points.

abstract quadrature_info(self)[source]

Returns a tuple of dicts (O, N) where O[quadrature_type] is a list of orders which are implemented for quadrature_type and N[quadrature_type] is a list of the corresponding numbers of integration points.

quadrature_types(self)[source]
class pymor.discretizers.builtin.grids.interfaces.Grid[source]

Bases: pymor.core.cache.CacheableObject

Topological grid with geometry where each codim-0 entity is affinely mapped to the same ReferenceElement.

The grid is completely determined via the subentity relation given by subentities and the embeddings given by embeddings. In addition, only size and reference_element have to be implemented.

cache_region = memory[source]
MAX_DOMAIN_WIDTH = 1000000000000.0[source]
abstract size(self, codim)[source]

The number of entities of codimension codim.

abstract subentities(self, codim, subentity_codim)[source]

retval[e,s] is the global index of the s-th codim-subentity_codim subentity of the codim-codim entity with global index e.

The ordering of subentities(0, subentity_codim)[e] has to correspond, w.r.t. the embedding of e, to the local ordering inside the reference element.

For codim > 0, we provide a default implementation by calculating the subentities of e as follows:

  1. Find the codim-1 parent entity e_0 of e with minimal global index

  2. Lookup the local indices of the subentities of e inside e_0 using the reference element.

  3. Map these local indices to global indices using subentities(codim - 1, subentity_codim).

This procedures assures that subentities(codim, subentity_codim)[e] has the right ordering w.r.t. the embedding determined by e_0, which agrees with what is returned by embeddings(codim)

_subentities(self, codim, subentity_codim)[source]
superentities(self, codim, superentity_codim)[source]

retval[e,s] is the global index of the s-th codim-superentity_codim superentity of the codim-codim entity with global index e.

retval[e] is sorted by global index.

The default implementation is to compute the result from subentities(superentity_codim, codim).

_superentities(self, codim, superentity_codim)[source]
superentity_indices(self, codim, superentity_codim)[source]

retval[e,s] is the local index of the codim-codim entity e in the codim-superentity_codim superentity superentities(codim, superentity_codim)[e,s].

_superentity_indices(self, codim, superentity_codim)[source]
_superentities_with_indices(self, codim, superentity_codim)[source]
neighbours(self, codim, neighbour_codim, intersection_codim=None)[source]

Maps entity index and local neighbour index to global neighbour index

retval[e,n] is the global index of the n-th codim-neighbour_codim entity of the codim-codim entity e that shares with e a subentity of codimension intersection_codim.

If intersection_codim == None, it is set to codim + 1 if codim == neighbour_codim and to min(codim, neighbour_codim) otherwise.

The default implementation is to compute the result from subentities(codim, intersection_codim) and superentities(intersection_codim, neihbour_codim).

_neighbours(self, codim, neighbour_codim, intersection_codim)[source]
boundary_mask(self, codim)[source]

retval[e] is true iff the codim-codim entity with global index e is a boundary entity.

By definition, a codim-1 entity is a boundary entity if it has only one codim-0 superentity. For codim != 1, a codim-codim entity is a boundary entity if it has a codim-1 sub/super-entity.

_boundary_mask(self, codim)[source]
boundaries(self, codim)[source]

Returns the global indices of all codim-codim boundary entities.

By definition, a codim-1 entity is a boundary entity if it has only one codim-0 superentity. For codim != 1, a codim-codim entity is a boundary entity if it has a codim-1 sub/super-entity.

_boundaries(self, codim)[source]
abstract reference_element(self, codim)[source]

The ReferenceElement of the codim-codim entities.

abstract embeddings(self, codim)[source]

Returns tuple (A, B) where A[e] and B[e] are the linear part and the translation part of the map from the reference element of e to e.

For codim > 0, we provide a default implementation by taking the embedding of the codim-1 parent entity e_0 of e with lowest global index and composing it with the subentity_embedding of e into e_0 determined by the reference element.

_embeddings(self, codim)[source]
jacobian_inverse_transposed(self, codim)[source]

retval[e] is the transposed (pseudo-)inverse of the Jacobian of embeddings(codim)[e].

_jacobian_inverse_transposed(self, codim)[source]
integration_elements(self, codim)[source]

retval[e] is given as sqrt(det(A^T*A)), where A = embeddings(codim)[0][e].

_integration_elements(self, codim)[source]
volumes(self, codim)[source]

retval[e] is the (dim-codim)-dimensional volume of the codim-codim entity with global index e.

_volumes(self, codim)[source]
volumes_inverse(self, codim)[source]

retval[e] = 1 / volumes(codim)[e].

_volumes_inverse(self, codim)[source]
unit_outer_normals(self)[source]

retval[e,i] is the unit outer normal to the i-th codim-1 subentity of the codim-0 entity with global index e.

_unit_outer_normals(self)[source]
centers(self, codim)[source]

retval[e] is the barycenter of the codim-codim entity with global index e.

_centers(self, codim)[source]
diameters(self, codim)[source]

retval[e] is the diameter of the codim-codim entity with global index e.

_diameters(self, codim)[source]
quadrature_points(self, codim, order=None, npoints=None, quadrature_type='default')[source]

retval[e] is an array of quadrature points in global coordinates for the codim-codim entity with global index e.

The quadrature is of order order or has npoints integration points. To integrate a function f over e one has to form

np.dot(f(quadrature_points(codim, order)[e]),
    reference_element(codim).quadrature(order)[1]) *
integration_elements(codim)[e].  # NOQA
_quadrature_points(self, codim, order, npoints, quadrature_type)[source]
bounding_box(self)[source]

Returns a (2, dim)-shaped array containing lower/upper bounding box coordinates.

_bounding_box(self)[source]
classmethod _check_domain(cls, domain)[source]
class pymor.discretizers.builtin.grids.interfaces.GridWithOrthogonalCenters[source]

Bases: Grid

Grid with an additional orthogonal_centers method.

abstract orthogonal_centers(self)[source]

retval[e] is a point inside the codim-0 entity with global index e such that the line segment from retval[e] to retval[e2] is always orthogonal to the codim-1 entity shared by the codim-0 entities with global index e and e2.

(This is mainly useful for gradient approximation in finite volume schemes.)

class pymor.discretizers.builtin.grids.interfaces.BoundaryInfo[source]

Bases: pymor.core.cache.CacheableObject

Provides boundary types for the boundaries of a given Grid.

For every boundary type and codimension a mask is provided, marking grid entities of the respective type and codimension by their global index.

boundary_types[source]

set of all boundary types the grid has.

cache_region = memory[source]
mask(self, boundary_type, codim)[source]

retval[i] is True if the codim-codim entity of global index i is associated to the boundary type boundary_type.

unique_boundary_type_mask(self, codim)[source]

retval[i] is True if the codim-codim entity of global index i is associated to one and only one boundary type.

no_boundary_type_mask(self, codim)[source]

retval[i] is True if the codim-codim entity of global index i is associated to no boundary type.

check_boundary_types(self, assert_unique_type=(1,), assert_some_type=())[source]
property has_dirichlet(self)[source]
property has_neumann(self)[source]
property has_robin(self)[source]
dirichlet_mask(self, codim)[source]
neumann_mask(self, codim)[source]
robin_mask(self, codim)[source]
_boundaries(self, boundary_type, codim)[source]
boundaries(self, boundary_type, codim)[source]
dirichlet_boundaries(self, codim)[source]
neumann_boundaries(self, codim)[source]
robin_boundaries(self, codim)[source]