pymor.discretizers.builtin package¶
Subpackages¶
Submodules¶
cg module¶
This module provides some operators for continuous finite element discretizations.
-
class
pymor.discretizers.builtin.cg.
AdvectionOperatorP1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear advection
Operator
for linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
v
is vector-valued.Parameters
- grid
The
Grid
for which to assemble the operator.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- advection_function
The
Function
v(x)
withshape_range = (grid.dim, )
. IfNone
, constant one is assumed.- advection_constant
The constant
c
. IfNone
,c
is set to one.- dirichlet_clear_columns
If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.- dirichlet_clear_diag
If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.- solver_options
The
solver_options
for the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.
AdvectionOperatorQ1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear advection
Operator
for bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
v
has to be vector-valued.Parameters
- grid
The
Grid
for which to assemble the operator.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- advection_function
The
Function
v(x)
withshape_range = (grid.dim, )
. IfNone
, constant one is assumed.- advection_constant
The constant
c
. IfNone
,c
is set to one.- dirichlet_clear_columns
If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.- dirichlet_clear_diag
If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.- solver_options
The
solver_options
for the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.
BoundaryDirichletFunctional
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear finite element functional for enforcing Dirichlet boundary values.
Parameters
- grid
Grid
for which to assemble the functional.- dirichlet_data
Function
providing the Dirichlet boundary values.- boundary_info
BoundaryInfo
determining the Dirichlet boundaries.- name
The name of the functional.
-
class
pymor.discretizers.builtin.cg.
BoundaryL2ProductFunctional
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear finite element functional representing the inner product with an L2-
Function
on the boundary.Parameters
- grid
Grid
for which to assemble the functional.- function
The
Function
with which to take the inner product.- boundary_type
The type of domain boundary (e.g. ‘neumann’) on which to assemble the functional. If
None
the functional is assembled over the whole boundary.- dirichlet_clear_dofs
If
True
, set dirichlet boundary DOFs to zero.- boundary_info
If
boundary_type
is specified ordirichlet_clear_dofs
isTrue
, theBoundaryInfo
determining which boundary entity belongs to which physical boundary.- name
The name of the functional.
-
class
pymor.discretizers.builtin.cg.
DiffusionOperatorP1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Diffusion
Operator
for linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
d
can be scalar- or matrix-valued.Parameters
- grid
The
Grid
for which to assemble the operator.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- diffusion_function
The
Function
d(x)
withshape_range == ()
orshape_range = (grid.dim, grid.dim)
. IfNone
, constant one is assumed.- diffusion_constant
The constant
c
. IfNone
,c
is set to one.- dirichlet_clear_columns
If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.- dirichlet_clear_diag
If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.- solver_options
The
solver_options
for the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.
DiffusionOperatorQ1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Diffusion
Operator
for bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
d
can be scalar- or matrix-valued.Parameters
- grid
The
Grid
for which to assemble the operator.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- diffusion_function
The
Function
d(x)
withshape_range == ()
orshape_range = (grid.dim, grid.dim)
. IfNone
, constant one is assumed.- diffusion_constant
The constant
c
. IfNone
,c
is set to one.- dirichlet_clear_columns
If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.- dirichlet_clear_diag
If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.- solver_options
The
solver_options
for the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.
InterpolationOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Vector-like Lagrange interpolation
Operator
for continuous finite element spaces.
-
class
pymor.discretizers.builtin.cg.
L2ProductFunctionalP1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear finite element functional representing the inner product with an L2-
Function
.Parameters
- grid
Grid
for which to assemble the functional.- function
The
Function
with which to take the inner product.- dirichlet_clear_dofs
If
True
, set dirichlet boundary DOFs to zero.- boundary_info
BoundaryInfo
determining the Dirichlet boundaries in casedirichlet_clear_dofs
is set toTrue
.- name
The name of the functional.
-
class
pymor.discretizers.builtin.cg.
L2ProductFunctionalQ1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Bilinear finite element functional representing the inner product with an L2-
Function
.Parameters
- grid
Grid
for which to assemble the functional.- function
The
Function
with which to take the inner product.- dirichlet_clear_dofs
If
True
, set dirichlet boundary DOFs to zero.- boundary_info
BoundaryInfo
determining the Dirichlet boundaries in casedirichlet_clear_dofs
is set toTrue
.- name
The name of the functional.
-
class
pymor.discretizers.builtin.cg.
L2ProductP1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Operator
representing the L2-product between linear finite element functions.Parameters
- grid
The
Grid
for which to assemble the product.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- dirichlet_clear_rows
If
True
, set the rows of the system matrix corresponding to Dirichlet boundary DOFs to zero.- dirichlet_clear_columns
If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero.- dirichlet_clear_diag
If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise, if eitherdirichlet_clear_rows
ordirichlet_clear_columns
isTrue
, the diagonal entries are set to one.- coefficient_function
Coefficient
Function
for product withshape_range == ()
. IfNone
, constant one is assumed.- solver_options
The
solver_options
for the operator.- name
The name of the product.
-
class
pymor.discretizers.builtin.cg.
L2ProductQ1
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Operator
representing the L2-product between bilinear finite element functions.Parameters
- grid
The
Grid
for which to assemble the product.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- dirichlet_clear_rows
If
True
, set the rows of the system matrix corresponding to Dirichlet boundary DOFs to zero.- dirichlet_clear_columns
If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero.- dirichlet_clear_diag
If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise, if eitherdirichlet_clear_rows
ordirichlet_clear_columns
isTrue
, the diagonal entries are set to one.- coefficient_function
Coefficient
Function
for product withshape_range == ()
. IfNone
, constant one is assumed.- solver_options
The
solver_options
for the operator.- name
The name of the product.
-
class
pymor.discretizers.builtin.cg.
RobinBoundaryOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Robin boundary
Operator
for linear finite elements.The operator represents the contribution of Robin boundary conditions to the stiffness matrix, where the boundary condition is supposed to be given in the form
-[ d(x) ∇u(x) ] ⋅ n(x) = c(x) (u(x) - g(x))
d
andn
are the diffusion function (seeDiffusionOperatorP1
) and the unit outer normal inx
, whilec
is the (scalar) Robin parameter function andg
is the (also scalar) Robin boundary value function.Parameters
- grid
The
Grid
over which to assemble the operator.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- robin_data
Tuple providing two
Functions
that represent the Robin parameter and boundary value function. IfNone
, the resulting operator is zero.- solver_options
The
solver_options
for the operator.- name
Name of the operator.
-
pymor.discretizers.builtin.cg.
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
InstationaryProblem
with aStationaryProblem
as stationary part using finite elements.Parameters
- analytical_problem
The
InstationaryProblem
to discretize.- diameter
If not
None
,diameter
is 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_default
is used.- grid_type
If not
None
, this parameter is forwarded todomain_discretizer
to specify the type of the generatedGrid
.- grid
Instead of using a domain discretizer, the
Grid
can also be passed directly using this parameter.- boundary_info
A
BoundaryInfo
specifying the boundary types of the grid boundary entities. Must be provided ifgrid
is 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-stepper
to be used bysolve
.- nt
If
time_stepper
is not specified, the number of time steps for implicit Euler time stepping.- preassemble
If
True
, preassemble all operators in the resultingModel
.
Returns
- m
The
Model
that has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid
.- boundary_info
The generated
BoundaryInfo
.- unassembled_m
In case
preassemble
isTrue
, the generatedModel
before preassembling operators.
-
pymor.discretizers.builtin.cg.
discretize_stationary_cg
(analytical_problem, diameter=None, domain_discretizer=None, grid_type=None, grid=None, boundary_info=None, preassemble=True)[source]¶ Discretizes a
StationaryProblem
using finite elements.Parameters
- analytical_problem
The
StationaryProblem
to discretize.- diameter
If not
None
,diameter
is 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_default
is used.- grid_type
If not
None
, this parameter is forwarded todomain_discretizer
to specify the type of the generatedGrid
.- grid
Instead of using a domain discretizer, the
Grid
can also be passed directly using this parameter.- boundary_info
A
BoundaryInfo
specifying the boundary types of the grid boundary entities. Must be provided ifgrid
is specified.- preassemble
If
True
, preassemble all operators in the resultingModel
.
Returns
- m
The
Model
that has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid
.- boundary_info
The generated
BoundaryInfo
.- unassembled_m
In case
preassemble
isTrue
, the generatedModel
before preassembling operators.
fv module¶
This module provides some operators for finite volume discretizations.
-
class
pymor.discretizers.builtin.fv.
DiffusionOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Finite Volume Diffusion
Operator
.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
Parameters
- grid
The
Grid
over which to assemble the operator.- boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- diffusion_function
The scalar-valued
Function
d(x)
. IfNone
, constant one is assumed.- diffusion_constant
The constant
c
. IfNone
,c
is set to one.- solver_options
The
solver_options
for the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.fv.
EngquistOsherFlux
(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.fv.NumericalConvectiveFlux
Engquist-Osher numerical flux.
If
f
is the analytical flux, andf'
its derivative, the Engquist-Osher flux is given by:F(U_in, U_out, normal, vol) = vol * [c^+(U_in, normal) + c^-(U_out, normal)] U_in c^+(U_in, normal) = f(0)⋅normal + ∫ max(f'(s)⋅normal, 0) ds s=0 U_out c^-(U_out, normal) = ∫ min(f'(s)⋅normal, 0) ds s=0
Parameters
Methods
evaluate_stage1
,evaluate_stage2
-
class
pymor.discretizers.builtin.fv.
L2Product
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Operator
representing the L2-product between finite volume functions.Parameters
- grid
The
Grid
for which to assemble the product.- solver_options
The
solver_options
for the operator.- name
The name of the product.
-
class
pymor.discretizers.builtin.fv.
L2ProductFunctional
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Finite volume functional representing the inner product with an L2-
Function
.Additionally, boundary conditions can be enforced by providing
dirichlet_data
andneumann_data
functions.Parameters
- grid
Grid
for which to assemble the functional.- function
The
Function
with which to take the inner product orNone
.- boundary_info
BoundaryInfo
determining the Dirichlet and Neumann boundaries orNone
. IfNone
, no boundary treatment is performed.- dirichlet_data
Function
providing the Dirichlet boundary values. IfNone
, constant-zero boundary is assumed.- diffusion_function
See
DiffusionOperator
. Has to be specified in casedirichlet_data
is given.- diffusion_constant
See
DiffusionOperator
. Has to be specified in casedirichlet_data
is given.- neumann_data
Function
providing the Neumann boundary values. IfNone
, constant-zero is assumed.- order
Order of the Gauss quadrature to use for numerical integration.
- name
The name of the functional.
-
class
pymor.discretizers.builtin.fv.
LaxFriedrichsFlux
(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.fv.NumericalConvectiveFlux
Lax-Friedrichs numerical flux.
If
f
is the analytical flux, the Lax-Friedrichs fluxF
is given by:F(U_in, U_out, normal, vol) = vol * [normal⋅(f(U_in) + f(U_out))/2 + (U_in - U_out)/(2*λ)]
Parameters
- flux
Function
defining the analytical fluxf
.- lxf_lambda
The stabilization parameter
λ
.
Methods
evaluate_stage1
,evaluate_stage2
-
class
pymor.discretizers.builtin.fv.
LinearAdvectionLaxFriedrichs
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear advection finite Volume
Operator
using Lax-Friedrichs flux.The operator is of the form
L(u, mu)(x) = ∇ ⋅ (v(x, mu)⋅u(x))
See
LaxFriedrichsFlux
for the definition of the Lax-Friedrichs flux.Parameters
- grid
Grid
over which to assemble the operator.- boundary_info
BoundaryInfo
determining the Dirichlet and Neumann boundaries.- velocity_field
Function
defining the velocity fieldv
.- lxf_lambda
The stabilization parameter
λ
.- solver_options
The
solver_options
for the operator.- name
The name of the operator.
-
class
pymor.discretizers.builtin.fv.
NonlinearAdvectionOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.interface.Operator
Nonlinear finite volume advection
Operator
.The operator is of the form
L(u, mu)(x) = ∇ ⋅ f(u(x), mu)
Parameters
- grid
Grid
for which to evaluate the operator.- boundary_info
BoundaryInfo
determining the Dirichlet and Neumann boundaries.- numerical_flux
The
NumericalConvectiveFlux
to use.- dirichlet_data
Function
providing the Dirichlet boundary values. IfNone
, constant-zero boundary is assumed.- solver_options
The
solver_options
for the operator.- name
The name of the operator.
Methods
apply
,jacobian
,restricted
,with_numerical_flux
apply2
,apply_adjoint
,apply_inverse
,apply_inverse_adjoint
,as_range_array
,as_source_array
,as_vector
,assemble
,d_mu
,pairwise_apply2
,__add__
,__matmul__
,__mul__
,__radd__
Attributes
parameters
,parameters_inherited
,parameters_internal
,parameters_own
,parametric
-
apply
(U, mu=None)[source]¶ Apply the operator to a
VectorArray
.Parameters
- U
VectorArray
of vectors to which the operator is applied.- mu
The
parameter values
for which to evaluate the operator.
Returns
VectorArray
of the operator evaluations.
-
jacobian
(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator
.Parameters
- U
Length 1
VectorArray
containing the vector for which to compute the Jacobian.- mu
The
parameter values
for which to compute the Jacobian.
Returns
Linear
Operator
representing the Jacobian.
-
restricted
(dofs)[source]¶ Restrict the operator range to a given set of degrees of freedom.
This method returns a restricted version
restricted_op
of the operator along with an arraysource_dofs
such that for anyVectorArray
U
inself.source
the following is true:self.apply(U, mu).dofs(dofs) == restricted_op.apply(NumpyVectorArray(U.dofs(source_dofs)), mu))
Such an operator is mainly useful for
empirical interpolation
where the evaluation of the original operator only needs to be known for few selected degrees of freedom. If the operator has a small stencil, only fewsource_dofs
will be needed to evaluate the restricted operator which can make its evaluation very fast compared to evaluating the original operator.Parameters
- dofs
One-dimensional
NumPy array
of degrees of freedom in the operatorrange
to which to restrict.
Returns
- restricted_op
The restricted operator as defined above. The operator will have
NumpyVectorSpace
(len(source_dofs))
assource
andNumpyVectorSpace
(len(dofs))
asrange
.- source_dofs
One-dimensional
NumPy array
of source degrees of freedom as defined above.
-
class
pymor.discretizers.builtin.fv.
NonlinearReactionOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.interface.Operator
Methods
apply2
,apply_adjoint
,apply_inverse
,apply_inverse_adjoint
,as_range_array
,as_source_array
,as_vector
,assemble
,d_mu
,pairwise_apply2
,restricted
,__add__
,__matmul__
,__mul__
,__radd__
Attributes
parameters
,parameters_inherited
,parameters_internal
,parameters_own
,parametric
-
apply
(U, ind=None, mu=None)[source]¶ Apply the operator to a
VectorArray
.Parameters
- U
VectorArray
of vectors to which the operator is applied.- mu
The
parameter values
for which to evaluate the operator.
Returns
VectorArray
of the operator evaluations.
-
jacobian
(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator
.Parameters
- U
Length 1
VectorArray
containing the vector for which to compute the Jacobian.- mu
The
parameter values
for which to compute the Jacobian.
Returns
Linear
Operator
representing the Jacobian.
-
-
class
pymor.discretizers.builtin.fv.
NumericalConvectiveFlux
(*args, **kwargs)[source]¶ Bases:
pymor.parameters.base.ParametricObject
Interface for numerical convective fluxes for finite volume schemes.
Numerical fluxes defined by this interfaces are functions of the form
F(U_inner, U_outer, unit_outer_normal, edge_volume, mu)
.- The flux evaluation is vectorized and happens in two stages:
evaluate_stage1
receives aNumPy array
U
of all values which appear asU_inner
orU_outer
for all edges the flux shall be evaluated at and returns atuple
ofNumPy arrays
each of the same length asU
.evaluate_stage2
receives the reorderedstage1_data
for each edge as well as the unit outer normal and the volume of the edges.stage1_data
is given as follows: IfR_l
isl
-th entry of thetuple
returned byevaluate_stage1
, thel
-th entryD_l
of of thestage1_data
tuple has the shape(num_edges, 2) + R_l.shape[1:]
. If for edgek
the valuesU_inner
andU_outer
are thei
-th andj
-th value in theU
array provided toevaluate_stage1
, we haveD_l[k, 0] == R_l[i], D_l[k, 1] == R_l[j].
evaluate_stage2
returns aNumPy array
of the flux evaluations for each edge.
Methods
evaluate_stage1
,evaluate_stage2
-
class
pymor.discretizers.builtin.fv.
ReactionOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Finite Volume reaction
Operator
.The operator is of the form
L(u, mu)(x) = c(x, mu)⋅u(x)
Parameters
- grid
The
Grid
for which to assemble the operator.- reaction_coefficient
The function ‘c’
- solver_options
The
solver_options
for the operator.- name
The name of the operator.
-
class
pymor.discretizers.builtin.fv.
SimplifiedEngquistOsherFlux
(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.fv.NumericalConvectiveFlux
Engquist-Osher numerical flux. Simplified Implementation for special case.
For the definition of the Engquist-Osher flux see
EngquistOsherFlux
. This class provides a faster and more accurate implementation for the special case thatf(0) == 0
and the derivative off
only changes sign at0
.Parameters
Methods
evaluate_stage1
,evaluate_stage2
-
pymor.discretizers.builtin.fv.
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]¶ Discretizes an
InstationaryProblem
with aStationaryProblem
as stationary part using the finite volume method.Parameters
- analytical_problem
The
InstationaryProblem
to discretize.- diameter
If not
None
,diameter
is 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_default
is used.- grid_type
If not
None
, this parameter is forwarded todomain_discretizer
to 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
Grid
can also be passed directly using this parameter.- boundary_info
A
BoundaryInfo
specifying the boundary types of the grid boundary entities. Must be provided ifgrid
is 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-stepper
to be used bysolve
.- nt
If
time_stepper
is not specified, the number of time steps for implicit Euler time stepping.- preassemble
If
True
, preassemble all operators in the resultingModel
.
Returns
- m
The
Model
that has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid
.- boundary_info
The generated
BoundaryInfo
.- unassembled_m
In case
preassemble
isTrue
, the generatedModel
before preassembling operators.
-
pymor.discretizers.builtin.fv.
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
StationaryProblem
using the finite volume method.Parameters
- analytical_problem
The
StationaryProblem
to discretize.- diameter
If not
None
,diameter
is 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_default
is used.- grid_type
If not
None
, this parameter is forwarded todomain_discretizer
to 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
Grid
can also be passed directly using this parameter.- boundary_info
A
BoundaryInfo
specifying the boundary types of the grid boundary entities. Must be provided ifgrid
is specified.- preassemble
If
True
, preassemble all operators in the resultingModel
.
Returns
- m
The
Model
that has been generated.- data
Dictionary with the following entries:
- grid
The generated
Grid
.- boundary_info
The generated
BoundaryInfo
.- unassembled_m
In case
preassemble
isTrue
, the generatedModel
before preassembling operators.
-
pymor.discretizers.builtin.fv.
nonlinear_advection_engquist_osher_operator
(grid, boundary_info, flux, flux_derivative, gausspoints=5, intervals=1, dirichlet_data=None, solver_options=None, name=None)[source]¶ Instantiate a
NonlinearAdvectionOperator
usingEngquistOsherFlux
.
-
pymor.discretizers.builtin.fv.
nonlinear_advection_lax_friedrichs_operator
(grid, boundary_info, flux, lxf_lambda=1.0, dirichlet_data=None, solver_options=None, name=None)[source]¶ Instantiate a
NonlinearAdvectionOperator
usingLaxFriedrichsFlux
.
-
pymor.discretizers.builtin.fv.
nonlinear_advection_simplified_engquist_osher_operator
(grid, boundary_info, flux, flux_derivative, dirichlet_data=None, solver_options=None, name=None)[source]¶ Instantiate a
NonlinearAdvectionOperator
usingSimplifiedEngquistOsherFlux
.
inplace module¶
inverse module¶
-
pymor.discretizers.builtin.inverse.
inv_transposed_two_by_two
(A)[source]¶ Efficiently compute the tranposed inverses of a
NumPy array
of 2x2-matrices| retval[i1,...,ik,m,n] = numpy.linalg.inv(A[i1,...,ik,:,:]).
-
pymor.discretizers.builtin.inverse.
inv_two_by_two
(A)[source]¶ Efficiently compute the inverses of a
NumPy array
of 2x2-matrices| retval[i1,...,ik,m,n] = numpy.linalg.inv(A[i1,...,ik,:,:]).
list module¶
-
class
pymor.discretizers.builtin.list.
ConvertToNumpyListVectorArrayRules
[source]¶ Bases:
pymor.algorithms.rules.RuleTable
Methods
append_rule
,apply
,apply_children
,breakpoint_for_name
,breakpoint_for_obj
,get_children
,insert_rule
,replace_children
Attributes
action_NumpyMatrixOperator
,action_recurse
,action_VectorArrayOperator
,rules
-
pymor.discretizers.builtin.list.
convert_to_numpy_list_vector_array
(obj)[source]¶ Use NumpyListVectorArrayMatrixOperator instead of NumpyMatrixOperator.
This simple function recursively converts
NumpyMatrixOperators
to correspondingNumpyListVectorArrayMatrixOperators
.Parameters
- obj
Either an
Operator
, e.g.NumpyMatrixOperator
or aLincombOperator
ofNumpyMatrixOperators
, or an entireModel
that is to be converted.
quadratures module¶
-
class
pymor.discretizers.builtin.quadratures.
GaussQuadratures
[source]¶ Bases:
object
Gauss quadrature on the interval [0, 1]
Methods
iter_quadrature
,maxpoints
,quadrature
Attributes
a
,order_map
,orders
,points
,weights