fem.tensor#

We assume here that a tensor space is the product of fem spaces whom basis are of compact support

Classes#

Inheritance diagram of psydac.fem.tensor

TensorFemSpace(domain_decomposition, *spaces)

Tensor-product Finite Element space V.

Details#

We assume here that a tensor space is the product of fem spaces whom basis are of compact support

class TensorFemSpace(domain_decomposition, *spaces, vector_space=None, cart=None, dtype=<class 'float'>)[source]#

Bases: FemSpace

Tensor-product Finite Element space V.

Parameters:
domain_decompositionpsydac.ddm.cart.DomainDecomposition
*spacespsydac.fem.splines.SplineSpace

1D finite element spaces.

vector_spacepsydac.linalg.stencil.StencilVectorSpace or None

The vector space to which the coefficients belong (optional).

cartpsydac.ddm.CartDecomposition or None

Object that contains all information about the Cartesian decomposition of a tensor-product grid of coefficients.

dtype{float, complex}

Data type of the coefficients.

Notes

For now we assume that this tensor-product space can ONLY be constructed from 1D spline spaces.

property ldim#

Parametric dimension.

property periodic#

Tuple of booleans: along each logical dimension, say if domain is periodic.

property domain_decomposition#
property mapping#

Mapping from logical coordinates ‘eta’ to physical coordinates ‘x’. If None, we assume identity mapping (hence x=eta).

property vector_space#

Returns the topological associated vector space.

property is_product#

Boolean flag that describes whether the space is a product space. If True, an element of this space can be decomposed into separate fields.

property interfaces#
property symbolic_space#

Symbolic space.

eval_field(field, *eta, weights=None)[source]#

Evaluate field at location(s) eta.

Parameters:
fieldFemField

Field object (element of FemSpace) to be evaluated.

etalist of float or numpy.ndarray

Evaluation point(s) in logical domain.

weightsStencilVector, optional

Weights of the basis functions, such that weights.space == field.coeffs.space.

Returns:
valuefloat or numpy.ndarray

Field value(s) at location(s) eta.

preprocess_regular_tensor_grid(grid, der=0, overlap=0)[source]#

Returns all the quantities needed to evaluate fields on a regular tensor-grid.

Parameters:
gridList of ndarray

List of 2D arrays representing each direction of the grid. Each of these arrays should have shape (ne_xi, nv_xi) where ne_xi is the number of cells in the domain in the direction xi and nv_xi is the number of evaluation points in the same direction.

derint, default=0

Number of derivatives of the basis functions to pre-compute.

overlapint

How much to overlap. Only used in the distributed context.

Returns:
degreetuple of int

Degree in each direction

global_basisList of ndarray

List of 4D arrays, one per direction, containing the values of the p+1 non-vanishing basis functions (and their derivatives) at each grid point. The array for direction xi has shape (ne_xi, der + 1, p+1, nv_xi).

global_spansList of ndarray

List of 1D arrays, one per direction, containing the index of the last non-vanishing basis function in each cell. The array for direction xi has shape (ne_xi,).

local_shapeList of tuple

Shape of what is local to this instance.

preprocess_irregular_tensor_grid(grid, der=0, overlap=0)[source]#

Returns all the quantities needed to evaluate fields on a regular tensor-grid.

Parameters:
gridList of ndarray

List of 1D arrays representing each direction of the grid.

derint, default=0

Number of derivatives of the basis functions to pre-compute.

overlapint

How much to overlap. Only used in the distributed context.

Returns:
padstuple of int

Padding in each direction

degreetuple of int

Degree in each direction

global_basisList of ndarray

List of 4D arrays, one per direction, containing the values of the p+1 non-vanishing basis functions (and their derivatives) at each grid point. The array for direction xi has shape (n_xi, p+1, der + 1).

global_spansList of ndarray

List of 1D arrays, one per direction, containing the index of the last non-vanishing basis function in each cell. The array for direction xi has shape (n_xi,).

cell_indexeslist of ndarray

List of 1D arrays, one per direction, containing the index of the cell in which the corresponding point in grid is.

local_shapeList of tuple

Shape of what is local to this instance.

eval_fields(grid, *fields, weights=None, npts_per_cell=None, overlap=0)[source]#

Evaluate one or several fields at the given location(s) grid.

Parameters:
gridList of ndarray

Grid on which to evaluate the fields

*fieldstuple of psydac.fem.basic.FemField

Fields to evaluate

weightspsydac.fem.basic.FemField or None, optional

Weights field.

npts_per_cell: int or tuple of int or None, optional

number of evaluation points in each cell. If an integer is given, then assume that it is the same in every direction.

overlapint

How much to overlap. Only used in the distributed context.

Returns:
List of ndarray of floats

List of the evaluated fields.

eval_fields_regular_tensor_grid(grid, *fields, weights=None, overlap=0)[source]#

Evaluate fields on a regular tensor grid

Parameters:
gridList of ndarray

