pymor.discretizers.builtin.grids.tria¶
Module Contents¶
Classes¶
Basic implementation of a triangular grid on a rectangular domain. |
- class pymor.discretizers.builtin.grids.tria.TriaGrid(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 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.
- 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.
- 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