pymor.discretizers.builtin.fv¶
This module provides some operators for finite volume discretizations.
Module Contents¶
- class pymor.discretizers.builtin.fv.BoundaryL2Functional(grid, function, boundary_type=None, boundary_info=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorFV 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.boundary_info – If
boundary_typeis specified, theBoundaryInfodetermining which boundary entity belongs to which physical boundary.name – The name of the functional.
- class pymor.discretizers.builtin.fv.DiffusionOperator(grid, boundary_info, diffusion_function=None, diffusion_constant=None, solver=None, name=None)[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 – The
Solverfor the operator.name – Name of the operator.
- class pymor.discretizers.builtin.fv.EngquistOsherFlux(flux, flux_derivative, gausspoints=5, intervals=1)[source]¶
Bases:
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=0- Parameters:
Methods
- class pymor.discretizers.builtin.fv.InterpolationOperator(grid, function, order=0)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorVector-like L^2-projection interpolation
Operatorfor finite volume spaces.- Parameters:
- class pymor.discretizers.builtin.fv.L2Functional(grid, function=None, boundary_info=None, dirichlet_data=None, diffusion_function=None, diffusion_constant=None, neumann_data=None, order=1, name=None)[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.L2Product(grid, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorOperatorrepresenting the L2-product between finite volume functions.- Parameters:
- class pymor.discretizers.builtin.fv.LaxFriedrichsFlux(flux, lxf_lambda=1.0)[source]¶
Bases:
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
- class pymor.discretizers.builtin.fv.LinearAdvectionLaxFriedrichsOperator(grid, boundary_info, velocity_field, lxf_lambda=1.0, solver=None, name=None)[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 – The
Solverfor the operator.name – The name of the operator.
- class pymor.discretizers.builtin.fv.NonlinearAdvectionOperator(grid, boundary_info, numerical_flux, dirichlet_data=None, jacobian_delta=1e-07, solver=None, name=None)[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 – The
Solverfor the operator.name – The name of the operator.
Methods
Apply the operator to a
VectorArray.Return the operator's Jacobian as a new
Operator.Restrict the operator range to a given set of degrees of freedom.
- 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:
|VectorArray| of the operator evaluations.
- jacobian(U, mu=None)[source]¶
Return the operator’s Jacobian as a new
Operator.If the operator has a
Solver, the JacobianOperatorwill be equipped with the solver’sjacobian_solver.- Parameters:
U – Length 1
VectorArraycontaining the vector for which to compute the Jacobian.mu – The
parameter valuesfor 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_opof the operator along with an arraysource_dofssuch that for anyVectorArrayUinself.sourcethe following is true:self.apply(U, mu).dofs(dofs) == restricted_op.apply(restricted_op.source.from_numpy(U.dofs(source_dofs)) mu).to_numpy()
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(grid, reaction_function, reaction_function_derivative=None, name=None)[source]¶
Bases:
pymor.operators.interface.OperatorFinite Volume nonlinear reaction
Operator.- Parameters:
grid – The
Gridfor which to assemble the operator.reaction_function – The reaction function.
reaction_function_derivative – The reaction function derivative.
name – The name of the operator.
Methods
Apply the operator to a
VectorArray.Return the operator's Jacobian as a new
Operator.- 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:
|VectorArray| of the operator evaluations.
- jacobian(U, mu=None)[source]¶
Return the operator’s Jacobian as a new
Operator.If the operator has a
Solver, the JacobianOperatorwill be equipped with the solver’sjacobian_solver.- Parameters:
U – Length 1
VectorArraycontaining the vector for which to compute the Jacobian.mu – The
parameter valuesfor which to compute the Jacobian.
- Returns:
Linear |Operator| representing the Jacobian.
- class pymor.discretizers.builtin.fv.NumericalConvectiveFlux[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
- class pymor.discretizers.builtin.fv.ReactionOperator(grid, reaction_coefficient, solver=None, name=None)[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:
- class pymor.discretizers.builtin.fv.SimplifiedEngquistOsherFlux(flux, flux_derivative)[source]¶
Bases:
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
- 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, solver=None)[source]¶
FV Discretization of an
InstationaryProblemwith aStationaryProblemas stationary part.- 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 bysolveofInstationaryModel.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.solver – The
Solverto be used.
- 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, solver=None)[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.solver – The
Solverto be used.
- 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, jacobian_delta=None, solver=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, jacobian_delta=None, solver=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, jacobian_delta=None, solver=None, name=None)[source]¶
Create a
NonlinearAdvectionOperatorusingSimplifiedEngquistOsherFlux.