| DOLFINx 0.9.0
    DOLFINx C++ interface | 
Linear algebra interface. More...
| Namespaces | |
| namespace | petsc | 
| PETSc linear algebra functions. | |
| Classes | |
| class | MatrixCSR | 
| Distributed sparse matrix.  More... | |
| class | SLEPcEigenSolver | 
| This class provides an eigenvalue solver for PETSc matrices. It is a wrapper for the SLEPc eigenvalue solver.  More... | |
| class | SparsityPattern | 
| class | Vector | 
| Concepts | |
| concept | MatSet | 
| Matrix accumulate/set concept for functions that can be used in assemblers to accumulate or set values in a matrix. | |
| Enumerations | |
| enum class | BlockMode : int { compact = 0 , expanded = 1 } | 
| Modes for representing block structured matrices.  More... | |
| enum class | Norm { l1 , l2 , linf , frobenius } | 
| Norm types. | |
| Functions | |
| template<class V > | |
| auto | inner_product (const V &a, const V &b) | 
| template<class V > | |
| auto | squared_norm (const V &a) | 
| template<class V > | |
| auto | norm (const V &x, Norm type=Norm::l2) | 
| template<class V > | |
| void | orthonormalize (std::vector< std::reference_wrapper< V > > basis) | 
| template<class V > | |
| bool | is_orthonormal (std::vector< std::reference_wrapper< const V > > basis, dolfinx::scalar_value_type_t< typename V::value_type > eps=std::numeric_limits< dolfinx::scalar_value_type_t< typename V::value_type > >::epsilon()) | 
| Test if basis is orthonormal. | |
Linear algebra interface.
Interface to linear algebra data structures and solvers.
| 
 | strong | 
| auto inner_product | ( | const V & | a, | 
| const V & | b ) | 
Compute the inner product of two vectors. The two vectors must have the same parallel layout
| a | A vector | 
| b | A vector | 
a^{H} b (a^{T} b if a and b are real) | bool is_orthonormal | ( | std::vector< std::reference_wrapper< const V > > | basis, | 
| dolfinx::scalar_value_type_t< typename V::value_type > | eps = std::numeric_limits< dolfinx::scalar_value_type_t<typename V::value_type>>::epsilon() ) | 
Test if basis is orthonormal.
Returns true if ||x_i - x_j|| - delta_{ij} < eps fro all i, j, and otherwise false.
| [in] | basis | Set of vectors to check. | 
| [in] | eps | Tolerance. | 
| auto norm | ( | const V & | x, | 
| Norm | type = Norm::l2 ) | 
| void orthonormalize | ( | std::vector< std::reference_wrapper< V > > | basis | ) | 
Orthonormalize a set of vectors
| [in,out] | basis | The set of vectors to orthonormalise. The vectors must have identical parallel layouts. The vectors are modified in-place. | 
| V | dolfinx::la::Vector | 
| auto squared_norm | ( | const V & | a | ) | 
Compute the squared L2 norm of vector