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.NumpyMatrixBasedOperatorLinear advection
Operatorfor linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
vhas to be vector-valued.Parameters
- grid
The
Gridfor which to assemble the operator.- boundary_info
BoundaryInfofor the treatment of Dirichlet boundary conditions.- advection_function
The
Functionv(x)withshape_range = (grid.dim, ). IfNone, constant one is assumed.- advection_constant
The constant
c. IfNone,cis 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_optionsfor the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.AdvectionOperatorQ1(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear advection
Operatorfor bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
vhas to be vector-valued.Parameters
- grid
The
Gridfor which to assemble the operator.- boundary_info
BoundaryInfofor the treatment of Dirichlet boundary conditions.- advection_function
The
Functionv(x)withshape_range = (grid.dim, ). IfNone, constant one is assumed.- advection_constant
The constant
c. IfNone,cis 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_optionsfor the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.BoundaryDirichletFunctional(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear finite element functional for enforcing Dirichlet boundary values.
Parameters
- grid
Gridfor which to assemble the functional.- dirichlet_data
Functionproviding the Dirichlet boundary values.- boundary_info
BoundaryInfodetermining the Dirichlet boundaries.- name
The name of the functional.
-
class
pymor.discretizers.builtin.cg.BoundaryL2ProductFunctional(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear finite element functional representing the inner product with an L2-
Functionon the boundary.Parameters
- grid
Gridfor which to assemble the functional.- function
The
Functionwith which to take the inner product.- boundary_type
The type of domain boundary (e.g. ‘neumann’) on which to assemble the functional. If
Nonethe functional is assembled over the whole boundary.- dirichlet_clear_dofs
If
True, set dirichlet boundary DOFs to zero.- boundary_info
If
boundary_typeis specified ordirichlet_clear_dofsisTrue, theBoundaryInfodetermining 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.NumpyMatrixBasedOperatorDiffusion
Operatorfor linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
dcan be scalar- or matrix-valued.Parameters
- grid
The
Gridfor which to assemble the operator.- boundary_info
BoundaryInfofor the treatment of Dirichlet boundary conditions.- diffusion_function
The
Functiond(x)withshape_range == ()orshape_range = (grid.dim, grid.dim). IfNone, constant one is assumed.- diffusion_constant
The constant
c. IfNone,cis 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_optionsfor the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.DiffusionOperatorQ1(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorDiffusion
Operatorfor bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
dcan be scalar- or matrix-valued.Parameters
- grid
The
Gridfor which to assemble the operator.- boundary_info
BoundaryInfofor the treatment of Dirichlet boundary conditions.- diffusion_function
The
Functiond(x)withshape_range == ()orshape_range = (grid.dim, grid.dim). IfNone, constant one is assumed.- diffusion_constant
The constant
c. IfNone,cis 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_optionsfor the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.cg.InterpolationOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorVector-like Lagrange interpolation
Operatorfor continuous finite element spaces.
-
class
pymor.discretizers.builtin.cg.L2ProductFunctionalP1(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear finite element functional representing the inner product with an L2-
Function.Parameters
- grid
Gridfor which to assemble the functional.- function
The
Functionwith which to take the inner product.- dirichlet_clear_dofs
If
True, set dirichlet boundary DOFs to zero.- boundary_info
BoundaryInfodetermining the Dirichlet boundaries in casedirichlet_clear_dofsis set toTrue.- name
The name of the functional.
-
class
pymor.discretizers.builtin.cg.L2ProductFunctionalQ1(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorBilinear finite element functional representing the inner product with an L2-
Function.Parameters
- grid
Gridfor which to assemble the functional.- function
The
Functionwith which to take the inner product.- dirichlet_clear_dofs
If
True, set dirichlet boundary DOFs to zero.- boundary_info
BoundaryInfodetermining the Dirichlet boundaries in casedirichlet_clear_dofsis set toTrue.- name
The name of the functional.
-
class
pymor.discretizers.builtin.cg.L2ProductP1(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorOperatorrepresenting the L2-product between linear finite element functions.Parameters
- grid
The
Gridfor which to assemble the product.- boundary_info
BoundaryInfofor 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_rowsordirichlet_clear_columnsisTrue, the diagonal entries are set to one.- coefficient_function
Coefficient
Functionfor product withshape_range == (). IfNone, constant one is assumed.- solver_options
The
solver_optionsfor the operator.- name
The name of the product.
-
class
pymor.discretizers.builtin.cg.L2ProductQ1(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorOperatorrepresenting the L2-product between bilinear finite element functions.Parameters
- grid
The
Gridfor which to assemble the product.- boundary_info
BoundaryInfofor 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_rowsordirichlet_clear_columnsisTrue, the diagonal entries are set to one.- coefficient_function
Coefficient
Functionfor product withshape_range == (). IfNone, constant one is assumed.- solver_options
The
solver_optionsfor the operator.- name
The name of the product.
-
class
pymor.discretizers.builtin.cg.RobinBoundaryOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorRobin boundary
Operatorfor 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))
dandnare the diffusion function (seeDiffusionOperatorP1) and the unit outer normal inx, whilecis the (scalar) Robin parameter function andgis the (also scalar) Robin boundary value function.Parameters
- grid
The
Gridover which to assemble the operator.- boundary_info
BoundaryInfofor the treatment of Dirichlet boundary conditions.- robin_data
Tuple providing two
Functionsthat represent the Robin parameter and boundary value function. IfNone, the resulting operator is zero.- solver_options
The
solver_optionsfor 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
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.cg.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.
fv module¶
This module provides some operators for finite volume discretizations.
-
class
pymor.discretizers.builtin.fv.DiffusionOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorFinite Volume Diffusion
Operator.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
Parameters
- grid
The
Gridover which to assemble the operator.- boundary_info
BoundaryInfofor the treatment of Dirichlet boundary conditions.- diffusion_function
The scalar-valued
Functiond(x). IfNone, constant one is assumed.- diffusion_constant
The constant
c. IfNone,cis set to one.- solver_options
The
solver_optionsfor the operator.- name
Name of the operator.
-
class
pymor.discretizers.builtin.fv.EngquistOsherFlux(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.fv.NumericalConvectiveFluxEngquist-Osher numerical flux.
If
fis 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=0Parameters
Methods
evaluate_stage1,evaluate_stage2
-
class
pymor.discretizers.builtin.fv.L2Product(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorOperatorrepresenting the L2-product between finite volume functions.Parameters
- grid
The
Gridfor which to assemble the product.- solver_options
The
solver_optionsfor the operator.- name
The name of the product.
-
class
pymor.discretizers.builtin.fv.L2ProductFunctional(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorFinite volume functional representing the inner product with an L2-
Function.Additionally, boundary conditions can be enforced by providing
dirichlet_dataandneumann_datafunctions.Parameters
- grid
Gridfor which to assemble the functional.- function
The
Functionwith which to take the inner product orNone.- boundary_info
BoundaryInfodetermining the Dirichlet and Neumann boundaries orNone. IfNone, no boundary treatment is performed.- dirichlet_data
Functionproviding the Dirichlet boundary values. IfNone, constant-zero boundary is assumed.- diffusion_function
See
DiffusionOperator. Has to be specified in casedirichlet_datais given.- diffusion_constant
See
DiffusionOperator. Has to be specified in casedirichlet_datais given.- neumann_data
Functionproviding 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.NumericalConvectiveFluxLax-Friedrichs numerical flux.
If
fis the analytical flux, the Lax-Friedrichs fluxFis 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
Functiondefining the analytical fluxf.- lxf_lambda
The stabilization parameter
λ.
Methods
evaluate_stage1,evaluate_stage2
-
class
pymor.discretizers.builtin.fv.LinearAdvectionLaxFriedrichsOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear advection finite Volume
Operatorusing Lax-Friedrichs flux.The operator is of the form
L(u, mu)(x) = ∇ ⋅ (v(x, mu)⋅u(x))
See
LaxFriedrichsFluxfor the definition of the Lax-Friedrichs flux.Parameters
- grid
Gridover which to assemble the operator.- boundary_info
BoundaryInfodetermining the Dirichlet and Neumann boundaries.- velocity_field
Functiondefining the velocity fieldv.- lxf_lambda
The stabilization parameter
λ.- solver_options
The
solver_optionsfor the operator.- name
The name of the operator.
-
class
pymor.discretizers.builtin.fv.NonlinearAdvectionOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.interface.OperatorNonlinear finite volume advection
Operator.The operator is of the form
L(u, mu)(x) = ∇ ⋅ f(u(x), mu)
Parameters
- grid
Gridfor which to evaluate the operator.- boundary_info
BoundaryInfodetermining the Dirichlet and Neumann boundaries.- numerical_flux
The
NumericalConvectiveFluxto use.- dirichlet_data
Functionproviding the Dirichlet boundary values. IfNone, constant-zero boundary is assumed.- solver_options
The
solver_optionsfor the operator.- name
The name of the operator.
Methods
apply,jacobian,restricted,with_numerical_fluxapply2,apply_adjoint,apply_inverse,apply_inverse_adjoint,as_range_array,as_source_array,as_vector,assemble,d_mu,pairwise_apply2,__matmul__Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
apply(U, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
parameter valuesfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
jacobian(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator.Parameters
- U
Length 1
VectorArraycontaining the vector for which to compute the Jacobian.- mu
The
parameter valuesfor which to compute the Jacobian.
Returns
Linear
Operatorrepresenting the Jacobian.
-
restricted(dofs)[source]¶ Restrict the operator range to a given set of degrees of freedom.
This method returns a restricted version
restricted_opof the operator along with an arraysource_dofssuch that for anyVectorArrayUinself.sourcethe 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 interpolationwhere 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_dofswill 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 arrayof degrees of freedom in the operatorrangeto which to restrict.
Returns
- restricted_op
The restricted operator as defined above. The operator will have
NumpyVectorSpace(len(source_dofs))assourceandNumpyVectorSpace(len(dofs))asrange.- source_dofs
One-dimensional
NumPy arrayof source degrees of freedom as defined above.
-
class
pymor.discretizers.builtin.fv.NonlinearReactionOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.interface.OperatorMethods
apply2,apply_adjoint,apply_inverse,apply_inverse_adjoint,as_range_array,as_source_array,as_vector,assemble,d_mu,pairwise_apply2,restricted,__matmul__Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
apply(U, ind=None, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
parameter valuesfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
jacobian(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator.Parameters
- U
Length 1
VectorArraycontaining the vector for which to compute the Jacobian.- mu
The
parameter valuesfor which to compute the Jacobian.
Returns
Linear
Operatorrepresenting the Jacobian.
-
-
class
pymor.discretizers.builtin.fv.NumericalConvectiveFlux(*args, **kwargs)[source]¶ Bases:
pymor.parameters.base.ParametricObjectInterface 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_stage1receives aNumPy arrayUof all values which appear asU_innerorU_outerfor all edges the flux shall be evaluated at and returns atupleofNumPy arrayseach of the same length asU.evaluate_stage2receives the reorderedstage1_datafor each edge as well as the unit outer normal and the volume of the edges.stage1_datais given as follows: IfR_lisl-th entry of thetuplereturned byevaluate_stage1, thel-th entryD_lof of thestage1_datatuple has the shape(num_edges, 2) + R_l.shape[1:]. If for edgekthe valuesU_innerandU_outerare thei-th andj-th value in theUarray provided toevaluate_stage1, we haveD_l[k, 0] == R_l[i], D_l[k, 1] == R_l[j].
evaluate_stage2returns aNumPy arrayof the flux evaluations for each edge.
Methods
evaluate_stage1,evaluate_stage2
-
class
pymor.discretizers.builtin.fv.ReactionOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorFinite Volume reaction
Operator.The operator is of the form
L(u, mu)(x) = c(x, mu)⋅u(x)
Parameters
- grid
The
Gridfor which to assemble the operator.- reaction_coefficient
The function ‘c’
- solver_options
The
solver_optionsfor the operator.- name
The name of the operator.
-
class
pymor.discretizers.builtin.fv.SimplifiedEngquistOsherFlux(*args, **kwargs)[source]¶ Bases:
pymor.discretizers.builtin.fv.NumericalConvectiveFluxEngquist-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) == 0and the derivative offonly 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
InstationaryProblemwith aStationaryProblemas stationary part using the finite volume method.Parameters
- 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.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
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.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
NonlinearAdvectionOperatorusingEngquistOsherFlux.
-
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
NonlinearAdvectionOperatorusingLaxFriedrichsFlux.
-
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
NonlinearAdvectionOperatorusingSimplifiedEngquistOsherFlux.
inplace module¶
inverse module¶
-
pymor.discretizers.builtin.inverse.inv_transposed_two_by_two(A)[source]¶ Efficiently compute the tranposed inverses of a
NumPy arrayof 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 arrayof 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.RuleTableMethods
append_rule,apply,apply_children,breakpoint_for_name,breakpoint_for_obj,get_children,insert_rule,replace_childrenAttributes
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
NumpyMatrixOperatorsto correspondingNumpyListVectorArrayMatrixOperators.Parameters
- obj
Either an
Operator, e.g.NumpyMatrixOperatoror aLincombOperatorofNumpyMatrixOperators, or an entireModelthat is to be converted.
quadratures module¶
-
class
pymor.discretizers.builtin.quadratures.GaussQuadratures[source]¶ Bases:
objectGauss quadrature on the interval [0, 1]
Methods
iter_quadrature,maxpoints,quadratureAttributes
a,order_map,orders,points,weights