|
DOLFINx 0.10.0
DOLFINx C++ interface
|
Functions supporting assembly of finite element fem::Form and fem::Expression. More...
#include "Function.h"#include "FunctionSpace.h"#include "assemble_expression_impl.h"#include "assemble_matrix_impl.h"#include "assemble_scalar_impl.h"#include "assemble_vector_impl.h"#include "pack.h"#include "traits.h"#include "utils.h"#include <algorithm>#include <basix/mdspan.hpp>#include <cstdint>#include <dolfinx/common/types.h>#include <memory>#include <optional>#include <span>#include <vector>

Go to the source code of this file.
Namespaces | |
| namespace | dolfinx |
| Top-level namespace. | |
| namespace | dolfinx::fem |
| Finite element method functionality. | |
Functions | |
| template<dolfinx::scalar T, std::floating_point U> | |
| void | tabulate_expression (std::span< T > values, const fem::Expression< T, U > &e, md::mdspan< const T, md::dextents< std::size_t, 2 > > coeffs, std::span< const T > constants, const mesh::Mesh< U > &mesh, fem::MDSpan2 auto entities, std::optional< std::pair< std::reference_wrapper< const FiniteElement< U > >, std::size_t > > element) |
| Evaluate an Expression on cells or facets. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| void | tabulate_expression (std::span< T > values, const fem::Expression< T, U > &e, const mesh::Mesh< U > &mesh, fem::MDSpan2 auto entities) |
| Evaluate an Expression on cells or facets. | |
| template<dolfinx::scalar T> | |
| std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > | make_coefficients_span (const std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > &coeffs) |
| Create a map of std::spans from a map of std::vectors. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| T | assemble_scalar (const Form< T, U > &M, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients) |
| Assemble functional into scalar. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| T | assemble_scalar (const Form< T, U > &M) |
| Assemble functional into scalar. | |
| template<typename V, std::floating_point U, dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type> requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T> | |
| void | assemble_vector (V &&b, const Form< T, U > &L, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients) |
| Assemble linear form into a vector. | |
| template<typename V, std::floating_point U, dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type> requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T> | |
| void | assemble_vector (V &&b, const Form< T, U > &L) |
| Assemble linear form into a vector. | |
| template<typename V, std::floating_point U = scalar_value_t<typename std::remove_cvref_t<V>::value_type>, dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type> requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T> | |
| void | apply_lifting (V &&b, std::vector< std::optional< std::reference_wrapper< const Form< T, U > > > > a, const std::vector< std::span< const T > > &constants, const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > > &coeffs, const std::vector< std::vector< std::reference_wrapper< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T alpha) |
| Modify the right-hand side vector to account for constraints (Dirichlet boundary condition constraints). This modification is known as 'lifting'. | |
| template<typename V, std::floating_point U = scalar_value_t<typename std::remove_cvref_t<V>::value_type>, dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type> requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T> | |
| void | apply_lifting (V &&b, std::vector< std::optional< std::reference_wrapper< const Form< T, U > > > > a, const std::vector< std::vector< std::reference_wrapper< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T alpha) |
| Modify the right-hand side vector to account for constraints (Dirichlet boundary conditions constraints). This modification is known as 'lifting'. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| void | assemble_matrix (la::MatSet< T > auto mat_add, const Form< T, U > &a, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients, std::span< const std::int8_t > dof_marker0, std::span< const std::int8_t > dof_marker1) |
| Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| void | assemble_matrix (auto mat_add, const Form< T, U > &a, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients, const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > &bcs) |
| Assemble bilinear form into a matrix. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| void | assemble_matrix (auto mat_add, const Form< T, U > &a, const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > &bcs) |
| Assemble bilinear form into a matrix. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| void | assemble_matrix (auto mat_add, const Form< T, U > &a, std::span< const std::int8_t > dof_marker0, std::span< const std::int8_t > dof_marker1) |
| Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix. | |
| template<dolfinx::scalar T> | |
| void | set_diagonal (auto set_fn, std::span< const std::int32_t > rows, T diagonal=1.0) |
| Sets a value to the diagonal of a matrix for specified rows. | |
| template<dolfinx::scalar T, std::floating_point U> | |
| void | set_diagonal (auto set_fn, const FunctionSpace< U > &V, const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > &bcs, T diagonal=1.0) |
| Sets a value to the diagonal of the matrix for rows with a Dirichlet boundary conditions applied. | |
Functions supporting assembly of finite element fem::Form and fem::Expression.