regpy.vecsps.tensor_bases¶
Classes¶
Consider an evaluation domain given as Tensor product \(D_1\otimes \dots\otimes D_n\) with \(D_1,\dots,D_n\) being \(n\) |
Functions¶
|
Implements a tensor basis of Chebyshev polynomials for product spaces. It requires that |
|
Implements a tensor basis of Legendre polynomials for product spaces. It requires that |
|
Implements a B-Spline basis in an arbitrary Dimension (given by dim) |
Module Contents¶
- class regpy.vecsps.tensor_bases.TensorBasis(coef_domain, eval_domain, bases, dtype=float)[source]¶
Bases:
regpy.operators.OperatorConsider an evaluation domain given as Tensor product \(D_1\otimes \dots\otimes D_n\) with \(D_1,\dots,D_n\) being \(n\) regpy.vecsps.VectorSpace’s and a tensor in the coefficients domain \(V_1\otimes \dots\otimes V_m\) then we define an operator mapping coefficients to some function f: eval_domain -> dtype:
\[f(d_1,...,d_n) = \sum_{k_1=0}^{N_1-1} ... \sum_{k_n=0}^{N_n-1} c_{k_1,...k_n} b^1_{k_1}(x_1) .... b^n_{k_n}(x_n).\]So that the operator TensorBasis maps the coefficient tensor \(c = (c_\{k_1,....k_n\})\) to the tensor of function values :math`(f(x))_{x in eval_domain}`
- Parameters:
eval_domain (regpy.vecsps.Prod) – an instance of the class regpy.vecsps.Prod in vector spaces of size where each D_i has size M_i
coef_domain (regpy.vecsps.Prod) – an instance of the class regpy.vecsps.Prod in vector spaces of size where each V_i has size N_i
bases ([ np.ndarray, ... ]) –
a list of matrices \([B_1,..., B_n]\) where the matrix \(B_l\) of size \(M_l \times N_l\) and contains the function values of the basis \(\{b^l_0, b^l_{M_l-1}\}\) of the l-th coordinate:
\[B_l = (b^l_{k}(x_{l,j}))_{j=0:M_l-1, k=0:N_l-1}\]
- dtype¶
dtype of the vector spaces.
- ndim¶
dimension of the coef_domain.
- bases¶
List of all the bases transforms as a list of np.ndarray’s
- regpy.vecsps.tensor_bases.chebyshev_basis(coef_domain, eval_domain, dtype=float)[source]¶
Implements a tensor basis of Chebyshev polynomials for product spaces. It requires that both coef_domain and eval_domain has the same dimension.
- Parameters:
coef_domain (regpy.vecsps.Prod) – Coefficients of tensor products of Chebychev polynomials.
eval_domain (regpy.vecsps.Prod) – Tensor product of GridFcts instances on which the Chebyshev polynomial are evaluated.
dtype (np.dtype, optional) – type of the underlying spaces, by default float
- Returns:
A bases transform from coefficients of Chebychev polynomial to their evaluation.
- Return type:
TesnorBasis
- regpy.vecsps.tensor_bases.legendre_basis(coef_domain, eval_domain, dtype=float)[source]¶
Implements a tensor basis of Legendre polynomials for product spaces. It requires that both coef_domain and eval_domain has the same dimension.
- Parameters:
coef_domain (regpy.vecsps.Prod) – Coefficients of tensor products of Legendre polynomials.
eval_domain (regpy.vecsps.Prod) – Tensor product of GridFcts instances on which the Legendre polynomial are evaluated.
dtype (np.dtype, optional) – type of the underlying spaces, by default float
- Returns:
A bases transform from coefficients of Legendre polynomial to their evaluation.
- Return type:
TesnorBasis
- regpy.vecsps.tensor_bases.bspline_basis(k, t, dim=1, add_points=10)[source]¶
Implements a B-Spline basis in an arbitrary Dimension (given by dim) the splines are generated via BSpline from scipy.interpolate. In each dimension it uses the knots given in t to generate a B-Spline Basis. The evaluation domain is a refined grid determined by the point added between points given by add_points:
np.linspace(t[0],t[-1],t.size*add_points)
Note, that to do that accurately construct Splines, we use the key extrapolate=False and extend the original knot points given in t by additionally 2k points with equidistant distance to T.
By this construction the splines will be zero at the boundary.
- Parameters:
k (integer) – Smoothness degree of the used Splines.
t (np.ndarray) – One dimensional array of knot points for evaluating the Splines.
dim (int, optional) – Dimension of the resulting spaces, by default 1
add_points (int, optional) – Number of points to be added between each Knot for evaluation, by default 10
- Returns:
A base transform from coefficients of Splines to evaluation on a grid constructed from a refined grid of the given evaluation knots.
- Return type: