pymor.discretizers.builtin.grids.rect¶
Module Contents¶
Classes¶
Basic implementation of a rectangular |
- class pymor.discretizers.builtin.grids.rect.RectGrid(num_intervals=(2, 2), domain=([0, 0], [1, 1]), identify_left_right=False, identify_bottom_top=False)[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.
- subentities(self, 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)
- embeddings(self, 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.
- bounding_box(self)[source]¶
Returns a
(2, dim)-shaped array containing lower/upper bounding box coordinates.
- structured_to_global(self, 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.
- global_to_structured(self, 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)).
- vertex_coordinates(self, 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]])
- orthogonal_centers(self)[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.)
- visualize(self, 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