pymor.discretizers.builtin.grids package¶
Submodules¶
_unstructured module¶
boundaryinfos module¶
-
class
pymor.discretizers.builtin.grids.boundaryinfos.AllDirichletBoundaryInfo(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.BoundaryInfoBoundaryInfowhere the boundary type ‘dirichlet’ is attached to each boundary entity.Methods
boundaries,check_boundary_types,dirichlet_boundaries,dirichlet_mask,neumann_boundaries,neumann_mask,no_boundary_type_mask,robin_boundaries,robin_mask,unique_boundary_type_maskAttributes
boundary_types,cache_region,has_dirichlet,has_neumann,has_robin
-
class
pymor.discretizers.builtin.grids.boundaryinfos.EmptyBoundaryInfo(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.BoundaryInfoBoundaryInfowith no boundary types attached to any boundary.Methods
boundaries,check_boundary_types,dirichlet_boundaries,dirichlet_mask,neumann_boundaries,neumann_mask,no_boundary_type_mask,robin_boundaries,robin_mask,unique_boundary_type_maskAttributes
boundary_types,cache_region,has_dirichlet,has_neumann,has_robin
-
class
pymor.discretizers.builtin.grids.boundaryinfos.GenericBoundaryInfo(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.BoundaryInfoGeneric
BoundaryInfostoring entity masks per boundary type.Methods
boundaries,check_boundary_types,dirichlet_boundaries,dirichlet_mask,neumann_boundaries,neumann_mask,no_boundary_type_mask,robin_boundaries,robin_mask,unique_boundary_type_maskAttributes
boundary_types,cache_region,has_dirichlet,has_neumann,has_robin-
classmethod
from_indicators(grid, indicators, assert_unique_type=None, assert_some_type=None)[source]¶ Create
BoundaryInfofrom indicator functions.Parameters
- grid
The
Gridto which theBoundaryInfois associated.- indicators
Dict where each key is a boundary type and the corresponding value is a boolean valued function defined on the analytical domain which indicates if a point belongs to a boundary of the given boundary type (the indicator functions must be vectorized).
-
classmethod
constructions module¶
-
pymor.discretizers.builtin.grids.constructions.flatten_grid(grid)[source]¶ This method is used by our visualizers to render n-dimensional grids which cannot be embedded into R^n by duplicating vertices which would have to be mapped to multiple points at once (think of grids on rectangular domains with identified edges).
Parameters
- grid
The
Gridto flatten.
Returns
- subentities
The
subentities(0, grid.dim)relation for the flattened grid.- coordinates
The coordinates of the codim-
grid.dimentities.- entity_map
Maps the indices of the codim-
grid.dimentities of the flattened grid to the indices of the corresponding entities in the original grid.
gmsh module¶
-
pymor.discretizers.builtin.grids.gmsh.load_gmsh(filename)[source]¶ Parse a Gmsh file and create a corresponding
GmshGridandGmshBoundaryInfo.Parameters
- filename
Path of the Gmsh MSH-file.
Returns
- grid
The generated
GmshGrid.- boundary_info
The generated
GmshBoundaryInfo.
interfaces module¶
-
class
pymor.discretizers.builtin.grids.interfaces.BoundaryInfo(*args, **kwargs)[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
boundaries,check_boundary_types,dirichlet_boundaries,dirichlet_mask,mask,neumann_boundaries,neumann_mask,no_boundary_type_mask,robin_boundaries,robin_mask,unique_boundary_type_maskAttributes
boundary_types,cache_region,has_dirichlet,has_neumann,has_robin-
boundary_types¶ set of all boundary types the grid has.
-
mask(boundary_type, codim)[source]¶ retval[i] is
Trueif the codim-codimentity of global indexiis associated to the boundary typeboundary_type.
-
-
class
pymor.discretizers.builtin.grids.interfaces.Grid(*args, **kwargs)[source]¶ Bases:
pymor.core.cache.CacheableObjectTopological 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
Attributes
-
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]¶ 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]¶ 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]¶ retval[e]is the transposed (pseudo-)inverse of the Jacobian ofembeddings(codim)[e].
-
neighbours(codim, neighbour_codim, intersection_codim=None)[source]¶ retval[e,n]is the global index of then-th codim-neighbour_codimentitiy of the codim-codimentityethat shares withea subentity of codimensionintersection_codim.If
intersection_codim == None, it is set tocodim + 1ifcodim == neighbour_codimand tomin(codim, neighbour_codim)otherwise.The default implementation is to compute the result from
subentities(codim, intersection_codim)andsuperentities(intersection_codim, neihbour_codim).
-
quadrature_points(codim, order=None, npoints=None, quadrature_type='default')[source]¶ 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]¶ 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]¶ 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]¶ retval[e,s]is the local index of the codim-codimentityein the codim-superentity_codimsuperentitysuperentities(codim, superentity_codim)[e,s].
-
unit_outer_normals()[source]¶ retval[e,i]is the unit outer normal to the i-th codim-1 subentity of the codim-0 entitiy with global indexe.
-
-
class
pymor.discretizers.builtin.grids.interfaces.GridWithOrthogonalCenters(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.GridGridwith an additionalorthogonal_centersmethod.Methods
Attributes
-
abstract
orthogonal_centers()[source]¶ 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.)
-
abstract
-
class
pymor.discretizers.builtin.grids.interfaces.ReferenceElement(*args, **kwargs)[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
center,mapped_diameter,quadrature,quadrature_info,quadrature_types,size,sub_reference_element,subentities,subentity_embedding,unit_outer_normals,__call__Attributes
-
dim¶ The dimension of the reference element
-
volume¶ The volume of the reference element
-
abstract
mapped_diameter(A)[source]¶ The diameter of the reference element after transforming it with the matrix
A(vectorized).
-
abstract
quadrature(order=None, npoints=None, quadrature_type='default')[source]¶ 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]¶ 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]¶ 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]¶ 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.
-
oned module¶
-
class
pymor.discretizers.builtin.grids.oned.OnedGrid(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.GridWithOrthogonalCentersOne-dimensional
Gridon an interval.Parameters
- domain
Tuple
(left, right)containing the left and right boundary of the domain.- num_intervals
The number of codim-0 entities.
Methods
Attributes
dim,reference_element-
bounding_box()[source]¶ returns a
(2, dim)-shaped array containing lower/upper bounding box coordinates.
-
embeddings(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 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.
-
orthogonal_centers()[source]¶ 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.)
-
subentities(codim, subentity_codim)[source]¶ 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)
-
visualize(U, codim=2, **kwargs)[source]¶ Visualize scalar data associated to the grid as a patch plot.
Parameters
- U
NumPy arrayof the data to visualize. IfU.dim == 2 and len(U) > 1, the data is visualized as a time series of plots. Alternatively, a tuple ofNumPy arrayscan be provided, in which case a subplot is created for each entry of the tuple. The lengths of all arrays have to agree.- codim
The codimension of the entities the data in
Uis attached to (either 0 or 2).- kwargs
See
visualize_patch
rect module¶
-
class
pymor.discretizers.builtin.grids.rect.RectGrid(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.GridWithOrthogonalCentersBasic implementation of a rectangular
Gridon a rectangular domain.The global face, edge and vertex indices are given as follows
x1 ^ | 6--10---7--11---8 | | | 3 2 4 3 5 | | | 3---8---4---9---5 | | | 0 0 1 1 2 | | | 0---6---1---7---2 --> x0
Parameters
- num_intervals
Tuple
(n0, n1)determining a grid withn0xn1codim-0 entities.- domain
Tuple
(ll, ur)wherelldefines the lower left andurthe upper right corner of the domain.- identify_left_right
If
True, the left and right boundaries are identified, i.e. the left-most codim-0 entities become neighbors of the right-most codim-0 entities.- identify_bottom_top
If
True, the bottom and top boundaries are identified, i.e. the bottom-most codim-0 entities become neighbors of the top-most codim-0 entities.
Methods
Attributes
dim,reference_element-
bounding_box()[source]¶ returns a
(2, dim)-shaped array containing lower/upper bounding box coordinates.
-
embeddings(codim=0)[source]¶ 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.
-
global_to_structured(codim)[source]¶ Returns an array which maps global codim-
codimindices to structured indices.I.e. if
GTS = global_to_structured(codim)andSTG = structured_to_global(codim), thenSTG[GTS[:, 0], GTS[:, 1]] == numpy.arange(size(codim)).
-
orthogonal_centers()[source]¶ 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.)
-
structured_to_global(codim)[source]¶ Returns an
NumPy arraywhich maps structured indices to global codim-codimindices.In other words
structured_to_global(codim)[i, j]is the global index of the i-th in x0-direction and j-th in x1-direction codim-codimentity of the grid.
-
subentities(codim, subentity_codim)[source]¶ 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)
-
vertex_coordinates(dim)[source]¶ Returns an array of the x_dim coordinates of the grid vertices.
I.e.
centers(2)[structured_to_global(2)[i, j]] == np.array([vertex_coordinates(0)[i], vertex_coordinates(1)[j]])
-
visualize(U, codim=2, **kwargs)[source]¶ Visualize scalar data associated to the grid as a patch plot.
Parameters
- U
NumPy arrayof the data to visualize. IfU.dim == 2 and len(U) > 1, the data is visualized as a time series of plots. Alternatively, a tuple ofNumPy arrayscan be provided, in which case a subplot is created for each entry of the tuple. The lengths of all arrays have to agree.- codim
The codimension of the entities the data in
Uis attached to (either 0 or 2).- kwargs
See
visualize_patch
referenceelements module¶
-
class
pymor.discretizers.builtin.grids.referenceelements.Line(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.ReferenceElementMethods
center,mapped_diameter,quadrature,quadrature_info,size,sub_reference_element,subentities,subentity_embedding,unit_outer_normalsquadrature_types,__call__Attributes
-
mapped_diameter(A)[source]¶ The diameter of the reference element after transforming it with the matrix
A(vectorized).
-
quadrature(order=None, npoints=None, quadrature_type='default')[source]¶ Returns tuple
(P, W)wherePis an array of quadrature points with corresponding weightsW.The quadrature is of order
orderor hasnpointsintegration points.
-
quadrature_info()[source]¶ 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.
-
subentities(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_codimsubentity of thei-th codim-codimsubentity of the reference element.
-
subentity_embedding(subentity_codim)[source]¶ 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.
-
-
class
pymor.discretizers.builtin.grids.referenceelements.Point(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.ReferenceElementMethods
center,mapped_diameter,quadrature,quadrature_info,size,sub_reference_element,subentities,subentity_embedding,unit_outer_normalsquadrature_types,__call__Attributes
-
mapped_diameter(A)[source]¶ The diameter of the reference element after transforming it with the matrix
A(vectorized).
-
quadrature(order=None, npoints=None, quadrature_type='default')[source]¶ Returns tuple
(P, W)wherePis an array of quadrature points with corresponding weightsW.The quadrature is of order
orderor hasnpointsintegration points.
-
quadrature_info()[source]¶ 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.
-
subentities(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_codimsubentity of thei-th codim-codimsubentity of the reference element.
-
subentity_embedding(subentity_codim)[source]¶ 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.
-
-
class
pymor.discretizers.builtin.grids.referenceelements.Square(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.ReferenceElementMethods
center,mapped_diameter,quadrature,quadrature_info,size,sub_reference_element,subentities,subentity_embedding,unit_outer_normalsquadrature_types,__call__Attributes
-
mapped_diameter(A)[source]¶ The diameter of the reference element after transforming it with the matrix
A(vectorized).
-
quadrature(order=None, npoints=None, quadrature_type='default')[source]¶ Returns tuple
(P, W)wherePis an array of quadrature points with corresponding weightsW.The quadrature is of order
orderor hasnpointsintegration points.
-
quadrature_info()[source]¶ 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.
-
subentities(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_codimsubentity of thei-th codim-codimsubentity of the reference element.
-
subentity_embedding(subentity_codim)[source]¶ 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.
-
-
class
pymor.discretizers.builtin.grids.referenceelements.Triangle(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.ReferenceElementMethods
center,mapped_diameter,quadrature,quadrature_info,size,sub_reference_element,subentities,subentity_embedding,unit_outer_normalsquadrature_types,__call__Attributes
-
mapped_diameter(A)[source]¶ The diameter of the reference element after transforming it with the matrix
A(vectorized).
-
quadrature(order=None, npoints=None, quadrature_type='default')[source]¶ Returns tuple
(P, W)wherePis an array of quadrature points with corresponding weightsW.The quadrature is of order
orderor hasnpointsintegration points.
-
quadrature_info()[source]¶ 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.
-
subentities(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_codimsubentity of thei-th codim-codimsubentity of the reference element.
-
subentity_embedding(subentity_codim)[source]¶ 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.
-
subgrid module¶
-
class
pymor.discretizers.builtin.grids.subgrid.SubGrid(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.GridA subgrid of a
Grid.Given a
Gridand a list of codim-0 entities we construct the minimal subgrid of the grid, containing all the given entities.Parameters
- parent_grid
Gridof which a subgrid is to be created.- parent_entities
NumPy arrayof global indices of the codim-0 entities which are to be contained in the subgrid.
Methods
Attributes
parent_grid,reference_element-
parent_grid¶ The
Gridfrom which the subgrid was constructed.Subgridonly stores aweakrefto the grid, so accessing this property might returnNoneif the original grid has been destroyed.
-
embeddings(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 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.
-
indices_from_parent_indices(ind, codim)[source]¶ Maps a
NumPy arrayof indicies of codim-codimentites of the parent grid to indicies of the subgrid.Raises
- ValueError
Not all provided indices correspond to entities contained in the subgrid.
-
parent_indices(codim)[source]¶ retval[e]is the index of thee-th codim-codimentity in the parent grid.
-
subentities(codim, subentity_codim)[source]¶ 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)
-
pymor.discretizers.builtin.grids.subgrid.make_sub_grid_boundary_info(sub_grid, parent_grid, parent_grid_boundary_info, new_boundary_type=None)[source]¶ Derives a
BoundaryInfofor aSubGrid.Parameters
- sub_grid
The
SubGridfor which aBoundaryInfois created.- parent_grid
The parent
Grid.- parent_grid_boundary_info
The
BoundaryInfoof the parentGridfrom which to derive theBoundaryInfo- new_boundary_type
The boundary type which is assigned to the new boundaries of
subgrid. IfNone, no boundary type is assigned.
Returns
BoundaryInfoassociated with sub_grid.
tria module¶
-
class
pymor.discretizers.builtin.grids.tria.TriaGrid(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.GridWithOrthogonalCentersBasic implementation of a triangular grid on a rectangular domain.
The global face, edge and vertex indices are given as follows
6---------10----------7---------11----------8 | \ / | \ / | | 22 10 18 | 23 11 19 | | \ / | \ / | 3 14 11 6 4 15 12 7 5 | / \ | / \ | | 14 2 26 | 15 3 27 | | / \ | / \ | 3----------8----------4----------9----------5 | \ / | \ / | | 20 8 16 | 21 9 17 | | \ / | \ / | 0 12 9 4 1 13 10 5 2 | / \ | / \ | | 12 0 24 | 13 1 25 | | / \ | / \ | 0----------6----------1----------7----------2
Parameters
- num_intervals
Tuple
(n0, n1)determining a grid withn0xn1codim-0 entities.- domain
Tuple
(ll, ur)wherelldefines the lower left andurthe upper right corner of the domain.- identify_left_right
If
True, the left and right boundaries are identified, i.e. the left-most codim-0 entities become neighbors of the right-most codim-0 entities.- identify_bottom_top
If
True, the bottom and top boundaries are identified, i.e. the bottom-most codim-0 entities become neighbors of the top-most codim-0 entities.
Methods
Attributes
dim,reference_element-
bounding_box()[source]¶ returns a
(2, dim)-shaped array containing lower/upper bounding box coordinates.
-
embeddings(codim=0)[source]¶ 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.
-
orthogonal_centers()[source]¶ 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.)
-
subentities(codim, subentity_codim)[source]¶ 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)
-
visualize(U, codim=2, **kwargs)[source]¶ Visualize scalar data associated to the grid as a patch plot.
Parameters
- U
NumPy arrayof the data to visualize. IfU.dim == 2 and len(U) > 1, the data is visualized as a time series of plots. Alternatively, a tuple ofNumPy arrayscan be provided, in which case a subplot is created for each entry of the tuple. The lengths of all arrays have to agree.- codim
The codimension of the entities the data in
Uis attached to (either 0 or 2).- kwargs
See
visualize_patch
unstructured module¶
-
class
pymor.discretizers.builtin.grids.unstructured.UnstructuredTriangleGrid(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.grids.interfaces.GridA generic unstructured, triangular grid.
Use
from_verticesto instantiate the grid from vertex coordinates and connectivity data.Methods
Attributes
dim,reference_element-
embeddings(codim=0)[source]¶ 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.
-
classmethod
from_vertices(vertices, faces)[source]¶ Instantiate grid from vertex coordinates and connectivity data.
Parameters
- vertices
A (num_vertices, 2)-shaped
NumPy arraycontaining the coordinates of all vertices in the grid. The row numbers in the array will be the global indices of the given vertices (codim 2 entities).- faces
A (num_faces, 3)-shaped
NumPy arraycontaining the global indices of the vertices which define a given triangle in the grid. The row numbers in the array will be the global indices of the given triangles (codim 0 entities).
-
subentities(codim=0, subentity_codim=None)[source]¶ 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)
-
visualize(U, codim=2, **kwargs)[source]¶ Visualize scalar data associated to the grid as a patch plot.
Parameters
- U
NumPy arrayof the data to visualize. IfU.dim == 2 and len(U) > 1, the data is visualized as a time series of plots. Alternatively, a tuple ofNumPy arrayscan be provided, in which case a subplot is created for each entry of the tuple. The lengths of all arrays have to agree.- codim
The codimension of the entities the data in
Uis attached to (either 0 or 2).- kwargs
See
visualize_patch
-
vtkio module¶
-
pymor.discretizers.builtin.grids.vtkio.write_vtk(grid, data, filename_base, codim=2, binary_vtk=True, last_step=None)[source]¶ Output grid-associated data in (legacy) vtk format
Parameters
- grid
A
Gridwith triangular or rectilinear reference element.- data
VectorArraywith either cell (ie one datapoint per codim 0 entity) or vertex (ie one datapoint per codim 2 entity) data in each array element.- codim
the codimension associated with the data
- filename_base
common component for output files in timeseries
- binary_vtk
if false, output files contain human readable inline ascii data, else appended binary
- last_step
if set must be <= len(data) to restrict output of timeseries