pymor.discretizers.builtin¶
Subpackages¶
Submodules¶
Package Contents¶
Classes¶
One-dimensional |
|
Basic implementation of a rectangular |
|
Basic implementation of a triangular grid on a rectangular domain. |
Functions¶
Discretizes a |
|
Discretizes an |
|
Discretizes a |
|
FV Discretization of an |
|
Parse a Gmsh file and create a corresponding |
- pymor.discretizers.builtin.discretize_stationary_cg(analytical_problem, diameter=None, domain_discretizer=None, grid_type=None, grid=None, boundary_info=None, preassemble=True, mu_energy_product=None)[source]¶
Discretizes a
StationaryProblemusing finite elements.Parameters
- analytical_problem
The
StationaryProblemto discretize.- diameter
If not
None,diameteris passed as an argument to thedomain_discretizer.- domain_discretizer
Discretizer to be used for discretizing the analytical domain. This has to be a function
domain_discretizer(domain_description, diameter, ...). IfNone,discretize_domain_defaultis used.- grid_type
If not
None, this parameter is forwarded todomain_discretizerto specify the type of the generatedGrid.- grid
Instead of using a domain discretizer, the
Gridcan also be passed directly using this parameter.- boundary_info
A
BoundaryInfospecifying the boundary types of the grid boundary entities. Must be provided ifgridis specified.- preassemble
If
True, preassemble all operators in the resultingModel.- mu_energy_product
If not
None,parameter valuesfor which to assemble the symmetric part of theOperatorof the resultingModelfom(ignoring the advection part). Thus, assuming no advection and a symmetric diffusion tensor,fom.products['energy']is equal tofom.operator.assemble(mu), except for the fact that the former has cleared Dirichlet rows and columns, while the latter only has cleared Dirichlet rows).
Returns
- m
The
Modelthat has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid.- boundary_info
The generated
BoundaryInfo.- unassembled_m
In case
preassembleisTrue, the generatedModelbefore preassembling operators.
- pymor.discretizers.builtin.discretize_instationary_cg(analytical_problem, diameter=None, domain_discretizer=None, grid_type=None, grid=None, boundary_info=None, num_values=None, time_stepper=None, nt=None, preassemble=True)[source]¶
Discretizes an
InstationaryProblemwith aStationaryProblemas stationary part using finite elements.Parameters
- analytical_problem
The
InstationaryProblemto discretize.- diameter
If not
None,diameteris passed as an argument to thedomain_discretizer.- domain_discretizer
Discretizer to be used for discretizing the analytical domain. This has to be a function
domain_discretizer(domain_description, diameter, ...). IfNone,discretize_domain_defaultis used.- grid_type
If not
None, this parameter is forwarded todomain_discretizerto specify the type of the generatedGrid.- grid
Instead of using a domain discretizer, the
Gridcan also be passed directly using this parameter.- boundary_info
A
BoundaryInfospecifying the boundary types of the grid boundary entities. Must be provided ifgridis specified.- num_values
The number of returned vectors of the solution trajectory. If
None, each intermediate vector that is calculated is returned.- time_stepper
The
time-stepperto be used bysolve.- nt
If
time_stepperis not specified, the number of time steps for implicit Euler time stepping.- preassemble
If
True, preassemble all operators in the resultingModel.
Returns
- m
The
Modelthat has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid.- boundary_info
The generated
BoundaryInfo.- unassembled_m
In case
preassembleisTrue, the generatedModelbefore preassembling operators.
- pymor.discretizers.builtin.discretize_stationary_fv(analytical_problem, diameter=None, domain_discretizer=None, grid_type=None, num_flux='lax_friedrichs', lxf_lambda=1.0, eo_gausspoints=5, eo_intervals=1, grid=None, boundary_info=None, preassemble=True)[source]¶
Discretizes a
StationaryProblemusing the finite volume method.Parameters
- analytical_problem
The
StationaryProblemto discretize.- diameter
If not
None,diameteris passed as an argument to thedomain_discretizer.- domain_discretizer
Discretizer to be used for discretizing the analytical domain. This has to be a function
domain_discretizer(domain_description, diameter, ...). IfNone,discretize_domain_defaultis used.- grid_type
If not
None, this parameter is forwarded todomain_discretizerto specify the type of the generatedGrid.- num_flux
The numerical flux to use in the finite volume formulation. Allowed values are
'lax_friedrichs','engquist_osher','simplified_engquist_osher'(seepymor.discretizers.builtin.fv).- lxf_lambda
The stabilization parameter for the Lax-Friedrichs numerical flux (ignored, if different flux is chosen).
- eo_gausspoints
Number of Gauss points for the Engquist-Osher numerical flux (ignored, if different flux is chosen).
- eo_intervals
Number of sub-intervals to use for integration when using Engquist-Osher numerical flux (ignored, if different flux is chosen).
- grid
Instead of using a domain discretizer, the
Gridcan also be passed directly using this parameter.- boundary_info
A
BoundaryInfospecifying the boundary types of the grid boundary entities. Must be provided ifgridis specified.- preassemble
If
True, preassemble all operators in the resultingModel.
Returns
- m
The
Modelthat has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid.- boundary_info
The generated
BoundaryInfo.- unassembled_m
In case
preassembleisTrue, the generatedModelbefore preassembling operators.
- pymor.discretizers.builtin.discretize_instationary_fv(analytical_problem, diameter=None, domain_discretizer=None, grid_type=None, num_flux='lax_friedrichs', lxf_lambda=1.0, eo_gausspoints=5, eo_intervals=1, grid=None, boundary_info=None, num_values=None, time_stepper=None, nt=None, preassemble=True)[source]¶
FV Discretization of an
InstationaryProblemwith aStationaryProblemas stationary partParameters
- analytical_problem
The
InstationaryProblemto discretize.- diameter
If not
None,diameteris passed to thedomain_discretizer.- domain_discretizer
Discretizer to be used for discretizing the analytical domain. This has to be a function
domain_discretizer(domain_description, diameter, ...). If further arguments should be passed to the discretizer, usefunctools.partial. IfNone,discretize_domain_defaultis used.- grid_type
If not
None, this parameter is forwarded todomain_discretizerto specify the type of the generatedGrid.- num_flux
The numerical flux to use in the finite volume formulation. Allowed values are
'lax_friedrichs','engquist_osher','simplified_engquist_osher'(seepymor.discretizers.builtin.fv).- lxf_lambda
The stabilization parameter for the Lax-Friedrichs numerical flux (ignored, if different flux is chosen).
- eo_gausspoints
Number of Gauss points for the Engquist-Osher numerical flux (ignored, if different flux is chosen).
- eo_intervals
Number of sub-intervals to use for integration when using Engquist-Osher numerical flux (ignored, if different flux is chosen).
- grid
Instead of using a domain discretizer, the
Gridcan also be passed directly using this parameter.- boundary_info
A
BoundaryInfospecifying the boundary types of the grid boundary entities. Must be provided ifgridis specified.- num_values
The number of returned vectors of the solution trajectory. If
None, each intermediate vector that is calculated is returned.- time_stepper
The
time-stepperto be used bysolve.- nt
If
time_stepperis not specified, the number of time steps for implicit Euler time stepping.- preassemble
If
True, preassemble all operators in the resultingModel.
Returns
- m
The
Modelthat has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid.- boundary_info
The generated
BoundaryInfo.- unassembled_m
In case
preassembleisTrue, the generatedModelbefore preassembling operators.
- pymor.discretizers.builtin.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.
- class pymor.discretizers.builtin.OnedGrid(domain=(0, 1), num_intervals=4, identify_left_right=False)[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.
- 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)[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=1, **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 1).- kwargs
See
visualize
- class pymor.discretizers.builtin.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
- class pymor.discretizers.builtin.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