pymor.discretizers.builtin.grids.interfaces

Module Contents

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]
boundaries(boundary_type, codim)[source]
check_boundary_types(assert_unique_type=(1,), assert_some_type=())[source]
dirichlet_boundaries(codim)[source]
dirichlet_mask(codim)[source]
mask(boundary_type, codim)[source]

Return mask.

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

neumann_boundaries(codim)[source]
neumann_mask(codim)[source]
no_boundary_type_mask(codim)[source]

Return no boundary type mask.

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

robin_boundaries(codim)[source]
robin_mask(codim)[source]
unique_boundary_type_mask(codim)[source]

Return unique boundary type mask.

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

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

Bases: pymor.core.cache.CacheableObject

Affine grid.

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.

Methods

boundaries

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

boundary_mask

Return boundary mask.

bounding_box

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

centers

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

diameters

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

embeddings

Return embeddings.

integration_elements

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

jacobian_inverse_transposed

Return the inverse of transposed Jacobian.

neighbors

Maps entity index and local neighbor index to global neighbor index.

quadrature_points

Returns the quadrature points.

reference_element

The ReferenceElement of the codim-codim entities.

size

The number of entities of codimension codim.

subentities

Return subentities.

superentities

Return superentities.

superentity_indices

Return indices of superentities.

unit_outer_normals

Return unit outer normals.

volumes

Return volumes.

volumes_inverse

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

MAX_DOMAIN_RATIO = 1000000.0[source]
MAX_DOMAIN_WIDTH = 10000000000.0[source]
MIN_DOMAIN_WIDTH = 1e-10[source]
cache_region = 'memory'[source]
boundaries(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.

boundary_mask(codim)[source]

Return boundary mask.

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.

bounding_box()[source]

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

centers(codim)[source]

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

diameters(codim)[source]

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

abstract embeddings(codim)[source]

Return embeddings.

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.

integration_elements(codim)[source]

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

jacobian_inverse_transposed(codim)[source]

Return the inverse of transposed Jacobian.

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

neighbors(codim, neighbor_codim, intersection_codim=None)[source]

Maps entity index and local neighbor index to global neighbor index.

retval[e,n] is the global index of the n-th codim-neighbor_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 == neighbor_codim and to min(codim, neighbor_codim) otherwise.

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

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

Returns the quadrature points.

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
abstract reference_element(codim)[source]

The ReferenceElement of the codim-codim entities.

abstract size(codim)[source]

The number of entities of codimension codim.

abstract subentities(codim, subentity_codim)[source]

Return subentities.

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)

superentities(codim, superentity_codim)[source]

Return superentities.

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).

superentity_indices(codim, superentity_codim)[source]

Return indices of superentities.

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].

unit_outer_normals()[source]

Return unit outer normals.

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

volumes(codim)[source]

Return volumes.

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

volumes_inverse(codim)[source]

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

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

Bases: Grid

Grid with an additional orthogonal_centers method.

Methods

orthogonal_centers

Return orthogonal centers.

abstract orthogonal_centers()[source]

Return orthogonal centers.

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.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

Methods

center

Coordinates of the barycenter.

mapped_diameter

Return the diameter after transformation.

quadrature

Return quadrature points and weights.

quadrature_info

Return quadrature information.

quadrature_types

size

Number of subentities of codimension codim.

sub_reference_element

Returns the reference element of the codim-codim subentities.

subentities

Return subentities.

subentity_embedding

Return subentity embedding.

unit_outer_normals

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

cache_region = 'memory'[source]
abstract center()[source]

Coordinates of the barycenter.

abstract mapped_diameter(A)[source]

Return the diameter after transformation.

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

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

Return quadrature points and weights.

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()[source]

Return quadrature information.

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()[source]
abstract size(codim)[source]

Number of subentities of codimension codim.

abstract sub_reference_element(codim)[source]

Returns the reference element of the codim-codim subentities.

abstract subentities(codim, subentity_codim)[source]

Return subentities.

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(subentity_codim)[source]

Return subentity embedding.

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.

abstract unit_outer_normals()[source]

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