DOLFINx
0.1.0
DOLFINx C++ interface
|
Finite element method functionality. More...
Namespaces | |
sparsitybuild | |
Functions to build sparsity patterns from degree-of-freedom maps. | |
Classes | |
class | Constant |
A constant value which can be attached to a Form. Constants may be scalar (rank 0), vector (rank 1), or tensor valued. More... | |
class | CoordinateElement |
This class manages coordinate mappings for isoparametric cells. More... | |
class | DirichletBC |
Interface for setting (strong) Dirichlet boundary conditions. More... | |
class | DofMap |
Degree-of-freedom map. More... | |
class | ElementDofLayout |
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh entity. This class also handles sub-space dofs, which are views into the parent dofs. More... | |
class | Expression |
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell. This class closely follows the concept of a UFC Expression. More... | |
class | FiniteElement |
Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis. More... | |
class | Form |
Class for variational forms. More... | |
class | Function |
This class represents a function \( u_h \) in a finite element function space \( V_h \), given by. More... | |
class | FunctionSpace |
This class represents a finite element function space defined by a mesh, a finite element, and a local-to-global map of the degrees of freedom (dofmap). More... | |
Enumerations | |
enum | IntegralType : std::int8_t { cell = 0, exterior_facet = 1, interior_facet = 2, vertex = 3 } |
Type of integral. | |
Functions | |
template<typename T > | |
T | assemble_scalar (const Form< T > &M) |
Assemble functional into scalar. Caller is responsible for accumulation across processes. More... | |
template<typename T > | |
void | assemble_vector (xtl::span< T > b, const Form< T > &L) |
Assemble linear form into a vector. More... | |
template<typename T > | |
void | apply_lifting (xtl::span< T > b, const std::vector< std::shared_ptr< const Form< T >>> &a, const std::vector< std::vector< std::shared_ptr< const DirichletBC< T >>>> &bcs1, const std::vector< xtl::span< const T >> &x0, double scale) |
Modify b such that: More... | |
template<typename T > | |
void | assemble_matrix (const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> &mat_add, const Form< T > &a, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs) |
Assemble bilinear form into a matrix. More... | |
template<typename T > | |
void | assemble_matrix (const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> &mat_add, const Form< T > &a, const std::vector< bool > &dof_marker0, const std::vector< bool > &dof_marker1) |
Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix. More... | |
template<typename T > | |
void | set_diagonal (const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> &set_fn, const xtl::span< const std::int32_t > &rows, T diagonal=1.0) |
Sets a value to the diagonal of a matrix for specified rows. It is typically called after assembly. The assembly function zeroes Dirichlet rows and columns. For block matrices, this function should normally be called only on the diagonal blocks, i.e. blocks for which the test and trial spaces are the same. More... | |
template<typename T > | |
void | set_diagonal (const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> &set_fn, const fem::FunctionSpace &V, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs, T diagonal=1.0) |
Sets a value to the diagonal of the matrix for rows with a Dirichlet boundary conditions applied. This function is typically called after assembly. The assembly function zeroes Dirichlet rows and columns. This function adds the value only to rows that are locally owned, and therefore does not create a need for parallel communication. For block matrices, this function should normally be called only on the diagonal blocks, i.e. blocks for which the test and trial spaces are the same. More... | |
template<typename T > | |
void | set_bc (xtl::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs, const xtl::span< const T > &x0, double scale=1.0) |
Set bc values in owned (local) part of the vector, multiplied by 'scale'. The vectors b and x0 must have the same local size. The bcs should be on (sub-)spaces of the form L that b represents. | |
template<typename T > | |
void | set_bc (xtl::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs, double scale=1.0) |
Set bc values in owned (local) part of the vector, multiplied by 'scale'. The bcs should be on (sub-)spaces of the form L that b represents. | |
template<typename T > | |
std::vector< std::vector< std::shared_ptr< const fem::DirichletBC< T > > > > | bcs_rows (const std::vector< const Form< T > * > &L, const std::vector< std::shared_ptr< const fem::DirichletBC< T >>> &bcs) |
Arrange boundary conditions by block. More... | |
template<typename T > | |
std::vector< std::vector< std::vector< std::shared_ptr< const fem::DirichletBC< T > > > > > | bcs_cols (const std::vector< std::vector< std::shared_ptr< const Form< T >>>> &a, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs) |
Arrange boundary conditions by block. More... | |
std::array< std::vector< std::int32_t >, 2 > | locate_dofs_topological (const std::array< std::reference_wrapper< const fem::FunctionSpace >, 2 > &V, const int dim, const xtl::span< const std::int32_t > &entities, bool remote=true) |
Find degrees-of-freedom which belong to the provided mesh entities (topological). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex. More... | |
std::vector< std::int32_t > | locate_dofs_topological (const fem::FunctionSpace &V, const int dim, const xtl::span< const std::int32_t > &entities, bool remote=true) |
Find degrees-of-freedom which belong to the provided mesh entities (topological). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex. More... | |
std::array< std::vector< std::int32_t >, 2 > | locate_dofs_geometrical (const std::array< std::reference_wrapper< const fem::FunctionSpace >, 2 > &V, const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &marker_fn) |
Finds degrees of freedom whose geometric coordinate is true for the provided marking function. More... | |
std::vector< std::int32_t > | locate_dofs_geometrical (const fem::FunctionSpace &V, const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &marker_fn) |
Finds degrees of freedom whose geometric coordinate is true for the provided marking function. More... | |
la::SparsityPattern | create_sparsity_discrete_gradient (const fem::FunctionSpace &V0, const fem::FunctionSpace &V1) |
template<typename T > | |
void | assemble_discrete_gradient (const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> &mat_set, const fem::FunctionSpace &V0, const fem::FunctionSpace &V1) |
graph::AdjacencyList< std::int32_t > | transpose_dofmap (const graph::AdjacencyList< std::int32_t > &dofmap, std::int32_t num_cells) |
Create an adjacency list that maps a global index (process-wise) to the 'unassembled' cell-wise contributions. It is built from the usual (cell, local index) -> global index dof map. An 'unassembled' vector is the stacked cell contributions, ordered by cell index. More... | |
std::tuple< std::shared_ptr< common::IndexMap >, int, graph::AdjacencyList< std::int32_t > > | build_dofmap_data (MPI_Comm comm, const mesh::Topology &topology, const ElementDofLayout &element_dof_layout) |
Build dofmap data for an element on a mesh topology. More... | |
template<typename T > | |
void | eval (array2d< T > &values, const fem::Expression< T > &e, const xtl::span< const std::int32_t > &active_cells) |
Evaluate a UFC expression. More... | |
std::array< std::vector< std::shared_ptr< const FunctionSpace > >, 2 > | common_function_spaces (const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace >, 2 >>> &V) |
Extract FunctionSpaces for (0) rows blocks and (1) columns blocks from a rectangular array of (test, trial) space pairs. The test space must be the same for each row and the trial spaces must be the same for each column. Raises an exception if there is an inconsistency. e.g. if each form in row i does not have the same test space then an exception is raised. More... | |
xt::xtensor< double, 2 > | interpolation_coords (const fem::FiniteElement &element, const mesh::Mesh &mesh, const xtl::span< const std::int32_t > &cells) |
Compute the evaluation points in the physical space at which an expression should be computed to interpolate it in a finite elemenet space. More... | |
template<typename T > | |
void | interpolate (Function< T > &u, const Function< T > &v) |
Interpolate a finite element Function (on possibly non-matching meshes) in another finite element space. More... | |
template<typename T > | |
void | interpolate (Function< T > &u, const std::function< xt::xarray< T >(const xt::xtensor< double, 2 > &)> &f, const xt::xtensor< double, 2 > &x, const xtl::span< const std::int32_t > &cells) |
Interpolate an expression in a finite element space. More... | |
template<typename T > | |
void | interpolate_c (Function< T > &u, const std::function< void(xt::xarray< T > &, const xt::xtensor< double, 2 > &)> &f, const xt::xtensor< double, 2 > &x, const xtl::span< const std::int32_t > &cells) |
Interpolate an expression f(x) More... | |
Mat | create_matrix (const Form< PetscScalar > &a, const std::string &type=std::string()) |
Create a matrix. More... | |
Mat | create_matrix_block (const std::vector< std::vector< const fem::Form< PetscScalar > * >> &a, const std::string &type=std::string()) |
Initialise a monolithic matrix for an array of bilinear forms. More... | |
Mat | create_matrix_nest (const std::vector< std::vector< const fem::Form< PetscScalar > * >> &a, const std::vector< std::vector< std::string >> &types) |
Create nested (MatNest) matrix. More... | |
Vec | create_vector_block (const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int >> &maps) |
Initialise monolithic vector. Vector is not zeroed. More... | |
Vec | create_vector_nest (const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int >> &maps) |
Create nested (VecNest) vector. Vector is not zeroed. | |
void | assemble_vector_petsc (Vec b, const Form< PetscScalar > &L) |
Assemble linear form into an already allocated PETSc vector. Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End. More... | |
void | apply_lifting_petsc (Vec b, const std::vector< std::shared_ptr< const Form< PetscScalar >>> &a, const std::vector< std::vector< std::shared_ptr< const DirichletBC< PetscScalar >>>> &bcs1, const std::vector< Vec > &x0, double scale) |
Modify b such that: More... | |
void | set_bc_petsc (Vec b, const std::vector< std::shared_ptr< const DirichletBC< PetscScalar >>> &bcs, const Vec x0, double scale=1.0) |
Set bc values in owned (local) part of the PETScVector, multiplied by 'scale'. The vectors b and x0 must have the same local size. The bcs should be on (sub-)spaces of the form L that b represents. | |
template<typename T > | |
std::vector< std::vector< std::array< std::shared_ptr< const fem::FunctionSpace >, 2 > > > | extract_function_spaces (const std::vector< std::vector< const fem::Form< T > * >> &a) |
Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array of forms. More... | |
template<typename T > | |
la::SparsityPattern | create_sparsity_pattern (const Form< T > &a) |
Create a sparsity pattern for a given form. The pattern is not finalised, i.e. the caller is responsible for calling SparsityPattern::assemble. More... | |
la::SparsityPattern | create_sparsity_pattern (const mesh::Topology &topology, const std::array< const std::reference_wrapper< const fem::DofMap >, 2 > &dofmaps, const std::set< IntegralType > &integrals) |
Create a sparsity pattern for a given form. The pattern is not finalised, i.e. the caller is responsible for calling SparsityPattern::assemble. | |
ElementDofLayout | create_element_dof_layout (const ufc_dofmap &dofmap, const mesh::CellType cell_type, const std::vector< int > &parent_map={}) |
Create an ElementDofLayout from a ufc_dofmap. | |
DofMap | create_dofmap (MPI_Comm comm, const ufc_dofmap &dofmap, mesh::Topology &topology) |
Create dof map on mesh from a ufc_dofmap. More... | |
std::vector< std::string > | get_coefficient_names (const ufc_form &ufc_form) |
Get the name of each coefficient in a UFC form. More... | |
std::vector< std::string > | get_constant_names (const ufc_form &ufc_form) |
Get the name of each constant in a UFC form. More... | |
template<typename T > | |
Form< T > | create_form (const ufc_form &ufc_form, const std::vector< std::shared_ptr< const fem::FunctionSpace >> &spaces, const std::vector< std::shared_ptr< const fem::Function< T >>> &coefficients, const std::vector< std::shared_ptr< const fem::Constant< T >>> &constants, const std::map< IntegralType, const mesh::MeshTags< int > * > &subdomains, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr) |
Create a Form from UFC input. More... | |
template<typename T > | |
Form< T > | create_form (const ufc_form &ufc_form, const std::vector< std::shared_ptr< const fem::FunctionSpace >> &spaces, const std::map< std::string, std::shared_ptr< const fem::Function< T >>> &coefficients, const std::map< std::string, std::shared_ptr< const fem::Constant< T >>> &constants, const std::map< IntegralType, const mesh::MeshTags< int > * > &subdomains, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr) |
Create a Form from UFC input. More... | |
template<typename T > | |
std::shared_ptr< Form< T > > | create_form (ufc_form *(*fptr)(), const std::vector< std::shared_ptr< const fem::FunctionSpace >> &spaces, const std::map< std::string, std::shared_ptr< const fem::Function< T >>> &coefficients, const std::map< std::string, std::shared_ptr< const fem::Constant< T >>> &constants, const std::map< IntegralType, const mesh::MeshTags< int > * > &subdomains, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr) |
Create a Form using a factory function that returns a pointer to a ufc_form. More... | |
std::shared_ptr< fem::FunctionSpace > | create_functionspace (ufc_function_space *(*fptr)(const char *), const std::string function_name, std::shared_ptr< mesh::Mesh > mesh) |
Create FunctionSpace from UFC. More... | |
template<typename U > | |
array2d< typename U::scalar_type > | pack_coefficients (const U &u) |
Pack coefficients of u of generic type U ready for assembly. | |
template<typename U > | |
std::vector< typename U::scalar_type > | pack_constants (const U &u) |
Pack constants of u of generic type U ready for assembly. | |
Finite element method functionality.
Classes and algorithms for finite element method spaces and operations.
void dolfinx::fem::apply_lifting | ( | xtl::span< T > | b, |
const std::vector< std::shared_ptr< const Form< T >>> & | a, | ||
const std::vector< std::vector< std::shared_ptr< const DirichletBC< T >>>> & | bcs1, | ||
const std::vector< xtl::span< const T >> & | x0, | ||
double | scale | ||
) |
Modify b such that:
b <- b - scale * A_j (g_j - x0_j)
where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs1 are on the trial spaces V_j. The forms in [a] must have the same test space as L (from which b was built), but the trial space may differ. If x0 is not supplied, then it is treated as zero.
Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.
void dolfinx::fem::apply_lifting_petsc | ( | Vec | b, |
const std::vector< std::shared_ptr< const Form< PetscScalar >>> & | a, | ||
const std::vector< std::vector< std::shared_ptr< const DirichletBC< PetscScalar >>>> & | bcs1, | ||
const std::vector< Vec > & | x0, | ||
double | scale | ||
) |
Modify b such that:
b <- b - scale * A_j (g_j - x0_j)
where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs1 are on the trial spaces V_j. The forms in [a] must have the same test space as L (from which b was built), but the trial space may differ. If x0 is not supplied, then it is treated as zero.
Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.
void dolfinx::fem::assemble_discrete_gradient | ( | const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> & | mat_set, |
const fem::FunctionSpace & | V0, | ||
const fem::FunctionSpace & | V1 | ||
) |
Build the discrete gradient operator A that takes a \(w \in H^1\) (P1, nodal Lagrange) to \(v \in H(curl)\) (lowest order Nedelec), i.e. v = Aw. V0 is the H(curl) space, and V1 is the P1 Lagrange space.
[in] | mat_set | A function (or lambda capture) to set values in a matrix |
[in] | V0 | A H(curl) space |
[in] | V1 | A P1 Lagrange space |
void dolfinx::fem::assemble_matrix | ( | const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> & | mat_add, |
const Form< T > & | a, | ||
const std::vector< bool > & | dof_marker0, | ||
const std::vector< bool > & | dof_marker1 | ||
) |
Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix.
[in] | mat_add | The function for adding values into the matrix |
[in] | a | The bilinear form to assemble |
[in] | dof_marker0 | Boundary condition markers for the rows. If bc[i] is true then rows i in A will be zeroed. The index i is a local index. |
[in] | dof_marker1 | Boundary condition markers for the columns. If bc[i] is true then rows i in A will be zeroed. The index i is a local index. |
void dolfinx::fem::assemble_matrix | ( | const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> & | mat_add, |
const Form< T > & | a, | ||
const std::vector< std::shared_ptr< const DirichletBC< T >>> & | bcs | ||
) |
Assemble bilinear form into a matrix.
[in] | mat_add | The function for adding values into the matrix |
[in] | a | The bilinear from to assemble |
[in] | bcs | Boundary conditions to apply. For boundary condition dofs the row and column are zeroed. The diagonal entry is not set. |
T dolfinx::fem::assemble_scalar | ( | const Form< T > & | M | ) |
Assemble functional into scalar. Caller is responsible for accumulation across processes.
[in] | M | The form (functional) to assemble |
void dolfinx::fem::assemble_vector | ( | xtl::span< T > | b, |
const Form< T > & | L | ||
) |
Assemble linear form into a vector.
[in,out] | b | The vector to be assembled. It will not be zeroed before assembly. |
[in] | L | The linear forms to assemble into b |
void dolfinx::fem::assemble_vector_petsc | ( | Vec | b, |
const Form< PetscScalar > & | L | ||
) |
Assemble linear form into an already allocated PETSc vector. Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.
[in,out] | b | The PETsc vector to assemble the form into. The vector must already be initialised with the correct size. The process-local contribution of the form is assembled into this vector. It is not zeroed before assembly. |
[in] | L | The linear form to assemble |
std::vector< std::vector<std::vector<std::shared_ptr<const fem::DirichletBC<T> > > > > dolfinx::fem::bcs_cols | ( | const std::vector< std::vector< std::shared_ptr< const Form< T >>>> & | a, |
const std::vector< std::shared_ptr< const DirichletBC< T >>> & | bcs | ||
) |
Arrange boundary conditions by block.
[in] | a | Biinear forms for each block |
[in] | bcs | Boundary conditions |
std::vector<std::vector<std::shared_ptr<const fem::DirichletBC<T> > > > dolfinx::fem::bcs_rows | ( | const std::vector< const Form< T > * > & | L, |
const std::vector< std::shared_ptr< const fem::DirichletBC< T >>> & | bcs | ||
) |
Arrange boundary conditions by block.
[in] | L | Linear forms for each block |
[in] | bcs | Boundary conditions |
std::tuple< std::shared_ptr< common::IndexMap >, int, graph::AdjacencyList< std::int32_t > > dolfinx::fem::build_dofmap_data | ( | MPI_Comm | comm, |
const mesh::Topology & | topology, | ||
const ElementDofLayout & | element_dof_layout | ||
) |
Build dofmap data for an element on a mesh topology.
[in] | comm | MPI communicator |
[in] | topology | The mesh topology |
[in] | element_dof_layout | The element dof layout for the function space |
std::array< std::vector< std::shared_ptr< const fem::FunctionSpace > >, 2 > dolfinx::fem::common_function_spaces | ( | const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace >, 2 >>> & | V | ) |
Extract FunctionSpaces for (0) rows blocks and (1) columns blocks from a rectangular array of (test, trial) space pairs. The test space must be the same for each row and the trial spaces must be the same for each column. Raises an exception if there is an inconsistency. e.g. if each form in row i does not have the same test space then an exception is raised.
[in] | V | Vector function spaces for (0) each row block and (1) each column block |
fem::DofMap dolfinx::fem::create_dofmap | ( | MPI_Comm | comm, |
const ufc_dofmap & | dofmap, | ||
mesh::Topology & | topology | ||
) |
Create dof map on mesh from a ufc_dofmap.
[in] | comm | MPI communicator |
[in] | dofmap | The ufc_dofmap |
[in] | topology | The mesh topology |
Form<T> dolfinx::fem::create_form | ( | const ufc_form & | ufc_form, |
const std::vector< std::shared_ptr< const fem::FunctionSpace >> & | spaces, | ||
const std::map< std::string, std::shared_ptr< const fem::Function< T >>> & | coefficients, | ||
const std::map< std::string, std::shared_ptr< const fem::Constant< T >>> & | constants, | ||
const std::map< IntegralType, const mesh::MeshTags< int > * > & | subdomains, | ||
const std::shared_ptr< const mesh::Mesh > & | mesh = nullptr |
||
) |
Create a Form from UFC input.
[in] | ufc_form | The UFC form |
[in] | spaces | The function spaces for the Form arguments |
[in] | coefficients | Coefficient fields in the form (by name) |
[in] | constants | Spatial constants in the form (by name) |
[in] | subdomains | Subdomain makers |
[in] | mesh | The mesh of the domain. This is required if the form has no arguments, e.g. a functional. |
Form<T> dolfinx::fem::create_form | ( | const ufc_form & | ufc_form, |
const std::vector< std::shared_ptr< const fem::FunctionSpace >> & | spaces, | ||
const std::vector< std::shared_ptr< const fem::Function< T >>> & | coefficients, | ||
const std::vector< std::shared_ptr< const fem::Constant< T >>> & | constants, | ||
const std::map< IntegralType, const mesh::MeshTags< int > * > & | subdomains, | ||
const std::shared_ptr< const mesh::Mesh > & | mesh = nullptr |
||
) |
Create a Form from UFC input.
[in] | ufc_form | The UFC form |
[in] | spaces | Vector of function spaces |
[in] | coefficients | Coefficient fields in the form |
[in] | constants | Spatial constants in the form |
[in] | subdomains | Subdomain markers |
[in] | mesh | The mesh of the domain |
std::shared_ptr<Form<T> > dolfinx::fem::create_form | ( | ufc_form *(*)() | fptr, |
const std::vector< std::shared_ptr< const fem::FunctionSpace >> & | spaces, | ||
const std::map< std::string, std::shared_ptr< const fem::Function< T >>> & | coefficients, | ||
const std::map< std::string, std::shared_ptr< const fem::Constant< T >>> & | constants, | ||
const std::map< IntegralType, const mesh::MeshTags< int > * > & | subdomains, | ||
const std::shared_ptr< const mesh::Mesh > & | mesh = nullptr |
||
) |
Create a Form using a factory function that returns a pointer to a ufc_form.
[in] | fptr | pointer to a function returning a pointer to ufc_form |
[in] | spaces | The function spaces for the Form arguments |
[in] | coefficients | Coefficient fields in the form (by name) |
[in] | constants | Spatial constants in the form (by name) |
[in] | subdomains | Subdomain markers |
[in] | mesh | The mesh of the domain. This is required if the form has no arguments, e.g. a functional. |
std::shared_ptr< fem::FunctionSpace > dolfinx::fem::create_functionspace | ( | ufc_function_space *(*)(const char *) | fptr, |
const std::string | function_name, | ||
std::shared_ptr< mesh::Mesh > | mesh | ||
) |
Create FunctionSpace from UFC.
[in] | fptr | Function Pointer to a ufc_function_space_create function |
[in] | function_name | Name of a function whose function space to create. Function name is the name of Python variable for ufl.Coefficient, ufl.TrialFunction or ufl.TestFunction as defined in the UFL file. |
[in] | mesh | Mesh |
Mat dolfinx::fem::create_matrix | ( | const Form< PetscScalar > & | a, |
const std::string & | type = std::string() |
||
) |
Create a matrix.
[in] | a | A bilinear form |
[in] | type | The PETSc matrix type to create |
Mat dolfinx::fem::create_matrix_block | ( | const std::vector< std::vector< const fem::Form< PetscScalar > * >> & | a, |
const std::string & | type = std::string() |
||
) |
Initialise a monolithic matrix for an array of bilinear forms.
[in] | a | Rectangular array of bilinear forms. The a(i, j) form will correspond to the (i, j) block in the returned matrix |
[in] | type | The type of PETSc Mat. If empty the PETSc default is used. |
Mat dolfinx::fem::create_matrix_nest | ( | const std::vector< std::vector< const fem::Form< PetscScalar > * >> & | a, |
const std::vector< std::vector< std::string >> & | types | ||
) |
Create nested (MatNest) matrix.
The caller is responsible for destroying the Mat object
la::SparsityPattern dolfinx::fem::create_sparsity_discrete_gradient | ( | const fem::FunctionSpace & | V0, |
const fem::FunctionSpace & | V1 | ||
) |
Build the sparsity for the discrete gradient operator A that takes a \(w \in H^1\) (P1, nodal Lagrange) to \(v \in H(curl)\) (lowest order Nedelec), i.e. v = Aw. V0 is the H(curl) space, and V1 is the P1 Lagrange space.
[in] | V0 | A H(curl) space |
[in] | V1 | A P1 Lagrange space |
la::SparsityPattern dolfinx::fem::create_sparsity_pattern | ( | const Form< T > & | a | ) |
Create a sparsity pattern for a given form. The pattern is not finalised, i.e. the caller is responsible for calling SparsityPattern::assemble.
[in] | a | A bilinear form |
Vec dolfinx::fem::create_vector_block | ( | const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int >> & | maps | ) |
Initialise monolithic vector. Vector is not zeroed.
The caller is responsible for destroying the Mat object
void dolfinx::fem::eval | ( | array2d< T > & | values, |
const fem::Expression< T > & | e, | ||
const xtl::span< const std::int32_t > & | active_cells | ||
) |
Evaluate a UFC expression.
[out] | values | An array to evaluate the expression into |
[in] | e | The expression to evaluate |
[in] | active_cells | The cells on which to evaluate the expression |
std::vector< std::vector<std::array<std::shared_ptr<const fem::FunctionSpace>, 2> > > dolfinx::fem::extract_function_spaces | ( | const std::vector< std::vector< const fem::Form< T > * >> & | a | ) |
Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array of forms.
[in] | a | A rectangular block on bilinear forms |
a
with a pair of function spaces in each array entry. If a form is null, then the returned function space pair is (null, null). std::vector< std::string > dolfinx::fem::get_coefficient_names | ( | const ufc_form & | ufc_form | ) |
Get the name of each coefficient in a UFC form.
[in] | ufc_form | The UFC form return The name of each coefficient |
std::vector< std::string > dolfinx::fem::get_constant_names | ( | const ufc_form & | ufc_form | ) |
Get the name of each constant in a UFC form.
[in] | ufc_form | The UFC form return The name of each constant |
Interpolate a finite element Function (on possibly non-matching meshes) in another finite element space.
[out] | u | The function to interpolate into |
[in] | v | The function to be interpolated |
void dolfinx::fem::interpolate | ( | Function< T > & | u, |
const std::function< xt::xarray< T >(const xt::xtensor< double, 2 > &)> & | f, | ||
const xt::xtensor< double, 2 > & | x, | ||
const xtl::span< const std::int32_t > & | cells | ||
) |
Interpolate an expression in a finite element space.
[out] | u | The function to interpolate into |
[in] | f | The expression to be interpolated |
[in] | x | The points at which f should be evaluated, as computed by fem::interpolation_coords. The element used in fem::interpolation_coords should be the same element as associated with u. |
[in] | cells | Indices of the cells in the mesh on which to interpolate. Should be the same as the list used when calling fem::interpolation_coords. |
void dolfinx::fem::interpolate_c | ( | Function< T > & | u, |
const std::function< void(xt::xarray< T > &, const xt::xtensor< double, 2 > &)> & | f, | ||
const xt::xtensor< double, 2 > & | x, | ||
const xtl::span< const std::int32_t > & | cells | ||
) |
Interpolate an expression f(x)
[out] | u | The function to interpolate into |
[in] | f | The expression to be interpolated |
[in] | x | The points at which should be evaluated, as computed by fem::interpolation_coords |
[in] | cells | Indices of the cells in the mesh on which to interpolate. Should be the same as the list used when calling fem::interpolation_coords. |
xt::xtensor< double, 2 > dolfinx::fem::interpolation_coords | ( | const fem::FiniteElement & | element, |
const mesh::Mesh & | mesh, | ||
const xtl::span< const std::int32_t > & | cells | ||
) |
Compute the evaluation points in the physical space at which an expression should be computed to interpolate it in a finite elemenet space.
[in] | element | The element to be interpolated into |
[in] | mesh | The domain |
[in] | cells | Indices of the cells in the mesh to compute interpolation coordinates for |
std::vector< std::int32_t > dolfinx::fem::locate_dofs_geometrical | ( | const fem::FunctionSpace & | V, |
const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> & | marker_fn | ||
) |
Finds degrees of freedom whose geometric coordinate is true for the provided marking function.
[in] | V | The function (sub)space on which degrees of freedom will be located. |
[in] | marker_fn | Function marking tabulated degrees of freedom |
std::array< std::vector< std::int32_t >, 2 > dolfinx::fem::locate_dofs_geometrical | ( | const std::array< std::reference_wrapper< const fem::FunctionSpace >, 2 > & | V, |
const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> & | marker_fn | ||
) |
Finds degrees of freedom whose geometric coordinate is true for the provided marking function.
[in] | V | The function (sub)space(s) on which degrees of freedom will be located. The spaces must share the same mesh and element type. |
[in] | marker_fn | Function marking tabulated degrees of freedom |
std::vector< std::int32_t > dolfinx::fem::locate_dofs_topological | ( | const fem::FunctionSpace & | V, |
const int | dim, | ||
const xtl::span< const std::int32_t > & | entities, | ||
bool | remote = true |
||
) |
Find degrees-of-freedom which belong to the provided mesh entities (topological). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex.
[in] | V | The function (sub)space on which degrees-of-freedom (DOFs) will be located. |
[in] | dim | Topological dimension of mesh entities on which degrees-of-freedom will be located |
[in] | entities | Indices of mesh entities. All DOFs associated with the closure of these indices will be returned |
[in] | remote | True to return also "remotely located" degree-of-freedom indices. Remotely located degree-of-freedom indices are local/owned by the current process, but which the current process cannot identify because it does not recognize mesh entity as a marked. For example, a boundary condition dof at a vertex where this process does not have the associated boundary facet. This commonly occurs with partitioned meshes. |
std::array< std::vector< std::int32_t >, 2 > dolfinx::fem::locate_dofs_topological | ( | const std::array< std::reference_wrapper< const fem::FunctionSpace >, 2 > & | V, |
const int | dim, | ||
const xtl::span< const std::int32_t > & | entities, | ||
bool | remote = true |
||
) |
Find degrees-of-freedom which belong to the provided mesh entities (topological). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex.
[in] | V | The function (sub)spaces on which degrees-of-freedom (DOFs) will be located. The spaces must share the same mesh and element type. |
[in] | dim | Topological dimension of mesh entities on which degrees-of-freedom will be located |
[in] | entities | Indices of mesh entities. All DOFs associated with the closure of these indices will be returned |
[in] | remote | True to return also "remotely located" degree-of-freedom indices. Remotely located degree-of-freedom indices are local/owned by the current process, but which the current process cannot identify because it does not recognize mesh entity as a marked. For example, a boundary condition dof at a vertex where this process does not have the associated boundary facet. This commonly occurs with partitioned meshes. |
void dolfinx::fem::set_diagonal | ( | const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> & | set_fn, |
const fem::FunctionSpace & | V, | ||
const std::vector< std::shared_ptr< const DirichletBC< T >>> & | bcs, | ||
T | diagonal = 1.0 |
||
) |
Sets a value to the diagonal of the matrix for rows with a Dirichlet boundary conditions applied. This function is typically called after assembly. The assembly function zeroes Dirichlet rows and columns. This function adds the value only to rows that are locally owned, and therefore does not create a need for parallel communication. For block matrices, this function should normally be called only on the diagonal blocks, i.e. blocks for which the test and trial spaces are the same.
[in] | set_fn | The function for setting values to a matrix |
[in] | V | The function space for the rows and columns of the matrix. It is used to extract only the Dirichlet boundary conditions that are define on V or subspaces of V. |
[in] | bcs | The Dirichlet boundary condtions |
[in] | diagonal | The value to add to the diagonal for rows with a boundary condition applied |
void dolfinx::fem::set_diagonal | ( | const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const T *)> & | set_fn, |
const xtl::span< const std::int32_t > & | rows, | ||
T | diagonal = 1.0 |
||
) |
Sets a value to the diagonal of a matrix for specified rows. It is typically called after assembly. The assembly function zeroes Dirichlet rows and columns. For block matrices, this function should normally be called only on the diagonal blocks, i.e. blocks for which the test and trial spaces are the same.
[in] | set_fn | The function for setting values to a matrix |
[in] | rows | The row blocks, in local indices, for which to add a value to the diagonal |
[in] | diagonal | The value to add to the diagonal for the specified rows |
graph::AdjacencyList< std::int32_t > dolfinx::fem::transpose_dofmap | ( | const graph::AdjacencyList< std::int32_t > & | dofmap, |
std::int32_t | num_cells | ||
) |
Create an adjacency list that maps a global index (process-wise) to the 'unassembled' cell-wise contributions. It is built from the usual (cell, local index) -> global index dof map. An 'unassembled' vector is the stacked cell contributions, ordered by cell index.
If the usual dof map is:
Cell: 0 1 2 3 Global index: [ [0, 3, 5], [3, 2, 4], [4, 3, 2], [2, 1, 0]]
the 'transpose' dof map will be:
Global index: 0 1 2 3 4 5 Unassembled index: [ [0, 11], [10], [4, 8, 9], [1, 3, 7], [5, 6], [2] ]
[in] | dofmap | The standard dof map that for each cell (node) gives the global (process-wise) index of each local (cell-wise) index. |
[in] | num_cells | The number of cells (nodes) in dofmap to consider. The first num_cells are used. This is argument is typically used to exclude ghost cell contributions. |