pymor.discretizers.builtin.grids.interfaces¶
Module Contents¶
- class pymor.discretizers.builtin.grids.interfaces.BoundaryInfo[source]¶
Bases:
pymor.core.cache.CacheableObjectProvides 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.
Methods
Return mask.
Return no boundary type mask.
Return unique boundary type mask.
- mask(boundary_type, codim)[source]¶
Return mask.
retval[i] is
Trueif the codim-codimentity of global indexiis associated to the boundary typeboundary_type.
- class pymor.discretizers.builtin.grids.interfaces.Grid[source]¶
Bases:
pymor.core.cache.CacheableObjectAffine 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
subentitiesand the embeddings given byembeddings. In addition, onlysizeandreference_elementhave to be implemented.Methods
Returns the global indices of all codim-
codimboundary entities.Return boundary mask.
Returns a
(2, dim)-shaped array containing lower/upper bounding box coordinates.retval[e]is the barycenter of the codim-codimentity with global indexe.retval[e]is the diameter of the codim-codimentity with global indexe.Return embeddings.
retval[e]is given assqrt(det(A^T*A)), whereA = embeddings(codim)[0][e].Return the inverse of transposed Jacobian.
Maps entity index and local neighbor index to global neighbor index.
Returns the quadrature points.
The
ReferenceElementof the codim-codimentities.The number of entities of codimension
codim.Return subentities.
Return superentities.
Return indices of superentities.
Return unit outer normals.
Return volumes.
retval[e] = 1 / volumes(codim)[e].- boundaries(codim)[source]¶
Returns the global indices of all codim-
codimboundary entities.By definition, a codim-1 entity is a boundary entity if it has only one codim-0 superentity. For
codim != 1, a codim-codimentity 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-codimentity with global indexeis 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-codimentity 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.
- abstract embeddings(codim)[source]¶
Return embeddings.
Returns tuple
(A, B)whereA[e]andB[e]are the linear part and the translation part of the map from the reference element ofetoe.For
codim > 0, we provide a default implementation by taking the embedding of the codim-1 parent entitye_0ofewith lowest global index and composing it with the subentity_embedding ofeintoe_0determined by the reference element.
- integration_elements(codim)[source]¶
retval[e]is given assqrt(det(A^T*A)), whereA = 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 ofembeddings(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 then-th codim-neighbor_codimentity of the codim-codimentityethat shares withea subentity of codimensionintersection_codim.If
intersection_codim == None, it is set tocodim + 1ifcodim == neighbor_codimand tomin(codim, neighbor_codim)otherwise.The default implementation is to compute the result from
subentities(codim, intersection_codim)andsuperentities(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-codimentity with global indexe.The quadrature is of order
orderor hasnpointsintegration points. To integrate a functionfovereone has to formnp.dot(f(quadrature_points(codim, order)[e]), reference_element(codim).quadrature(order)[1]) * integration_elements(codim)[e]. # NOQA
- abstract reference_element(codim)[source]¶
The
ReferenceElementof the codim-codimentities.
- abstract subentities(codim, subentity_codim)[source]¶
Return subentities.
retval[e,s]is the global index of thes-th codim-subentity_codimsubentity of the codim-codimentity 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 ofeas follows:Find the
codim-1parent entitye_0ofewith minimal global indexLookup the local indices of the subentities of
einsidee_0using 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(codim, superentity_codim)[source]¶
Return superentities.
retval[e,s]is the global index of thes-th codim-superentity_codimsuperentity of the codim-codimentity 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(codim, superentity_codim)[source]¶
Return indices of superentities.
retval[e,s]is the local index of the codim-codimentityein the codim-superentity_codimsuperentitysuperentities(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 indexe.
- class pymor.discretizers.builtin.grids.interfaces.GridWithOrthogonalCenters[source]¶
Bases:
GridGridwith an additionalorthogonal_centersmethod.Methods
Return orthogonal centers.
- abstract orthogonal_centers()[source]¶
Return orthogonal centers.
retval[e]is a point inside the codim-0 entity with global indexesuch that the line segment fromretval[e]toretval[e2]is always orthogonal to the codim-1 entity shared by the codim-0 entities with global indexeande2.(This is mainly useful for gradient approximation in finite volume schemes.)
- class pymor.discretizers.builtin.grids.interfaces.ReferenceElement[source]¶
Bases:
pymor.core.cache.CacheableObjectDefines 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.
Methods
Coordinates of the barycenter.
Return the diameter after transformation.
Return quadrature points and weights.
Return quadrature information.
Number of subentities of codimension
codim.Returns the reference element of the codim-
codimsubentities.Return subentities.
Return subentity embedding.
retval[e]is the unit outer-normal vector to the codim-1 subentity with indexe.- 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)wherePis an array of quadrature points with corresponding weightsW.The quadrature is of order
orderor hasnpointsintegration points.
- abstract quadrature_info()[source]¶
Return quadrature information.
Returns a tuple of dicts
(O, N)whereO[quadrature_type]is a list of orders which are implemented forquadrature_typeandN[quadrature_type]is a list of the corresponding numbers of integration points.
- abstract sub_reference_element(codim)[source]¶
Returns the reference element of the codim-
codimsubentities.
- 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 thej-th codim-subentity_codimsubentity of thei-th codim-codimsubentity 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_codimsubentities 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.