Math namespace
Module for general mathematical operations.
Classes
-
template<class Operator, class Vector>class LinearSystem
- Represents a linear system of equations.
-
template<class MatrixScalar, class VectorScalar>class LinearSystem<Math::Matrix<MatrixScalar>, Math::Vector<VectorScalar>>
- Linear system specialization for dense matrices.
-
template<class MatrixScalar, class VectorScalar>class LinearSystem<Math::SparseMatrix<MatrixScalar>, Math::Vector<VectorScalar>>
- Linear system specialization for sparse matrices.
- class LinearSystemBase
- Base class for linear systems of the form .
- class Rad
- Represents an angle in radians.
-
template<class Derived, class T>class Unit
- Base class for units using CRTP.
Typedefs
-
template<class ScalarType>using Matrix = Eigen::Matrix<ScalarType, Eigen::Dynamic, Eigen::Dynamic>
- Dynamic-size dense matrix type.
- using RealMatrix = Matrix<Real>
- Real-valued dense matrix.
- using ComplexMatrix = Matrix<Complex>
- Complex-valued dense matrix.
-
template<class ScalarType>using SpatialMatrix = Eigen::Matrix<ScalarType, Eigen::Dynamic, Eigen::Dynamic, 0, RODIN_
MAXIMAL_ SPACE_ DIMENSION, RODIN_ MAXIMAL_ SPACE_ DIMENSION> - Spatial matrix with bounded maximum dimensions.
-
using PointMatrix = Eigen::Matrix<Real, Eigen::Dynamic, Eigen::Dynamic, 0, RODIN_
MAXIMAL_ SPACE_ DIMENSION, Eigen::Dynamic> - Point matrix with bounded row dimension.
-
template<class ScalarType, size_t Rows, size_t Cols>using FixedSizeMatrix = Eigen::Matrix<ScalarType, Rows, Cols>
- Fixed-size matrix type.
-
template<class ScalarType>using SparseMatrix = Eigen::SparseMatrix<ScalarType>
- Sparse matrix type.
-
template<class ScalarType>using Vector = Eigen::VectorX<ScalarType>
- Dynamic-size dense vector type.
- using ComplexVector = Vector<Complex>
- Dynamic-size complex-valued vector.
-
template<class ScalarType>using SpatialVector = Eigen::Matrix<ScalarType, Eigen::Dynamic, 1, 0, RODIN_
MAXIMAL_ SPACE_ DIMENSION, 1> - Spatial vector with bounded maximum size.
- using SpatialPoint = SpatialVector<Real>
- Real-valued spatial vector for point coordinates.
-
template<class ScalarType, size_t Size>using FixedSizeVector = Eigen::Vector<ScalarType, Size>
- Fixed-size vector type.
-
template<class ScalarType>using Vector2 = FixedSizeVector<ScalarType, 2>
- 2D fixed-size vector.
-
template<class ScalarType>using Vector3 = FixedSizeVector<ScalarType, 3>
- 3D fixed-size vector.
-
template<class ScalarType>using Vector4 = FixedSizeVector<ScalarType, 4>
- 4D fixed-size vector.
-
template<class ScalarType>using Vector8 = FixedSizeVector<ScalarType, 8>
- 8D fixed-size vector.
-
template<class ScalarType>using Vector16 = FixedSizeVector<ScalarType, 16>
- 16D fixed-size vector.
-
template<class ScalarType>using Vector32 = FixedSizeVector<ScalarType, 32>
- 32D fixed-size vector.
-
template<class ScalarType>using Vector64 = FixedSizeVector<ScalarType, 64>
- 64D fixed-size vector.
-
template<class ScalarType>using Vector128 = FixedSizeVector<ScalarType, 128>
- 128D fixed-size vector.
Functions
-
template<class T>auto abs(const T& x) -> auto constexpr
- Computes the absolute value of a value.
-
template<class T>auto exp(const T& x) -> auto constexpr
- Computes the exponential function.
- auto conj(const Complex& x) -> Complex constexpr
- Computes the complex conjugate of a complex number.
-
template<class T>auto conj(const Eigen::MatrixBase<T>& x) -> auto constexpr
- Computes the conjugate of an Eigen matrix.
- auto conj(const Real& x) -> Real constexpr
- Identity conjugate for real numbers.
-
template<class Base, class Exponent>auto pow2(const Base& base) -> auto constexpr
- Computes the square of a value.
-
template<class Base, class Exponent>auto pow(const Base& base, const Exponent& exponent) -> auto constexpr
- Computes a power with arbitrary exponent.
-
template<class T>auto sqrt(const T& x) -> auto constexpr
- Computes the square root of a value.
-
template<class T>auto isNaN(const T& x) -> Boolean constexpr
- Determines if a floating point number is not-a-number (NaN).
- auto isNaN(const Complex& x) -> Boolean constexpr
- Determines if a complex number is not-a-number (NaN).
-
template<class T>auto isInf(const T& x) -> Boolean constexpr
- Determines if a floating point number is infinite.
-
template<class T>auto cos(const T& x) -> auto constexpr
- Computes the cosine function.
-
template<class T>auto cosh(const T& x) -> auto constexpr
- Computes the hyperbolic cosine function.
-
template<class T>auto sin(const T& x) -> auto constexpr
- Computes the sine function.
-
template<class T>auto sinh(const T& x) -> auto constexpr
- Computes the hyperbolic sine function.
-
template<class T>auto tan(const T& x) -> auto constexpr
- Computes the tangent function.
-
template<typename T>auto sgn(const T& x) -> T constexpr
- Computes the sign (signum) function.
-
template<class T>auto binom(const T& n, const T& k) -> T constexpr
- Computes the binomial coefficient.
-
template<class T>auto factorial(const T& n) -> T constexpr
- Computes the factorial function.
-
template<class T>auto permutation(const T& n, const T& k) -> T constexpr
- Computes the permutation function.
-
template<class T>auto nan() -> auto constexpr
- Returns a quiet NaN (not-a-number) value.
- auto nan() -> Complex constexpr
- Returns a complex NaN value.
-
template<class LHS, class RHS>auto sum(const LHS& lhs, const RHS& rhs) -> auto constexpr
- Computes the sum of two values.
-
template<class Operand>auto minus(const Operand& op) -> auto constexpr
- Computes the unary minus (negation).
-
template<class LHS, class RHS>auto minus(const LHS& lhs, const RHS& rhs) -> auto constexpr
- Computes the difference of two values.
-
template<class LHS, class RHS>auto mult(const LHS& lhs, const RHS& rhs) -> auto constexpr
- Computes the product of two values.
-
template<class LHS, class RHS>auto division(const LHS& lhs, const RHS& rhs) -> auto constexpr
- Computes the division of two values.
- auto dot(const Real& lhs, const Real& rhs) -> Real constexpr
- Computes the dot product of two real numbers.
- auto dot(const Complex& lhs, const Complex& rhs) -> Complex constexpr
- Computes the dot product of two complex numbers.
-
template<class LHSDerived, class RHSDerived>auto dot(const Eigen::MatrixBase<LHSDerived>& lhs, const Eigen::MatrixBase<RHSDerived>& rhs) -> auto constexpr
- Computes the dot product of two Eigen vectors/matrices.
-
template<class AScalarType, class YScalarType, class XScalarType>void axpy(SparseMatrix<YScalarType>& y, AScalarType alpha, const SparseMatrix<XScalarType>& x)
- Performs the AXPY operation on sparse matrices.
Typedef documentation
#include <Rodin/Math/Matrix.h>
template<class ScalarType>
using Rodin:: Math:: Matrix = Eigen::Matrix<ScalarType, Eigen::Dynamic, Eigen::Dynamic>
Dynamic-size dense matrix type.
| Template parameters | |
|---|---|
| ScalarType | The element type (e.g., Real, Complex) |
A dense matrix with both dimensions determined at runtime. The matrix is stored in column-major order (Fortran style) and supports standard linear algebra operations.
using Rodin:: Math:: RealMatrix = Matrix<Real>
#include <Rodin/Math/Matrix.h>
Real-valued dense matrix.
Convenience alias for Matrix<Real>.
using Rodin:: Math:: ComplexMatrix = Matrix<Complex>
#include <Rodin/Math/Matrix.h>
Complex-valued dense matrix.
Convenience alias for Matrix<Complex>, commonly used in eigenvalue problems, frequency-domain analysis, and other applications requiring complex arithmetic.
#include <Rodin/Math/Matrix.h>
template<class ScalarType>
using Rodin:: Math:: SpatialMatrix = Eigen::Matrix<ScalarType, Eigen::Dynamic, Eigen::Dynamic, 0, RODIN_ MAXIMAL_ SPACE_ DIMENSION, RODIN_ MAXIMAL_ SPACE_ DIMENSION>
Spatial matrix with bounded maximum dimensions.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A dynamic-size matrix with maximum dimensions bounded by RODIN_MAXIMAL_SPACE_DIMENSION. Used for geometric transformations, Jacobians, and other spatial operators to optimize memory allocation.
Typical uses:
- Jacobian matrices:
- Metric tensors
- Local coordinate transformations
using Rodin:: Math:: PointMatrix = Eigen::Matrix<Real, Eigen::Dynamic, Eigen::Dynamic, 0, RODIN_ MAXIMAL_ SPACE_ DIMENSION, Eigen::Dynamic>
#include <Rodin/Math/Matrix.h>
Point matrix with bounded row dimension.
A dynamic-size matrix where the number of rows is bounded by RODIN_MAXIMAL_SPACE_DIMENSION. Commonly used to store collections of spatial points as columns.
Example: storing vertex coordinates of a mesh element where each column represents one vertex's coordinates.
#include <Rodin/Math/Matrix.h>
template<class ScalarType, size_t Rows, size_t Cols>
using Rodin:: Math:: FixedSizeMatrix = Eigen::Matrix<ScalarType, Rows, Cols>
Fixed-size matrix type.
| Template parameters | |
|---|---|
| ScalarType | The element type |
| Rows | The number of rows (must be known at compile time) |
| Cols | The number of columns (must be known at compile time) |
A compile-time fixed-size matrix. Both dimensions are known at compile time, enabling optimizations and stack allocation.
Example uses:
- 2×2 rotation matrices
- 3×3 transformation matrices
- Small dense blocks in larger systems
#include <Rodin/Math/SparseMatrix.h>
template<class ScalarType>
using Rodin:: Math:: SparseMatrix = Eigen::SparseMatrix<ScalarType>
Sparse matrix type.
| Template parameters | |
|---|---|
| ScalarType | The element type (e.g., Real, Complex) |
A sparse matrix stores only the non-zero elements, making it memory-efficient for large matrices arising from finite element discretizations where most entries are zero. The matrix uses compressed column storage (CCS) format by default.
Storage Format
Eigen's SparseMatrix uses Compressed Column Storage (CCS), also known as Compressed Sparse Column (CSC) format, which stores:
- Non-zero values
- Row indices for each non-zero value
- Column pointers indicating where each column starts
Typical Applications
- Stiffness matrices: in
- Mass matrices: in
- Discretized differential operators
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector = Eigen::VectorX<ScalarType>
Dynamic-size dense vector type.
| Template parameters | |
|---|---|
| ScalarType | The element type (e.g., Real, Complex) |
A column vector with dynamic size determined at runtime. The vector supports standard linear algebra operations.
using Rodin:: Math:: ComplexVector = Vector<Complex>
#include <Rodin/Math/Vector.h>
Dynamic-size complex-valued vector.
Convenience alias for Vector<Complex>.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: SpatialVector = Eigen::Matrix<ScalarType, Eigen::Dynamic, 1, 0, RODIN_ MAXIMAL_ SPACE_ DIMENSION, 1>
Spatial vector with bounded maximum size.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A dynamic-size vector with maximum size bounded by RODIN_MAXIMAL_SPACE_DIMENSION. Used for geometric quantities in 2D or 3D space to optimize memory allocation.
using Rodin:: Math:: SpatialPoint = SpatialVector<Real>
#include <Rodin/Math/Vector.h>
Real-valued spatial vector for point coordinates.
Convenience alias for SpatialVector<Real>, commonly used to represent points in 2D or 3D space.
#include <Rodin/Math/Vector.h>
template<class ScalarType, size_t Size>
using Rodin:: Math:: FixedSizeVector = Eigen::Vector<ScalarType, Size>
Fixed-size vector type.
| Template parameters | |
|---|---|
| ScalarType | The element type |
| Size | The number of elements (must be known at compile time) |
A compile-time fixed-size column vector. The size is known at compile time, enabling optimizations and stack allocation.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector2 = FixedSizeVector<ScalarType, 2>
2D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 2 elements.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector3 = FixedSizeVector<ScalarType, 3>
3D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 3 elements.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector4 = FixedSizeVector<ScalarType, 4>
4D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 4 elements.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector8 = FixedSizeVector<ScalarType, 8>
8D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 8 elements.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector16 = FixedSizeVector<ScalarType, 16>
16D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 16 elements.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector32 = FixedSizeVector<ScalarType, 32>
32D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 32 elements.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector64 = FixedSizeVector<ScalarType, 64>
64D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 64 elements.
#include <Rodin/Math/Vector.h>
template<class ScalarType>
using Rodin:: Math:: Vector128 = FixedSizeVector<ScalarType, 128>
128D fixed-size vector.
| Template parameters | |
|---|---|
| ScalarType | The element type |
A vector with exactly 128 elements.
Function documentation
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: abs(const T& x) constexpr
Computes the absolute value of a value.
| Template parameters | |
|---|---|
| T | Type of value (Real, Complex, etc.) |
| Parameters | |
| x in | Value |
| Returns | Absolute value |
Returns , the absolute value (magnitude) of .
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: exp(const T& x) constexpr
Computes the exponential function.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Exponent |
| Returns | |
Returns where is Euler's number.
Complex Rodin:: Math:: conj(const Complex& x) constexpr
#include <Rodin/Math/Common.h>
Computes the complex conjugate of a complex number.
| Parameters | |
|---|---|
| x in | Complex number |
| Returns | Complex conjugate |
For , returns .
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: conj(const Eigen::MatrixBase<T>& x) constexpr
Computes the conjugate of an Eigen matrix.
| Template parameters | |
|---|---|
| T | Eigen matrix type |
| Parameters | |
| x in | Matrix |
| Returns | Conjugate matrix |
For matrices with complex entries, returns the element-wise conjugate.
Real Rodin:: Math:: conj(const Real& x) constexpr
#include <Rodin/Math/Common.h>
Identity conjugate for real numbers.
| Parameters | |
|---|---|
| x in | Real number |
| Returns | |
For real numbers, the conjugate is the number itself.
#include <Rodin/Math/Common.h>
template<class Base, class Exponent>
auto Rodin:: Math:: pow2(const Base& base) constexpr
Computes the square of a value.
| Template parameters | |
|---|---|
| Base | Type of the value |
| Exponent | Unused template parameter for compatibility |
| Parameters | |
| base in | Value to square |
| Returns | |
Returns efficiently without calling std::pow.
#include <Rodin/Math/Common.h>
template<class Base, class Exponent>
auto Rodin:: Math:: pow(const Base& base,
const Exponent& exponent) constexpr
Computes a power with arbitrary exponent.
| Template parameters | |
|---|---|
| Base | Type of the base |
| Exponent | Type of the exponent |
| Parameters | |
| base in | Base value |
| exponent in | Exponent value |
| Returns | |
Returns where is the base and is the exponent.
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: sqrt(const T& x) constexpr
Computes the square root of a value.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Value (must be non-negative for real types) |
| Returns | |
Returns , the principal square root of .
#include <Rodin/Math/Common.h>
template<class T>
Boolean Rodin:: Math:: isNaN(const T& x) constexpr
Determines if a floating point number is not-a-number (NaN).
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Value to check |
| Returns | True if is NaN, false otherwise |
Checks if is NaN (not a number), which typically results from undefined operations like or for reals.
Boolean Rodin:: Math:: isNaN(const Complex& x) constexpr
#include <Rodin/Math/Common.h>
Determines if a complex number is not-a-number (NaN).
| Parameters | |
|---|---|
| x in | Complex number to check |
| Returns | True if real or imaginary part is NaN, false otherwise |
Returns true if either the real or imaginary part is NaN.
#include <Rodin/Math/Common.h>
template<class T>
Boolean Rodin:: Math:: isInf(const T& x) constexpr
Determines if a floating point number is infinite.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Value to check |
| Returns | True if is or , false otherwise |
Checks if , which can result from overflow or operations like .
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: cos(const T& x) constexpr
Computes the cosine function.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Angle in radians |
| Returns | |
Returns where is in radians.
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: cosh(const T& x) constexpr
Computes the hyperbolic cosine function.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Value |
| Returns | |
Returns .
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: sin(const T& x) constexpr
Computes the sine function.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Angle in radians |
| Returns | |
Returns where is in radians.
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: sinh(const T& x) constexpr
Computes the hyperbolic sine function.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Value |
| Returns | |
Returns .
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: tan(const T& x) constexpr
Computes the tangent function.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Angle in radians |
| Returns | |
Returns where is in radians.
#include <Rodin/Math/Common.h>
template<typename T>
T Rodin:: Math:: sgn(const T& x) constexpr
Computes the sign (signum) function.
| Template parameters | |
|---|---|
| T | Type of value |
| Parameters | |
| x in | Value |
| Returns | Sign of : -1, 0, or +1 |
Returns:
#include <Rodin/Math/Common.h>
template<class T>
T Rodin:: Math:: binom(const T& n,
const T& k) constexpr
Computes the binomial coefficient.
| Template parameters | |
|---|---|
| T | Type of value (typically Integer) |
| Parameters | |
| n in | Total number of items |
| k in | Number of items to choose |
| Returns | |
Returns , the number of ways to choose items from items without regard to order.
#include <Rodin/Math/Common.h>
template<class T>
T Rodin:: Math:: factorial(const T& n) constexpr
Computes the factorial function.
| Template parameters | |
|---|---|
| T | Type of value (typically Integer) |
| Parameters | |
| n in | Non-negative integer |
| Returns | |
Returns .
#include <Rodin/Math/Common.h>
template<class T>
T Rodin:: Math:: permutation(const T& n,
const T& k) constexpr
Computes the permutation function.
| Template parameters | |
|---|---|
| T | Type of value (typically Integer) |
| Parameters | |
| n in | Total number of items |
| k in | Number of items to arrange |
| Returns | |
Returns , the number of ways to arrange items from items where order matters.
#include <Rodin/Math/Common.h>
template<class T>
auto Rodin:: Math:: nan() constexpr
Returns a quiet NaN (not-a-number) value.
| Template parameters | |
|---|---|
| T | Type of value (must be a floating point type) |
| Returns | Quiet NaN value of type T |
Returns a NaN value of the specified type, useful for indicating undefined or invalid results.
Complex Rodin:: Math:: nan() constexpr
#include <Rodin/Math/Common.h>
Returns a complex NaN value.
| Returns | Complex NaN value |
|---|
Returns a complex number where both real and imaginary parts are NaN.
#include <Rodin/Math/Common.h>
template<class LHS, class RHS>
auto Rodin:: Math:: sum(const LHS& lhs,
const RHS& rhs) constexpr
Computes the sum of two values.
| Template parameters | |
|---|---|
| LHS | Type of left-hand side |
| RHS | Type of right-hand side |
| Parameters | |
| lhs in | Left-hand side operand |
| rhs in | Right-hand side operand |
| Returns | |
Returns . This function is used internally for type deduction in the form language.
#include <Rodin/Math/Common.h>
template<class Operand>
auto Rodin:: Math:: minus(const Operand& op) constexpr
Computes the unary minus (negation).
| Template parameters | |
|---|---|
| Operand | Type of operand |
| Parameters | |
| op in | Operand to negate |
| Returns | |
Returns . This function is used internally for type deduction in the form language.
#include <Rodin/Math/Common.h>
template<class LHS, class RHS>
auto Rodin:: Math:: minus(const LHS& lhs,
const RHS& rhs) constexpr
Computes the difference of two values.
| Template parameters | |
|---|---|
| LHS | Type of left-hand side |
| RHS | Type of right-hand side |
| Parameters | |
| lhs in | Left-hand side operand |
| rhs in | Right-hand side operand |
| Returns | |
Returns . This function is used internally for type deduction in the form language.
#include <Rodin/Math/Common.h>
template<class LHS, class RHS>
auto Rodin:: Math:: mult(const LHS& lhs,
const RHS& rhs) constexpr
Computes the product of two values.
| Template parameters | |
|---|---|
| LHS | Type of left-hand side |
| RHS | Type of right-hand side |
| Parameters | |
| lhs in | Left-hand side operand |
| rhs in | Right-hand side operand |
| Returns | |
Returns . This function is used internally for type deduction in the form language.
#include <Rodin/Math/Common.h>
template<class LHS, class RHS>
auto Rodin:: Math:: division(const LHS& lhs,
const RHS& rhs) constexpr
Computes the division of two values.
| Template parameters | |
|---|---|
| LHS | Type of left-hand side |
| RHS | Type of right-hand side |
| Parameters | |
| lhs in | Left-hand side operand (numerator) |
| rhs in | Right-hand side operand (denominator) |
| Returns | |
Returns . This function is used internally for type deduction in the form language.
Real Rodin:: Math:: dot(const Real& lhs,
const Real& rhs) constexpr
#include <Rodin/Math/Common.h>
Computes the dot product of two real numbers.
| Parameters | |
|---|---|
| lhs in | Left-hand side operand |
| rhs in | Right-hand side operand |
| Returns | |
For real numbers, the dot product is simply multiplication: .
Complex Rodin:: Math:: dot(const Complex& lhs,
const Complex& rhs) constexpr
#include <Rodin/Math/Common.h>
Computes the dot product of two complex numbers.
| Parameters | |
|---|---|
| lhs in | Left-hand side complex number |
| rhs in | Right-hand side complex number |
| Returns | |
For complex numbers, the dot product uses conjugation: .
#include <Rodin/Math/Common.h>
template<class LHSDerived, class RHSDerived>
auto Rodin:: Math:: dot(const Eigen::MatrixBase<LHSDerived>& lhs,
const Eigen::MatrixBase<RHSDerived>& rhs) constexpr
Computes the dot product of two Eigen vectors/matrices.
| Template parameters | |
|---|---|
| LHSDerived | Eigen type of left-hand side |
| RHSDerived | Eigen type of right-hand side |
| Parameters | |
| lhs in | Left-hand side vector/matrix |
| rhs in | Right-hand side vector/matrix |
| Returns | Dot product (scalar) |
For vectors and , computes:
where is the complex conjugate (identity for real values).
#include <Rodin/Math/SparseMatrix.h>
template<class AScalarType, class YScalarType, class XScalarType>
void Rodin:: Math:: axpy(SparseMatrix<YScalarType>& y,
AScalarType alpha,
const SparseMatrix<XScalarType>& x)
Performs the AXPY operation on sparse matrices.
| Template parameters | |
|---|---|
| AScalarType | Scalar type for alpha |
| YScalarType | Scalar type for matrix y |
| XScalarType | Scalar type for matrix x |
| Parameters | |
| y in/out | Matrix to be updated (accumulated result) |
| alpha in | Scalar multiplier for x |
| x in | Matrix to be scaled and added to y |
Computes where and are sparse matrices and is a scalar.
This is the sparse matrix version of the Level 1 BLAS AXPY operation.