List of 2D arrays representing each direction of the grid. Each of these arrays should have shape (ne_xi, nv_xi) where ne is the number of cells in the domain in the direction xi and nv_xi is the number of evaluation points in the same direction.

*fieldstuple of psydac.fem.basic.FemField

Fields to evaluate on grid.

weightspsydac.fem.basic.FemField or None, optional

Weights to apply to our fields.

overlapint

How much to overlap. Only used in the distributed context.

Returns:
List of ndarray of float

Values of the fields on the regular tensor grid

eval_fields_irregular_tensor_grid(grid, *fields, weights=None, overlap=0)[source]#

Evaluate fields on a regular tensor grid

Parameters:
gridList of ndarray

List of 2D arrays representing each direction of the grid. Each of these arrays should have shape (ne_xi, nv_xi) where ne is the number of cells in the domain in the direction xi and nv_xi is the number of evaluation points in the same direction.

*fieldstuple of psydac.fem.basic.FemField

Fields to evaluate on grid.

weightspsydac.fem.basic.FemField or None, optional

Weights to apply to our fields.

overlapint

How much to overlap. Only used in the distributed context.

Returns:
List of ndarray of float

Values of the fields on the regular tensor grid

eval_field_gradient(field, *eta, weights=None)[source]#

Evaluate field gradient at location(s) eta.

Parameters:
fieldFemField

Field object (element of FemSpace) to be evaluated.

etalist of float or numpy.ndarray

Evaluation point(s) in logical domain.

weightsStencilVector, optional

Weights of the basis functions, such that weights.space == field.coeffs.space.

Returns:
valuefloat or numpy.ndarray

Value(s) of field gradient at location(s) eta.

integral(f, *, nquads=None)[source]#
property is_scalar#
property dtype#
property nbasis#
property knots#
property breaks#
property degree#
property multiplicity#
property pads#
property ncells#
property spaces#
get_assembly_grids(*nquads)[source]#

Return a tuple of FemAssemblyGrid objects (one for each direction).

These objects are local to the process, and contain all 1D information which is necessary for the correct assembly of the l.h.s. matrix and r.h.s. vector in a finite element method. This information includes the coordinates and weights of all quadrature points, as well as the values of the non-zero basis functions, and their derivatives, at such points.

The computed FemAssemblyGrid objects are stored in a dictionary in self with nquads as key, and are reused if a match is found.

Parameters:
*nquadsint

Number of quadrature points per cell, along each direction.

Returns:
tuple of FemAssemblyGrid

The 1D assembly grids along each direction.

property local_domain#

Logical domain local to the process, assuming the global domain is decomposed across processes without any overlapping.

This information is fundamental for avoiding double-counting when computing integrals over the global domain.

Returns:
element_startstuple of int

Start element index along each direction.

element_endstuple of int

End element index along each direction.

property global_element_starts#
property global_element_ends#
property eta_lims#

Eta limits of domain local to the process (for field evaluation).

Returns:
eta_limits: tuple of (2-tuple of float)

Along each dimension i, limits are given as (eta^i_{min}, eta^i_{max}).

init_interpolation()[source]#
init_histopolation()[source]#
compute_interpolant(values, field)[source]#

Compute field (i.e. update its spline coefficients) such that it interpolates a certain function \(f(x1,x2,..)\) at the Greville points.

Parameters:
valuesStencilVector

Function values \(f(x_i)\) at the n-dimensional tensor grid of Greville points \(x_i\), to be interpolated.

fieldFemField

Input/output argument: tensor spline that has to interpolate the given values.

reduce_grid(axes=(), knots=())[source]#

Create a new TensorFemSpace object with a coarser grid than the original one we do that by giving a new knot sequence in the desired dimension.

Parameters:
axesList of int

Dimensions where we want to coarsen the grid.

knotsList or tuple

New knot sequences in each dimension.

Returns:
VTensorFemSpace

New space with a coarser grid.

export_fields(filename, **fields)[source]#

Write spline coefficients of given fields to HDF5 file.

import_fields(filename, *field_names)[source]#

Load fields from HDF5 file containing spline coefficients.

Parameters:
filenamestr

Name of HDF5 input file.

field_nameslist of str

Names of the datasets with the required spline coefficients.

Returns:
fieldslist of FemSpace objects

Distributed fields, given in the same order of the names.

reduce_degree(axes, multiplicity=None, basis='B')[source]#
add_refined_space(ncells)[source]#

refine the space with new ncells and add it to the dictionary of refined_space

create_interface_space(axis, ext, cart)[source]#

Create a new interface fem space along a given axis and extremity.

Parameters:
axisint

The axis of the new Interface space.

ext: int

The extremity of the new Interface space. the values must be 1 or -1.

cart: CartDecomposition

The cart of the new space, needed in the parallel case.

get_refined_space(ncells)[source]#
set_refined_space(ncells, new_space)[source]#
plot_2d_decomposition(mapping=None, refine=10)[source]#