Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d6/dfa/namespacedolfinx_1_1la.html
DOLFINx  0.4.1
DOLFINx C++ interface
Classes | Enumerations | Functions
dolfinx::la Namespace Reference

Linear algebra interface. More...

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
 This class provides a sparsity pattern data structure that can be used to initialize sparse matrices. After assembly, column indices are always sorted in increasing order. Ghost entries are kept after assembly. More...
 
class  Vector
 Distributed vector. More...
 

Enumerations

enum class  Norm { l1 , l2 , linf , frobenius }
 Norm types.
 

Functions

template<typename T , class Allocator = std::allocator<T>>
inner_product (const Vector< T, Allocator > &a, const Vector< T, Allocator > &b)
 Compute the inner product of two vectors. The two vectors must have the same parallel layout. More...
 
template<typename T , typename U >
void orthonormalize (const xtl::span< Vector< T, U >> &basis, double tol=1.0e-10)
 Orthonormalize a set of vectors. More...
 
template<typename T , typename U >
bool is_orthonormal (const xtl::span< const Vector< T, U >> &basis, double tol=1.0e-10)
 Test if basis is orthonormal. More...
 

Detailed Description

Linear algebra interface.

Interface to linear algebra data structures and solvers

Function Documentation

◆ inner_product()

T dolfinx::la::inner_product ( const Vector< T, Allocator > &  a,
const Vector< T, Allocator > &  b 
)

Compute the inner product of two vectors. The two vectors must have the same parallel layout.

Note
Collective
Parameters
aA vector
bA vector
Returns
Returns a^{H} b (a^{T} b if a and b are real)

◆ is_orthonormal()

bool dolfinx::la::is_orthonormal ( const xtl::span< const Vector< T, U >> &  basis,
double  tol = 1.0e-10 
)

Test if basis is orthonormal.

Parameters
[in]basisThe set of vectors to check
[in]tolThe tolerance used to test for orthonormality
Returns
True is basis is orthonormal, otherwise false

◆ orthonormalize()

void dolfinx::la::orthonormalize ( const xtl::span< Vector< T, U >> &  basis,
double  tol = 1.0e-10 
)

Orthonormalize a set of vectors.

Parameters
[in,out]basisThe set of vectors to orthonormalise. The vectors must have identical parallel layouts. The vectors are modified in-place.
[in]tolThe tolerance used to detect a linear dependency