FEniCSStepExpansion#
- class cuqipy_fenics.geometry.FEniCSStepExpansion(geometry, num_steps_x=8, num_steps_y=None)#
A geometry class that parameterizes finite element functions on a 1D or 2D domain with piecewise constant functions. This parameterization, which we refer to as step expansion parameterization, is built on the given input geometry that specifies the computational mesh.
In the 1D case, the domain is divided into num_steps_x constant functions, each function support is of length L_x/num_steps_x, where L_x is the length of the domain in the x direction.
In the 2D case, the domain is divided into num_steps_x \(\times\) num_steps_y constant functions each has a support area of L_x/num_steps_x x L_y/num_steps_y, where L_y is the length of the domain in the y direction. The constant functions are arranged in a row-wise order, from the bottom to top of the domain. That is, the first num_steps_x functions are in the bottom row, the next num_steps_x functions are in the next row, and so on.
Note: for accurate results, the underlying mesh should be structured in which the number of nodes in in the x and y directions are multiples of num_steps_x and num_steps_y respectively.
- Parameters:
geometry (cuqipy_fenics.geometry.Geometry) – An input geometry on which the step expansion is built (the geometry must have a mesh attribute)
num_steps_x (int) – Number of step expansion terms in the x direction
num_steps_y (int, optional) – Number of step expansion terms in the y direction, only required for geometries with physical dimension 2
Example
import numpy as np from cuqipy_fenics.geometry import FEniCSStepExpansion, FEniCSContinuous from cuqi.distribution import Gaussian import dolfin as dl mesh = dl.UnitSquareMesh(32,32) V = dl.FunctionSpace(mesh, 'DG', 0) G_FEM = FEniCSContinuous(V, labels=['$\xi_1$', '$\xi_2$']) G_step = FEniCSStepExpansion(G_FEM, num_steps_x=4, num_steps_y=8) x = Gaussian(mean=np.zeros(G_step.par_dim), cov=np.eye(G_step.par_dim), geometry=G_step) samples = x.sample() samples.plot()
- __init__(geometry, num_steps_x=8, num_steps_y=None)#
Methods
__init__
(geometry[, num_steps_x, num_steps_y])fun2par
(f)The function to parameter map used to map function values back to parameters, if available.
fun2vec
(fun)Maps the function value (FEniCS object) to the corresponding vector representation of the function (ndarray of the function DOF values).
gradient
(direction, wrt)par2field
(p)Applies linear transformation of the parameters p to generate a realization of the step expansion
par2fun
(p)The parameter to function map used to map parameters to function values in e.g. plotting.
plot
(values[, is_par, plot_par])Plots a function over the set defined by the geometry object.
plot_envelope
(lo_values, hi_values[, ...])Plots an envelope from lower and upper bounds over the set defined by the geometry object.
vec2fun
(funvec)Maps the vector representation of the function (ndarray of the function DOF values) to the function value (FEniCS object).
Attributes
The dimension of the geometry (function space).
Flag to indicate whether the function value is an array.
The shape of the geometry (function space).
The dimension of the geometry (dimension of the vector representation of the function value).
The shape of the geometry (shape of the vector representation of the function value).
The dimension of the geometry (parameter space).
The shape of the parameter space.
Returns the physical dimension of the geometry, e.g. 1, 2 or 3.