pymor.analyticalproblems package¶
Submodules¶
burgers module¶
-
pymor.analyticalproblems.burgers.burgers_problem(v=1.0, circle=True, initial_data_type='sin', parameter_range=1.0, 2.0)[source]¶ One-dimensional Burgers-type problem.
The problem is to solve
∂_t u(x, t, μ) + ∂_x (v * u(x, t, μ)^μ) = 0 u(x, 0, μ) = u_0(x)for u with t in [0, 0.3] and x in [0, 2].
Parameters
- v
The velocity v.
- circle
If
True, impose periodic boundary conditions. Otherwise Dirichlet left, outflow right.- initial_data_type
Type of initial data (
'sin'or'bump').- parameter_range
The interval in which μ is allowed to vary.
-
pymor.analyticalproblems.burgers.burgers_problem_2d(vx=1.0, vy=1.0, torus=True, initial_data_type='sin', parameter_range=1.0, 2.0)[source]¶ Two-dimensional Burgers-type problem.
The problem is to solve
∂_t u(x, t, μ) + ∇ ⋅ (v * u(x, t, μ)^μ) = 0 u(x, 0, μ) = u_0(x)for u with t in [0, 0.3], x in [0, 2] x [0, 1].
Parameters
- vx
The x component of the velocity vector v.
- vy
The y component of the velocity vector v.
- torus
If
True, impose periodic boundary conditions. Otherwise, Dirichlet left and bottom, outflow top and right.- initial_data_type
Type of initial data (
'sin'or'bump').- parameter_range
The interval in which μ is allowed to vary.
domaindescriptions module¶
-
class
pymor.analyticalproblems.domaindescriptions.CircleDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.DomainDescriptionDescribes a domain with the topology of a circle, i.e. a line with identified end points.
Parameters
- domain
List [x_l, x_r] providing the left and right endpoint.
Attributes
boundary_types,has_dirichlet,has_neumann,has_robin-
domain¶
-
class
pymor.analyticalproblems.domaindescriptions.CircularSectorDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.PolygonalDomainDescribes a circular sector domain of variable radius.
Parameters
- angle
The angle between 0 and 2*pi of the circular sector.
- radius
The radius of the circular sector.
- arc
The boundary type of the arc.
- radii
The boundary type of the two radii.
- num_points
The number of points of the polygonal chain approximating the circular boundary.
Attributes
angle,arc,num_points,radii,radiusboundary_types,has_dirichlet,has_neumann,has_robin-
angle¶
-
radius¶
-
arc¶
-
radii¶
-
num_points¶
-
class
pymor.analyticalproblems.domaindescriptions.CylindricalDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.DomainDescriptionDescribes a cylindrical domain.
Boundary types can be associated edgewise.
Parameters
- domain
List of two points defining the lower-left and upper-right corner of the domain. The left and right edge are identified.
- top
The boundary type of the top edge.
- bottom
The boundary type of the bottom edge.
Attributes
bottom,diameter,dim,domain,height,lower_left,top,upper_right,volume,widthboundary_types,has_dirichlet,has_neumann,has_robin-
domain¶
-
top¶
-
bottom¶
-
class
pymor.analyticalproblems.domaindescriptions.DiscDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.PolygonalDomainDescribes a disc domain of variable radius.
Parameters
- radius
The radius of the disc.
- boundary
The boundary type of the boundary.
- num_points
The number of points of the polygonal chain approximating the boundary.
Attributes
boundary_types,has_dirichlet,has_neumann,has_robin-
radius¶
-
boundary¶
-
num_points¶
-
class
pymor.analyticalproblems.domaindescriptions.DomainDescription(*args, **kwargs)[source]¶ Bases:
pymor.core.base.ImmutableObjectDescribes a geometric domain along with its boundary.
Attributes
boundary_types,dim,has_dirichlet,has_neumann,has_robin-
dim¶ The dimension of the domain
-
boundary_types¶ Set of boundary types the domain has.
-
-
class
pymor.analyticalproblems.domaindescriptions.LineDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.DomainDescriptionDescribes an interval domain.
Boundary types can be associated edgewise.
Parameters
- domain
List [x_l, x_r] providing the left and right endpoint.
- left
The boundary type of the left endpoint.
- right
The boundary type of the right endpoint.
Attributes
boundary_types,has_dirichlet,has_neumann,has_robin-
domain¶
-
left¶
-
right¶
-
class
pymor.analyticalproblems.domaindescriptions.PolygonalDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.DomainDescriptionDescribes a domain with a polygonal boundary and polygonal holes inside the domain.
Parameters
- points
List of points [x_0, x_1] that describe the polygonal chain that bounds the domain.
- boundary_types
Either a dictionary
{boundary_type: [i_0, ...], boundary_type: [j_0, ...], ...}withi_0, ...being the ids of boundary segments for a given boundary type (0is the line connecting point0to1,1is the line connecting point1to2etc.), or a function that returns the boundary type for a given coordinate.- holes
List of lists of points that describe the polygonal chains that bound the holes inside the domain.
Attributes
boundary_types,has_dirichlet,has_neumann,has_robin-
points¶
-
boundary_types¶
-
holes¶
-
class
pymor.analyticalproblems.domaindescriptions.RectDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.DomainDescriptionDescribes a rectangular domain.
Boundary types can be associated edgewise.
Parameters
- domain
List of two points defining the lower-left and upper-right corner of the domain.
- left
The boundary type of the left edge.
- right
The boundary type of the right edge.
- top
The boundary type of the top edge.
- bottom
The boundary type of the bottom edge.
Attributes
bottom,diameter,dim,domain,height,left,lower_left,right,top,upper_right,volume,widthboundary_types,has_dirichlet,has_neumann,has_robin-
domain¶
-
left¶
-
right¶
-
top¶
-
bottom¶
-
class
pymor.analyticalproblems.domaindescriptions.TorusDomain(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.domaindescriptions.DomainDescriptionDescribes a domain with the topology of a torus.
Parameters
- domain
List of two points defining the lower-left and upper-right corner of the domain. The left and right edge are identified, as well as the bottom and top edge
Attributes
diameter,dim,domain,height,lower_left,upper_right,volume,widthboundary_types,has_dirichlet,has_neumann,has_robin-
domain¶
elliptic module¶
-
class
pymor.analyticalproblems.elliptic.StationaryProblem(*args, **kwargs)[source]¶ Bases:
pymor.parameters.base.ParametricObjectLinear elliptic problem description.
The problem consists in solving
- ∇ ⋅ [d(x, μ) ∇ u(x, μ)] + ∇ ⋅ [f(x, u(x, μ), μ)] + c(x, u(x, μ), μ) = f(x, μ)
for u.
Parameters
- domain
A
DomainDescriptionof the domain the problem is posed on.- rhs
The
Functionf(x, μ).rhs.dim_domainhas to agree with the dimension ofdomain, whereasrhs.shape_rangehas to be().- diffusion
The
Functiond(x, μ) withshape_rangeof either()or(dim domain, dim domain).- advection
The
Functionf, only depending on x, withshape_rangeof(dim domain,).- nonlinear_advection
The
Functionf, only depending on u, withshape_rangeof(dim domain,).- nonlinear_advection_derivative
The derivative of f, only depending on u, with respect to u.
- reaction
The
Functionc, only depending on x, withshape_rangeof().- nonlinear_reaction
The
Functionc, only depending on u, withshape_rangeof().- nonlinear_reaction_derivative
The derivative of the
Functionc, only depending on u, withshape_rangeof().- dirichlet_data
Functionproviding the Dirichlet boundary values.- neumann_data
Functionproviding the Neumann boundary values.- robin_data
Tuple of two
Functionsproviding the Robin parameter and boundary values.- outputs
Tuple of additional output functionals to assemble. Each value must be a tuple of the form
(functional_type, data)wherefunctional_typeis a string defining the type of functional to assemble anddatais aFunctionholding the corresponding coefficient function. Currently implementedfunctional_typesare:- l2
Evaluate the l2-product with the given data function.
- l2_boundary
Evaluate the l2-product with the given data function on the boundary.
- parameter_ranges
Ranges of interest for the
Parametersof the problem.- name
Name of the problem.
Attributes
-
domain¶
-
rhs¶
-
diffusion¶
-
advection¶
-
nonlinear_advection¶
-
nonlinear_advection_derivative¶
-
reaction¶
-
nonlinear_reaction¶
-
nonlinear_reaction_derivative¶
-
dirichlet_data¶
-
neumann_data¶
-
robin_data¶
-
outputs¶
functions module¶
-
class
pymor.analyticalproblems.functions.BitmapFunction(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.functions.FunctionDefine a 2D
Functionvia a grayscale image.Parameters
- filename
Path of the image representing the function.
- bounding_box
Lower left and upper right coordinates of the domain of the function.
- range
A pixel of value p is mapped to
(p / 255.) * range[1] + range[0].
Methods
Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
evaluate(x, mu=None)[source]¶ Evaluate the function for given argument
xandparameter valuesmu.
-
class
pymor.analyticalproblems.functions.ConstantFunction(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.functions.FunctionA constant
Functionf: R^d -> R^shape(c), f(x) = c
Parameters
- value
The constant c.
- dim_domain
The dimension d.
- name
The name of the function.
Methods
Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
evaluate(x, mu=None)[source]¶ Evaluate the function for given argument
xandparameter valuesmu.
-
class
pymor.analyticalproblems.functions.ExpressionFunction(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.functions.GenericFunctionTurns a Python expression given as a string into a
Function.Some
NumPyarithmetic functions like ‘sin’, ‘log’, ‘min’ are supported. For a full list see thefunctionsclass attribute.Warning
evalis used to evaluate the given expression. Using this class with expression strings from untrusted sources will cause mayhem and destruction!Parameters
- expression
A Python expression of one variable
xand a parametermugiven as a string.- dim_domain
The dimension of the domain.
- shape_range
The shape of the values returned by the expression.
- parameters
The
Parametersthe expression accepts.- values
Dictionary of additional constants that can be used in
expressionwith their corresponding value.- name
The name of the function.
Methods
Attributes
functionsparameters,parameters_inherited,parameters_internal,parameters_own,parametric
-
class
pymor.analyticalproblems.functions.Function(*args, **kwargs)[source]¶ Bases:
pymor.parameters.base.ParametricObjectInterface for
Parameterdependent analytical functions.Every
Functionis a map of the formf(μ): Ω ⊆ R^d -> R^(shape_range)
The returned values are
NumPy arraysof arbitrary (but fixed) shape. Note that NumPy distinguishes between one-dimensional arrays of length 1 (with shape(1,)) and zero-dimensional scalar arrays (with shape()). In pyMOR, we usually expect scalar-valued functions to haveshape_range == ().While the function might raise an error if it is evaluated for an argument not in the domain Ω, the exact behavior is left undefined.
Functions are vectorized in the sense, that if
x.ndim == k, thenf(x, μ)[i0, i1, ..., i(k-2)] == f(x[i0, i1, ..., i(k-2)], μ).
In particular,
f(x, μ).shape == x.shape[:-1] + shape_range.Methods
Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
dim_domain¶ The dimension d > 0.
-
shape_range¶ The shape of the function values.
-
abstract
evaluate(x, mu=None)[source]¶ Evaluate the function for given argument
xandparameter valuesmu.
-
-
class
pymor.analyticalproblems.functions.GenericFunction(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.functions.FunctionWrapper making an arbitrary Python function between
NumPy arraysa properFunction.Note that a
GenericFunctioncan only bepickledif the function it is wrapping can be pickled (cf.dumps_function). For this reason, it is usually preferable to useExpressionFunctioninstead ofGenericFunction.Parameters
- mapping
The function to wrap. If
parametersisNone, the function is of the formmapping(x). Ifparametersis notNone, the function has to have the signaturemapping(x, mu). Moreover, the function is expected to be vectorized, i.e.:mapping(x).shape == x.shape[:-1] + shape_range.
- dim_domain
The dimension of the domain.
- shape_range
The shape of the values returned by the mapping.
- parameters
The
Parametersthe mapping accepts.- name
The name of the function.
Methods
Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
evaluate(x, mu=None)[source]¶ Evaluate the function for given argument
xandparameter valuesmu.
-
class
pymor.analyticalproblems.functions.LincombFunction(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.functions.FunctionA
Functionrepresenting a linear combination ofFunctions.The linear coefficients can be provided either as scalars or as
ParameterFunctionals.Parameters
- functions
List of
Functionswhose linear combination is formed.- coefficients
A list of linear coefficients. A linear coefficient can either be a fixed number or a
ParameterFunctional.- name
Name of the function.
Methods
Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
functions¶
-
coefficients¶
-
evaluate(x, mu=None)[source]¶ Evaluate the function for given argument
xandparameter valuesmu.
-
evaluate_coefficients(mu)[source]¶ Compute the linear coefficients for given
parameter valuesmu.
-
class
pymor.analyticalproblems.functions.ProductFunction(*args, **kwargs)[source]¶ Bases:
pymor.analyticalproblems.functions.FunctionA
Functionrepresenting a product ofFunctions.Parameters
- functions
List of
Functionswhose product is formed.- name
Name of the function.
Methods
Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
functions¶
-
evaluate(x, mu=None)[source]¶ Evaluate the function for given argument
xandparameter valuesmu.
helmholtz module¶
-
pymor.analyticalproblems.helmholtz.helmholtz_problem(domain=RectDomain(domain=array([[0, 0], [1, 1]])), rhs=None, parameter_range=0.0, 100.0, dirichlet_data=None, neumann_data=None)[source]¶ Helmholtz equation problem.
This problem is to solve the Helmholtz equation
- ∆ u(x, k) - k^2 u(x, k) = f(x, k)
on a given domain.
Parameters
- domain
A
DomainDescriptionof the domain the problem is posed on.- rhs
The
Functionf(x, μ).- parameter_range
A tuple
(k_min, k_max)describing the interval in which k is allowd to vary.- dirichlet_data
Functionproviding the Dirichlet boundary values.- neumann_data
Functionproviding the Neumann boundary values.
instationary module¶
-
class
pymor.analyticalproblems.instationary.InstationaryProblem(*args, **kwargs)[source]¶ Bases:
pymor.parameters.base.ParametricObjectInstationary problem description.
This class describes an instationary problem of the form
| ∂_t u(x, t, μ) + A(u(x, t, μ), t, μ) = f(x, t, μ), | u(x, 0, μ) = u_0(x, μ)
where A, f are given by the problem’s
stationary_partand t is allowed to vary in the interval [0, T].Parameters
- stationary_part
The stationary part of the problem.
- initial_data
Functionproviding the initial values u_0.- T
The final time T.
- parameter_ranges
Ranges of interest for the
Parametersof the problem.- name
Name of the problem.
Methods
with_stationary_partAttributes
parameter_ranges,parameter_space,stationary_part,Tparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
T¶
-
stationary_part¶
-
parameter_ranges¶
-
name¶
text module¶
thermalblock module¶
-
pymor.analyticalproblems.thermalblock.thermal_block_problem(num_blocks=3, 3, parameter_range=0.1, 1)[source]¶ Analytical description of a 2D ‘thermal block’ diffusion problem.
The problem is to solve the elliptic equation
- ∇ ⋅ [ d(x, μ) ∇ u(x, μ) ] = f(x, μ)
on the domain [0,1]^2 with Dirichlet zero boundary values. The domain is partitioned into nx x ny blocks and the diffusion function d(x, μ) is constant on each such block i with value μ_i.
---------------------------- | | | | | μ_4 | μ_5 | μ_6 | | | | | |--------------------------- | | | | | μ_1 | μ_2 | μ_3 | | | | | ----------------------------
Parameters
- num_blocks
The tuple
(nx, ny)- parameter_range
A tuple
(μ_min, μ_max). EachParametercomponent μ_i is allowed to lie in the interval [μ_min, μ_max].