pymor.discretizers.builtin.grids.interfaces
¶
Module Contents¶
Classes¶
Defines a reference element. |
|
Topological grid with geometry where each codim-0 entity is affinely mapped to the same |
|
|
|
Provides boundary types for the boundaries of a given |
- 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.
- 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 thej
-th codim-subentity_codim
subentity of thei
-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)
andsub_reference_element(subentity_codim - 1).subentity_embedding(1)
choosing always the superentity with smallest index.
- abstract sub_reference_element(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 indexe
.
- 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)
whereP
is an array of quadrature points with corresponding weightsW
.The quadrature is of order
order
or hasnpoints
integration points.
- 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 byembeddings
. In addition, onlysize
andreference_element
have to be implemented.- abstract subentities(self, codim, subentity_codim)[source]¶
retval[e,s]
is the global index of thes
-th codim-subentity_codim
subentity of the codim-codim
entity with global indexe
.The ordering of
subentities(0, subentity_codim)[e]
has to correspond, w.r.t. the embedding ofe
, to the local ordering inside the reference element.For
codim > 0
, we provide a default implementation by calculating the subentities ofe
as follows:Find the
codim-1
parent entitye_0
ofe
with minimal global indexLookup the local indices of the subentities of
e
insidee_0
using the reference element.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 bye_0
, which agrees with what is returned byembeddings(codim)
- superentities(self, codim, superentity_codim)[source]¶
retval[e,s]
is the global index of thes
-th codim-superentity_codim
superentity of the codim-codim
entity with global indexe
.retval[e]
is sorted by global index.The default implementation is to compute the result from
subentities(superentity_codim, codim)
.
- superentity_indices(self, codim, superentity_codim)[source]¶
retval[e,s]
is the local index of the codim-codim
entitye
in the codim-superentity_codim
superentitysuperentities(codim, superentity_codim)[e,s].
- 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 then
-th codim-neighbour_codim
entity of the codim-codim
entitye
that shares withe
a subentity of codimensionintersection_codim
.If
intersection_codim == None
, it is set tocodim + 1
ifcodim == neighbour_codim
and tomin(codim, neighbour_codim)
otherwise.The default implementation is to compute the result from
subentities(codim, intersection_codim)
andsuperentities(intersection_codim, neihbour_codim)
.
- boundary_mask(self, codim)[source]¶
retval[e]
is true iff the codim-codim
entity with global indexe
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.
- 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.
- abstract reference_element(self, codim)[source]¶
The
ReferenceElement
of the codim-codim
entities.
- abstract embeddings(self, codim)[source]¶
Returns tuple
(A, B)
whereA[e]
andB[e]
are the linear part and the translation part of the map from the reference element ofe
toe
.For
codim > 0
, we provide a default implementation by taking the embedding of the codim-1 parent entitye_0
ofe
with lowest global index and composing it with the subentity_embedding ofe
intoe_0
determined by the reference element.
- jacobian_inverse_transposed(self, codim)[source]¶
retval[e]
is the transposed (pseudo-)inverse of the Jacobian ofembeddings(codim)[e]
.
- integration_elements(self, codim)[source]¶
retval[e]
is given assqrt(det(A^T*A))
, whereA = embeddings(codim)[0][e]
.
- volumes(self, codim)[source]¶
retval[e]
is the (dim-codim
)-dimensional volume of the codim-codim
entity with global indexe
.
- 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 indexe
.
- centers(self, codim)[source]¶
retval[e]
is the barycenter of the codim-codim
entity with global indexe
.
- diameters(self, codim)[source]¶
retval[e]
is the diameter of the codim-codim
entity with global indexe
.
- 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 indexe
.The quadrature is of order
order
or hasnpoints
integration points. To integrate a functionf
overe
one has to formnp.dot(f(quadrature_points(codim, order)[e]), reference_element(codim).quadrature(order)[1]) * integration_elements(codim)[e]. # NOQA
- class pymor.discretizers.builtin.grids.interfaces.GridWithOrthogonalCenters[source]¶
Bases:
Grid
Grid
with an additionalorthogonal_centers
method.- abstract orthogonal_centers(self)[source]¶
retval[e]
is a point inside the codim-0 entity with global indexe
such that the line segment fromretval[e]
toretval[e2]
is always orthogonal to the codim-1 entity shared by the codim-0 entities with global indexe
ande2
.(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.
- mask(self, boundary_type, codim)[source]¶
retval[i] is
True
if the codim-codim
entity of global indexi
is associated to the boundary typeboundary_type
.
- unique_boundary_type_mask(self, codim)[source]¶
retval[i] is
True
if the codim-codim
entity of global indexi
is associated to one and only one boundary type.