CVM Class Library  8.1
This C++ class library encapsulates concepts of vector and different matrices including square, band, symmetric and hermitian ones in Euclidean space of real and complex numbers.
 All Classes Files Functions Variables Typedefs Friends Macros Pages
Public Member Functions | List of all members
basic_srmatrix< TR > Class Template Reference

End-user class encapsulating square matrix of real numbers. More...

#include <cvm.h>

Inheritance diagram for basic_srmatrix< TR >:
Inheritance graph
[legend]
Collaboration diagram for basic_srmatrix< TR >:
Collaboration graph
[legend]

Public Member Functions

 basic_srmatrix ()
 Default constructor.
 basic_srmatrix (tint nDim)
 Constructor.
 basic_srmatrix (TR *pd, tint nDim)
 Constructor.
 basic_srmatrix (const TR *pd, tint nDim)
 Constructor.
 basic_srmatrix (const basic_srmatrix &m)
 Copy constructor.
 basic_srmatrix (basic_srmatrix &&m) noexcept
 Move constructor.
 basic_srmatrix (const BaseRMatrix &m)
 Constructor.
 basic_srmatrix (const RVector &v)
 Constructor.
 basic_srmatrix (BaseRMatrix &m, tint nRow, tint nCol, tint nDim)
 Submatrix constructor.
type_proxy< TR, TR > operator() (tint nRow, tint nCol) throw (cvmexception)
 Reference to element (l-value)
TR operator() (tint nRow, tint nCol) const throw (cvmexception)
 Value of element (not l-value)
RVector operator() (tint nCol) throw (cvmexception)
 Column as l-value.
const RVector operator() (tint nCol) const throw (cvmexception)
 Column as not l-value.
RVector operator[] (tint nRow) throw (cvmexception)
 Row as l-value.
const RVector operator[] (tint nRow) const throw (cvmexception)
 Row as not l-value.
basic_srmatrixoperator= (const basic_srmatrix &m) throw (cvmexception)
 Assignment operator.
basic_srmatrixoperator= (basic_srmatrix &&m) throw (cvmexception)
 Move assignment operator.
basic_srmatrixassign (const RVector &v) throw (cvmexception)
 Vector (as array) assignment.
basic_srmatrixassign (const TR *pd)
 External array assignment.
basic_srmatrixassign (tint nRow, tint nCol, const BaseRMatrix &m) throw (cvmexception)
 Assignment to submatrix.
basic_srmatrixset (TR d)
 Sets all elements to one value.
basic_srmatrixresize (tint nNewDim) throw (cvmexception)
 Changes dimension.
basic_srmatrixoperator<< (const basic_srmatrix &m) throw (cvmexception)
 Matrix replacement.
basic_srmatrix operator+ (const basic_srmatrix &m) const throw (cvmexception)
 Addition operator.
basic_srmatrix operator- (const basic_srmatrix &m) const throw (cvmexception)
 Subtraction operator.
basic_srmatrixsum (const basic_srmatrix &m1, const basic_srmatrix &m2) throw (cvmexception)
 Sum of matrices.
basic_srmatrixdiff (const basic_srmatrix &m1, const basic_srmatrix &m2) throw (cvmexception)
 Difference of matrices.
basic_srmatrixoperator+= (const basic_srmatrix &m) throw (cvmexception)
 Increment operator.
basic_srmatrixoperator-= (const basic_srmatrix &m) throw (cvmexception)
 Decrement operator.
basic_srmatrix operator- () const
 Unary minus operator.
basic_srmatrixoperator++ ()
 Plus identity, prefix.
basic_srmatrix operator++ (int)
 Plus identity, postfix.
basic_srmatrixoperator-- ()
 Minus identity, prefix.
basic_srmatrix operator-- (int)
 Minus identity, postfix.
basic_srmatrix operator* (TR dMult) const
 Multiply by number operator.
basic_srmatrix operator/ (TR dDiv) const throw (cvmexception)
 Divide by number operator.
basic_srmatrixoperator*= (TR dMult)
 Multiply by number and assign.
basic_srmatrixoperator/= (TR dDiv)
 Divide by number and assign.
basic_srmatrixnormalize ()
 Matrix normalizer.
basic_srmatrix operator~ () const throw (cvmexception)
 Matrix transposition.
basic_srmatrixtranspose (const basic_srmatrix &m) throw (cvmexception)
 Matrix transposition.
basic_srmatrixtranspose ()
 Matrix transposition (in-place)
RVector operator* (const RVector &v) const throw (cvmexception)
 Matrix-vector product.
BaseRMatrix operator* (const BaseRMatrix &m) const throw (cvmexception)
 Matrix-matrix product.
basic_srmatrix operator* (const basic_srmatrix &m) const throw (cvmexception)
 Matrix-matrix product.
basic_srmatrixoperator*= (const basic_srmatrix &m) throw (cvmexception)
 Matrix-matrix product with assignment.
basic_srmatrixswap_rows (tint n1, tint n2) throw (cvmexception)
 Rows swap.
basic_srmatrixswap_cols (tint n1, tint n2) throw (cvmexception)
 Columns swap.
RVector solve (const RVector &vB, TR &dErr) const throw (cvmexception)
 Linear solver.
RVector solve_tran (const RVector &vB, TR &dErr) const throw (cvmexception)
 Linear solver (transposed)
RVector solve (const RVector &vB) const throw (cvmexception)
 Linear solver.
RVector solve_tran (const RVector &vB) const throw (cvmexception)
 Linear solver (transposed)
BaseRMatrix solve (const BaseRMatrix &mB, TR &dErr) const throw (cvmexception)
 Linear solver.
BaseRMatrix solve_tran (const BaseRMatrix &mB, TR &dErr) const throw (cvmexception)
 Linear solver (transposed)
BaseRMatrix solve (const BaseRMatrix &mB) const throw (cvmexception)
 Linear solver.
BaseRMatrix solve_tran (const BaseRMatrix &mB) const throw (cvmexception)
 Linear solver (transposed)
RVector operator% (const RVector &vB) const throw (cvmexception)
 Linear solver operator (transposed)
RVector operator/ (const RVector &vB) const throw (cvmexception)
 Linear solver operator.
RVector solve_lu (const basic_srmatrix &mLU, const tint *pPivots, const RVector &vB, TR &dErr) const throw (cvmexception)
 LU factorization based linear solver.
RVector solve_lu (const basic_srmatrix &mLU, const tint *pPivots, const RVector &vB) const throw (cvmexception)
 LU factorization based linear solver.
BaseRMatrix solve_lu (const basic_srmatrix &mLU, const tint *pPivots, const BaseRMatrix &mB, TR &dErr) const throw (cvmexception)
 LU factorization based linear solver.
BaseRMatrix solve_lu (const basic_srmatrix &mLU, const tint *pPivots, const BaseRMatrix &mB) const throw (cvmexception)
 LU factorization based linear solver.
TR det () const throw (cvmexception)
 Matrix determinant.
basic_srmatrixlow_up (const basic_srmatrix &m, tint *nPivots) throw (cvmexception)
 Low-up (LU) factorization.
basic_srmatrix low_up (tint *nPivots) const throw (cvmexception)
 Low-up (LU) factorization.
TR cond () const throw (cvmexception)
 Condition number reciprocal.
basic_srmatrixinv (const basic_srmatrix &m) throw (cvmexception)
 Matrix inversion.
basic_srmatrix inv () const throw (cvmexception)
 Matrix inversion.
basic_srmatrixexp (const basic_srmatrix &m, TR tol=basic_cvmMachSp< TR >()) throw (cvmexception)
 Matrix exponent.
basic_srmatrix exp (TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Matrix exponent.
basic_srmatrixpolynom (const basic_srmatrix &m, const RVector &v) throw (cvmexception)
 Matrix polynomial.
basic_srmatrix polynom (const RVector &v) const
 Matrix polynomial.
CVector eig (basic_scmatrix< TR, TC > &mEigVect, bool bRightVect=true) const throw (cvmexception)
 Eigenvalues and eigenvectors.
CVector eig () const throw (cvmexception)
 Eigenvalues.
basic_srmatrixcholesky (const basic_srsmatrix< TR > &m) throw (cvmexception)
 Cholesky factorization.
basic_srmatrixbunch_kaufman (const basic_srsmatrix< TR > &m, tint *nPivots) throw (cvmexception)
 Bunch-Kaufman factorization.
void qr (basic_srmatrix< TR > &mQ, basic_srmatrix< TR > &mR) const throw (cvmexception)
 QR factorization.
void lq (basic_srmatrix< TR > &mL, basic_srmatrix< TR > &mQ) const throw (cvmexception)
 LQ factorization.
void ql (basic_srmatrix< TR > &mQ, basic_srmatrix< TR > &mL) const throw (cvmexception)
 QL factorization.
void rq (basic_srmatrix< TR > &mR, basic_srmatrix< TR > &mQ) const throw (cvmexception)
 RQ factorization.
basic_srmatrixidentity ()
 Identity matrix.
basic_srmatrixvanish ()
 Set matrix to zero.
basic_srmatrixrandomize (TR dFrom, TR dTo)
 Randomizer.
- Public Member Functions inherited from basic_rmatrix< TR >
 basic_rmatrix ()
 Default constructor.
 basic_rmatrix (tint nM, tint nN)
 Constructor.
 basic_rmatrix (TR *pd, tint nM, tint nN)
 Constructor.
 basic_rmatrix (const TR *pd, tint nM, tint nN)
 Constructor.
 basic_rmatrix (const basic_rmatrix &m)
 Copy constructor.
 basic_rmatrix (basic_rmatrix &&m) noexcept
 Move constructor.
 basic_rmatrix (const RVector &v, bool bBeColumn=true)
 Constructor.
 basic_rmatrix (basic_rmatrix &m, tint nRow, tint nCol, tint nHeight, tint nWidth)
 Submatrix constructor.
RVector diag (tint nDiag) throw (cvmexception)
 Diagonal as l-value.
const RVector diag (tint nDiag) const throw (cvmexception)
 Diagonal (not l-value)
basic_rmatrixoperator= (const basic_rmatrix &m) throw (cvmexception)
 Assignment operator.
basic_rmatrixoperator= (basic_rmatrix &&m) throw (cvmexception)
 Move assignment operator.
basic_rmatrixassign (tint nRow, tint nCol, const basic_rmatrix &m) throw (cvmexception)
 Assignment to submatrix.
basic_rmatrixresize (tint nNewM, tint nNewN) throw (cvmexception)
 Changes dimensions.
bool operator== (const basic_rmatrix &m) const
 Matrix comparison.
bool operator!= (const basic_rmatrix &m) const
 Matrix comparison.
basic_rmatrixoperator<< (const basic_rmatrix &m) throw (cvmexception)
 Matrix replacement.
basic_rmatrix operator+ (const basic_rmatrix &m) const throw (cvmexception)
 Addition operator.
basic_rmatrix operator- (const basic_rmatrix &m) const throw (cvmexception)
 Subtraction operator.
basic_rmatrixsum (const basic_rmatrix &m1, const basic_rmatrix &m2) throw (cvmexception)
 Sum of matrices.
basic_rmatrixdiff (const basic_rmatrix &m1, const basic_rmatrix &m2) throw (cvmexception)
 Difference of matrices.
basic_rmatrixoperator+= (const basic_rmatrix &m) throw (cvmexception)
 Increment operator.
basic_rmatrixoperator-= (const basic_rmatrix &m) throw (cvmexception)
 Decrement operator.
basic_rmatrixtranspose (const basic_rmatrix &m) throw (cvmexception)
 Matrix transposition.
basic_rmatrix operator* (const basic_rmatrix &m) const throw (cvmexception)
 Matrix-matrix product.
basic_rmatrixmult (const basic_rmatrix &m1, const basic_rmatrix &m2) throw (cvmexception)
 Matrix-matrix product.
basic_rmatrixrank1update (const RVector &vCol, const RVector &vRow) throw (cvmexception)
 Rank-1 update.
basic_rmatrixsolve (const basic_srmatrix< TR > &mA, const basic_rmatrix &mB, TR &dErr) throw (cvmexception)
 Linear solver.
basic_rmatrixsolve_tran (const basic_srmatrix< TR > &mA, const basic_rmatrix &mB, TR &dErr) throw (cvmexception)
 Linear solver (transposed)
basic_rmatrixsolve (const basic_srmatrix< TR > &mA, const basic_rmatrix &mB) throw (cvmexception)
 Linear solver.
basic_rmatrixsolve_tran (const basic_srmatrix< TR > &mA, const basic_rmatrix &mB) throw (cvmexception)
 Linear solver (transposed)
basic_rmatrixsolve_lu (const basic_srmatrix< TR > &mA, const basic_srmatrix< TR > &mLU, const tint *pPivots, const basic_rmatrix &mB, TR &dErr) throw (cvmexception)
 LU factorization based linear solver.
basic_rmatrixsolve_lu (const basic_srmatrix< TR > &mA, const basic_srmatrix< TR > &mLU, const tint *pPivots, const basic_rmatrix &mB) throw (cvmexception)
 LU factorization based linear solver.
RVector svd () const throw (cvmexception)
 Singular value decomposition.
RVector svd (basic_srmatrix< TR > &mU, basic_srmatrix< TR > &mVH) const throw (cvmexception)
 Singular value decomposition.
basic_rmatrix pinv (TR threshold=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Pseudo (generalized) inversion.
basic_rmatrixpinv (const basic_rmatrix &mA, TR threshold=basic_cvmMachSp< TR >()) throw (cvmexception)
 Pseudo (generalized) inversion.
basic_rmatrix gels (bool transpose, const basic_rmatrix &mB, basic_rvector< TR > &vErr) const throw (cvmexception)
 Overdetermined or underdetermined linear solver.
basic_rmatrixgels (bool transpose, const basic_rmatrix &mA, const basic_rmatrix &mB, basic_rvector< TR > &vErr) throw (cvmexception)
 Overdetermined or underdetermined linear solver.
basic_rvector< TR > gels (bool transpose, const basic_rvector< TR > &vB, TR &dErr) const throw (cvmexception)
 Overdetermined or underdetermined linear solver.
basic_rmatrix gelsy (const basic_rmatrix &mB, tint &rank, TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Linear least squares problem.
basic_rmatrixgelsy (const basic_rmatrix &mA, const basic_rmatrix &mB, tint &rank, TR tol=basic_cvmMachSp< TR >()) throw (cvmexception)
 Linear least squares problem.
basic_rvector< TR > gelsy (const basic_rvector< TR > &vB, tint &rank, TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Linear least squares problem.
basic_rmatrix gelss (const basic_rmatrix &mB, basic_rvector< TR > &sv, tint &rank, TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Linear least squares problem.
basic_rmatrixgelss (const basic_rmatrix &mA, const basic_rmatrix &mB, basic_rvector< TR > &sv, tint &rank, TR tol=basic_cvmMachSp< TR >()) throw (cvmexception)
 Linear least squares problem.
basic_rvector< TR > gelss (const basic_rvector< TR > &vB, basic_rvector< TR > &sv, tint &rank, TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Linear least squares problem.
basic_rmatrix gelsd (const basic_rmatrix &mB, basic_rvector< TR > &sv, tint &rank, TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Linear least squares problem.
basic_rmatrixgelsd (const basic_rmatrix &mA, const basic_rmatrix &mB, basic_rvector< TR > &sv, tint &rank, TR tol=basic_cvmMachSp< TR >()) throw (cvmexception)
 Linear least squares problem.
basic_rvector< TR > gelsd (const basic_rvector< TR > &vB, basic_rvector< TR > &sv, tint &rank, TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Linear least squares problem.
tint rank (TR tol=basic_cvmMachSp< TR >()) const throw (cvmexception)
 Matrix rank.
void qr (basic_rmatrix< TR > &mQ, basic_srmatrix< TR > &mR) const throw (cvmexception)
 QR factorization ("economy" mode)
void qr (basic_srmatrix< TR > &mQ, basic_rmatrix< TR > &mR) const throw (cvmexception)
 QR factorization ("full" mode)
void lq (basic_srmatrix< TR > &mL, basic_rmatrix< TR > &mQ) const throw (cvmexception)
 LQ factorization ("economy" mode)
void lq (basic_rmatrix< TR > &mL, basic_srmatrix< TR > &mQ) const throw (cvmexception)
 LQ factorization ("full" mode)
void rq (basic_srmatrix< TR > &mR, basic_rmatrix< TR > &mQ) const throw (cvmexception)
 RQ factorization ("economy" mode)
void rq (basic_rmatrix< TR > &mR, basic_srmatrix< TR > &mQ) const throw (cvmexception)
 RQ factorization ("full" mode)
void ql (basic_rmatrix< TR > &mQ, basic_srmatrix< TR > &mL) const throw (cvmexception)
 QL factorization ("economy" mode)
void ql (basic_srmatrix< TR > &mQ, basic_rmatrix< TR > &mL) const throw (cvmexception)
 QL factorization ("full" mode)
basic_rmatrixger (TR alpha, const RVector &vCol, const RVector &vRow) throw (cvmexception)
 Rank-1 update matrix-vector operation.
basic_rmatrixgemm (const basic_rmatrix &m1, bool bTrans1, const basic_rmatrix &m2, bool bTrans2, TR dAlpha, TR dBeta) throw (cvmexception)
 Generic matrix-matrix operation.
basic_rmatrixsymm (bool bLeft, const basic_srsmatrix< TR > &ms, const basic_rmatrix &m, TR dAlpha, TR dBeta) throw (cvmexception)
 Generic symmetric matrix-matrix operation.
TR norm2 () const override
 2-norm
- Public Member Functions inherited from Matrix< TR, TR >
tint msize () const
 Number of rows.
tint nsize () const
 Number of columns.
tint ld () const
 Leading dimension.
tint rowofmax () const
 Row with maximum element.
tint rowofmin () const
 Row with minimum element.
tint colofmax () const
 Column with maximum element.
tint colofmin () const
 Column with minimum element.
TR norm1 () const override
 1-norm
TR norminf () const override
 Infinity norm.
- Public Member Functions inherited from basic_array< TR, TR >
 basic_array ()
 Default constructor.
 basic_array (tint nSize, bool bZeroMemory=true)
 Constructor.
 basic_array (TR *pd, tint nSize, tint nIncr=1)
 Constructor.
 basic_array (const TR *pd, tint nSize, tint nIncr=1)
 Constructor.
 basic_array (const TR *begin, const TR *end)
 Constructor.
 basic_array (const basic_array &a)
 Copy constructor.
 basic_array (basic_array &&a) noexcept
 Move constructor.
basic_arrayoperator= (const basic_array &a) throw (cvmexception)
 Assignment operator.
basic_arrayoperator= (basic_array &&a) throw (cvmexception)
 Move assignment operator.
virtual ~basic_array ()
 Destructor.
tint size () const
 Size (length) of array.
TR * get ()
 Pointer to data.
const TR * get () const
 Const pointer to data.
 operator TR * ()
 operator const TR * () const
void assign (size_type n, const TR &val) throw (cvmexception)
 (STL) assigns given value to n-th element (0-based)
void assign (const_iterator begin, const_iterator end) throw (cvmexception)
 (STL) assigns begin-end iteartor range to array
tint incr () const
 Increment between elements.
tint indofmax () const
 Index of element with maximum module.
tint indofmin () const
 Index of element with minimum module.
virtual TR norm () const
 Euclidean norm.
iterator begin ()
 (STL) iterator to begin
const_iterator begin () const
 (STL) const iterator to begin
iterator end ()
 (STL) iterator to end
const_iterator end () const
 (STL) const iterator to end
reverse_iterator rbegin ()
 (STL) iterator to begin reversed
const_reverse_iterator rbegin () const
 (STL) const iterator to begin reversed
reverse_iterator rend ()
 (STL) iterator to end reversed
const_reverse_iterator rend () const
 (STL) const iterator to end reversed
size_type max_size () const
 (STL) maximum possible size of array
size_type capacity () const
 (STL) current capacity of array, equal to size()
bool empty () const
 (STL) is array empty
reference front ()
 (STL) reference to first element
const_reference front () const
 (STL) const reference to first element
reference back ()
 (STL) reference to last element
const_reference back () const
 (STL) const reference to last element
void clear ()
 (STL) clears array, dealocates memory and sets size() to zero
void swap (basic_array &v) throw (cvmexception)
 (STL) swaps array values, throws cvmexception if sizes are different
reference at (size_type n) throw (cvmexception)
 (STL) returns reference to n-th element of array (0-based), throws cvmexception if n is out of boundaries
const_reference at (size_type n) const throw (cvmexception)
 (STL) returns const reference to n-th element of array (0-based), throws cvmexception if n is out of boundaries
void push_back (const TR &x) throw (cvmexception)
 (STL) pushes new value to the end of array
void pop_back () throw (cvmexception)
 (STL) removes last element from array
iterator insert (iterator position, const TR &x) throw (cvmexception)
 (STL) inserts element to given position in array
iterator erase (iterator position) throw (cvmexception)
 (STL) removes element from given position in array

Additional Inherited Members

- Public Types inherited from basic_array< TR, TR >
typedef TR value_type
 STL-specific value type definition.
typedef value_typepointer
 STL-specific value pointer definition.
typedef value_typeiterator
 STL-specific iterator definition.
typedef const value_typeconst_iterator
 STL-specific const iterator definition.
typedef const value_typeconst_pointer
 STL-specific const pointer definition.
typedef value_typereference
 STL-specific reference definition.
typedef const value_typeconst_reference
 STL-specific const reference definition.
typedef size_t size_type
 STL-specific size type definition.
typedef ptrdiff_t difference_type
 STL-specific difference type definition.
typedef std::reverse_iterator
< const_iterator
const_reverse_iterator
 STL-specific const reverse iterator definition.
typedef std::reverse_iterator
< iterator
reverse_iterator
 STL-specific reverse iterator definition.
- Protected Member Functions inherited from Matrix< TR, TR >
 Matrix ()
 Default constructor.
 Matrix (tint nM, tint nN, tint nLD, bool bZeroMemory)
 Constructor.
 Matrix (TR *pd, tint nM, tint nN, tint nLD, tint nSize)
 Constructor.
 Matrix (const TR *pd, tint nM, tint nN, tint nLD, tint nSize)
 Constructor.
 Matrix (const BaseArray &v, bool bBeColumn)
 Constructor.
 Matrix (Matrix &&m) noexcept
 Move constructor.
Matrixoperator= (Matrix &&m) throw (cvmexception)
 Move assignment operator.
- Protected Member Functions inherited from SqMatrix< TR, TR >
 SqMatrix ()
 internal constructor
virtual ~SqMatrix ()
- Protected Attributes inherited from Matrix< TR, TR >
tint mm
 Number of rows.
tint mn
 Number of columns.
tint mld
 Leading dimension.

Detailed Description

template<typename TR>
class basic_srmatrix< TR >

End-user class encapsulating square matrix of real numbers.

TR type stands for treal. Please use predefined srmatrix class in your applications.

See Also
Matrix

Definition at line 14940 of file cvm.h.

Constructor & Destructor Documentation

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( )
inline

Default constructor.

Creates empty square matrix. No memory gets allocated.

Example:
using namespace cvm;
std::cout << m.msize() << std::endl << m.nsize() << std::endl;
std::cout << m.size() << std::endl;
m.resize (3);
std::cout << m;
prints
0
0
0
0 0 0
0 0 0
0 0 0

Definition at line 14976 of file cvm.h.

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( tint  nDim)
inlineexplicit

Constructor.

Creates $n\times n$ srmatrix object where $n$ is passed in nDim parameter. Constructor sets all elements to zero. It throws cvmexception in case of non-positive size passed or memory allocation failure. Example:

using namespace cvm;
srmatrix m (4);
std::cout << m.msize() << std::endl << m.nsize()
<< std::endl << m.size() << std::endl << m;

prints

4
4
16
0 0 0 0
0 0 0 0
0 0 0 0
0 0 0 0
Parameters
[in]nDimNumber of rows and columns.

Definition at line 15006 of file cvm.h.

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( TR *  pd,
tint  nDim 
)
inline

Constructor.

Creates $n\times n$ srmatrix object where $n$ is passed in nDim parameter. It throws cvmexception in case of non-positive size passed. Unlike others, this constructor does not allocate memory. It just shares memory with array pointed to by pd (for matrices nIncr=1 is always satisfied).

Example:
using namespace cvm;
double a[] = {1., 1., 1., 1., 1., 1., 1., 1., 1.};
srmatrix m (a, 3);
m(1,1) = 5.;
std::cout << m << std::endl;
std::cout << a[0] << " " << a[1] << " " << a[2] << std::endl;
prints
5 1 1
1 1 1
1 1 1
5 1 1
See Also
http://cvmlib.com/faq.htm
basic_rmatrix::basic_rmatrix(TR*,tint,tint)
Parameters
[in]pdPointer to array to share memory with.
[in]nDimNumber of rows and columns.

Definition at line 15041 of file cvm.h.

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( const TR *  pd,
tint  nDim 
)
inline

Constructor.

Creates $n\times n$ srmatrix object where $n$ is passed in nDim parameter. Constructor throws cvmexception in case of non-positive sizes passed. This is const version, it allocates memory and copies every element (deep copy) from external array pointed to by pd parameter. It copies nDim*nDim elements total. Example:

using namespace cvm;
const double a[] = {1., 1., 1., 1., 1., 1., 1., 1., 1.};
srmatrix m (a, 3);
m(1,1) = 5.;
std::cout << m << std::endl;
std::cout << a[0] << " " << a[1] << " " << a[2] << std::endl;

prints

5 1 1
1 1 1
1 1 1
1 1 1
See Also
http://cvmlib.com/faq.htm
basic_rmatrix::basic_rmatrix(const TR*,tint,tint)
Parameters
[in]pdConst pointer to external array.
[in]nDimNumber of rows and columns.

Definition at line 15077 of file cvm.h.

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( const basic_srmatrix< TR > &  m)
inline

Copy constructor.

Creates srmatrix object as a copy of matrix m. It throws cvmexception in case of memory allocation failure.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m(a, 3);
srmatrix mc(m);
m(1,1) = 7.77;
std::cout << m << std::endl << mc;
prints
7.77 4 7
2 5 8
3 6 9
1 4 7
2 5 8
3 6 9
Parameters
[in]msrmatrix to copy from.

Definition at line 15110 of file cvm.h.

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( basic_srmatrix< TR > &&  m)
inline

Move constructor.

Implements move semantics introduced in new C++ standard. Moves data ownership from other matrix to newly created object. It's usually executed implicitly in cases like this:

rvector a(b + c);

or this

rvector a = b + c;

Here temporary result of calling b.operator+(c) will not be destroyed but rather moved to newly created object a.

Parameters
[in]mrvalue reference to other matrix.

Definition at line 15133 of file cvm.h.

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( const BaseRMatrix m)
inline

Constructor.

Creates srmatrix object as a copy of matrix m. It's assumed that $m\times n$ matrix m must have equal sizes, i.e. $m = n$ is satisfied. Constructor throws cvmexception if this is not true or in case of memory allocation failure. Please note that this constructor is not explicit.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m(a, 2, 3);
std::cout << m << std::endl;
m.resize(3, 3);
srmatrix ms (m);
std::cout << ms;
prints
1 3 5
2 4 6
1 3 5
2 4 6
0 0 0
Parameters
[in]mrmatrix to copy from.

Definition at line 15169 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( const RVector v)
inlineexplicit

Constructor.

Creates srmatrix object of size v.size() by v.size() and assigns vector v to its main diagonal. Constructor throws cvmexception in case of memory allocation failure.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5.};
rvector v(a, 5);
srmatrix m(v);
std::cout << m;
prints
1 0 0 0 0
0 2 0 0 0
0 0 3 0 0
0 0 0 4 0
0 0 0 0 5
Parameters
[in]vrvector to copy main diagonal from.

Definition at line 15201 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix< TR >::basic_srmatrix ( BaseRMatrix m,
tint  nRow,
tint  nCol,
tint  nDim 
)
inline

Submatrix constructor.

Creates srmatrix object as submatrix of matrix m. It means that the object created shares memory with some part of m. This part is defined by its upper left corner (parameters nRow and nCol, both are CVM0 based) and its dimension (parameter nDim).

Example:
using namespace cvm;
rmatrix m(4,5);
srmatrix subm(m, 2, 2, 2);
subm.set(1.);
std::cout << m;
prints
0 0 0 0 0
0 1 1 0 0
0 1 1 0 0
0 0 0 0 0
Parameters
[in]mParent rmatrix to attach to.
[in]nRowRow to start from (CVM0 based).
[in]nColColumn to start from (CVM0 based).
[in]nDimDimension of square submatrix.

Definition at line 15235 of file cvm.h.

Member Function Documentation

template<typename TR>
type_proxy<TR,TR> basic_srmatrix< TR >::operator() ( tint  nRow,
tint  nCol 
) throw (cvmexception)
inline

Reference to element (l-value)

Operator provides access to a particular element of calling matrix by its row and column index. Indexes passed are CVM0 based. It returns l-value in order to make possible write access to an element. Operator throws cvmexception if nRow or nCol is outside of boundaries.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6.};
const rmatrix m (a, 2, 3);
rmatrix ms(m);
std::cout << m(1,1) << " " << m(2,3) << std::endl << std::endl;
ms(2,2) = 7.77;
std::cout << ms;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1.00e+00 6.00e+00
1.00e+00 3.00e+00 5.00e+00
2.00e+00 7.77e+00 6.00e+00
See Also
http://cvmlib.com/faq.htm
Parameters
[in]nRowRow index (CVM0 based).
[in]nColColumn index (CVM0 based).
Returns
type_proxy Proxy to element (l-value).

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srbmatrix< TR >.

Definition at line 15241 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
TR basic_srmatrix< TR >::operator() ( tint  nRow,
tint  nCol 
) const throw (cvmexception)
inline

Value of element (not l-value)

Operator returns value of a particular element of calling matrix by its row and column index. Indexes passed are CVM0 based. Operator throws cvmexception if nRow or nCol is outside of boundaries.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6.};
const rmatrix m (a, 2, 3);
std::cout << m(1,1) << " " << m(2,3) << std::endl;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1.00e+00 6.00e+00
See Also
http://cvmlib.com/faq.htm
Parameters
[in]nRowRow index (CVM0 based).
[in]nColColumn index (CVM0 based).
Returns
TR treal value.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15248 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
RVector basic_srmatrix< TR >::operator() ( tint  nCol) throw (cvmexception)
inline

Column as l-value.

Operator provides access to nCol-th column (CVM0 based) of calling matrix by returning rvector sharing memory with it. Operator throws cvmexception if nCol is outside of boundaries.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6.};
const rmatrix m (a, 2, 3);
srmatrix ms(2);
std::cout << m(2) << std::endl;
ms(2) = m(3);
std::cout << ms;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
3.00e+00 4.00e+00
0.00e+00 5.00e+00
0.00e+00 6.00e+00
Parameters
[in]nColIndex of column (CVM0 based).
Returns
rvector Column as l-value.

Reimplemented from basic_rmatrix< TR >.

Definition at line 15256 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
const RVector basic_srmatrix< TR >::operator() ( tint  nCol) const throw (cvmexception)
inline

Column as not l-value.

Operator creates rvector object as a copy of nCol-th column (CVM0 based) of calling matrix. Operator throws cvmexception if nCol is outside of boundaries.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6.};
const rmatrix m (a, 2, 3);
std::cout << m(2) << std::endl;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
3.00e+00 4.00e+00
Parameters
[in]nColIndex of column (CVM0 based).
Returns
rvector Column value.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15263 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
RVector basic_srmatrix< TR >::operator[] ( tint  nRow) throw (cvmexception)
inline

Row as l-value.

Operator provides access to nRow-th row (CVM0 based) of calling matrix by returning rvector sharing memory with it. Operator throws cvmexception if nRow is outside of boundaries.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6.};
const rmatrix m (a, 2, 3);
srmatrix ms(3);
std::cout << m[1] << std::endl;
ms[1] = m[2];
std::cout << ms;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1.00e+00 3.00e+00 5.00e+00
2.00e+00 4.00e+00 6.00e+00
0.00e+00 0.00e+00 0.00e+00
0.00e+00 0.00e+00 0.00e+00
Parameters
[in]nRowIndex of row (CVM0 based).
Returns
rvector Row as l-value.

Reimplemented from basic_rmatrix< TR >.

Definition at line 15270 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
const RVector basic_srmatrix< TR >::operator[] ( tint  nRow) const throw (cvmexception)
inline

Row as not l-value.

Operator creates rvector object as a copy of nRow-th row (CVM0 based) of calling matrix. Operator throws cvmexception if nRow is outside of boundaries.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6.};
const rmatrix m (a, 2, 3);
std::cout << m[1] << std::endl;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1.00e+00 3.00e+00 5.00e+00
Parameters
[in]nRowIndex of row (CVM0 based).
Returns
rvector Row value.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15277 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator= ( const basic_srmatrix< TR > &  m) throw (cvmexception)
inline

Assignment operator.

Sets every element of calling srmatrix to be equal to appropriate element of smatrix m and returns reference to the matrix changed. Operator throws cvmexception in case of different matrix dimensions.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
const srmatrix m1(a, 3);
srmatrix m2(3);
m2 = m1;
std::cout << m2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1.00e+00 4.00e+00 7.00e+00
2.00e+00 5.00e+00 8.00e+00
3.00e+00 6.00e+00 9.00e+00
Parameters
[in]msrmatrix to assign from.
Returns
Reference to changed calling matrix.

Definition at line 15316 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator= ( basic_srmatrix< TR > &&  m) throw (cvmexception)
inline

Move assignment operator.

Implements move semantics introduced in new C++ standard. Moves data ownership from other matrix to calling object. It's usually executed implicitly in cases like this:

a = b + c;

Here temporary result of calling b.operator+(c) will not be destroyed but rather moved to calling object a.

Parameters
[in]mrvalue reference to other matrix.

Definition at line 15336 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::assign ( const RVector v) throw (cvmexception)
inline

Vector (as array) assignment.

Sets every element of calling matrix to be equal to appropriate element of rvector v as an array. Assignment is performed according to matrix storage (by columns). It's assumed that vector passed is long enough to fill calling matrix. Function throws cvmexception otherwise.

Example:
try {
const double a[] = {1., 2., 3., 4., 5., 6., 7., 8.};
const rvector v(a,8);
rmatrix m(2,3);
std::cout << v << std::endl;
m.assign(v);
std::cout << m;
} catch (std::exception& e) {
std::cout << "Exception: " << e.what () << std::endl;
}
prints
1 2 3 4 5 6 7 8
1 3 5
2 4 6
Parameters
[in]vrvector to assign.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srbmatrix< TR >.

Definition at line 15344 of file cvm.h.

Here is the call graph for this function:

Here is the caller graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::assign ( const TR *  pd)
inline

External array assignment.

Sets every element of calling matrix to be equal to appropriate element of array pointed to by parameter pd and returns reference to the matrix changed. Assignment is performed according to matrix storage (by columns). It's assumed that array passed is long enough to fill calling matrix.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
const double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m(2, 3);
m.assign(a);
std::cout << m;
prints
1.00e+000 3.00e+000 5.00e+000
2.00e+000 4.00e+000 6.00e+000
Parameters
[in]pdConst pointer to external array.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srbmatrix< TR >.

Definition at line 15352 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::assign ( tint  nRow,
tint  nCol,
const BaseRMatrix m 
) throw (cvmexception)
inline

Assignment to submatrix.

Sets submatrix of calling matrix beginning with row nRow and column nCol to matrix m and returns reference to the matrix changed. Function throws cvmexception if nRow or nCol are not positive or matrix m doesn't fit. Indexes are CVM0 based.

Example:
using namespace cvm;
srmatrix m1(5);
rmatrix m2(2,3);
m1.set(1.);
m2.set(2.);
m1.assign(2,3,m2);
std::cout << m1;
prints
1 1 1 1 1
1 1 2 2 2
1 1 2 2 2
1 1 1 1 1
1 1 1 1 1
Parameters
[in]nRowRow index (CVM0 based).
[in]nColColumn index (CVM0 based).
[in]mReference to rmatrix to assign.
Returns
Reference to changed calling matrix.

Definition at line 15389 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::set ( TR  d)
inline

Sets all elements to one value.

Sets every element of calling matrix to be equal to parameter d and returns reference to the matrix changed.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
rmatrix m(2, 3);
m.set(3.);
std::cout << m;
prints
3.00e+000 3.00e+000 3.00e+000
3.00e+000 3.00e+000 3.00e+000
Parameters
[in]dValue to set to.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15400 of file cvm.h.

Here is the caller graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::resize ( tint  nNewDim) throw (cvmexception)
inline

Changes dimension.

Changes dimension of calling square matrix to nNewDim and returns reference to the matrix changed. In case of increasing of its size, the matrix is filled up with zeroes. Function throws cvmexception in case of negative dimension passed or memory allocation failure.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4.};
srmatrix m(a, 2);
std::cout << m << std::endl;
m.resize (3);
std::cout << m;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1 3
2 4
1 3 0
2 4 0
0 0 0
Parameters
[in]nNewDimNew dimension.
Returns
Reference to changed calling matrix.

Reimplemented from basic_array< TR, TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15441 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator<< ( const basic_srmatrix< TR > &  m) throw (cvmexception)
inline

Matrix replacement.

Destroys calling matrix, creates a new one as a copy of m and returns reference to the matrix changed. Operator throws cvmexception in case of memory allocation failure.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
srmatrix m(3);
srmatrix mc(1);
m(1,2) = 1.;
m(2,3) = 2.;
std::cout << m << std::endl << mc << std::endl;
mc << m;
std::cout << mc;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
0.00e+00 1.00e+00 0.00e+00
0.00e+00 0.00e+00 2.00e+00
0.00e+00 0.00e+00 0.00e+00
0.00e+00
0.00e+00 1.00e+00 0.00e+00
0.00e+00 0.00e+00 2.00e+00
0.00e+00 0.00e+00 0.00e+00
See Also
operator =()
Parameters
[in]msrmatrix to replace by.
Returns
Reference to changed calling matrix.

Definition at line 15487 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator+ ( const basic_srmatrix< TR > &  m) const throw (cvmexception)
inline

Addition operator.

Creates object of type srmatrix as sum of calling matrix and matrix m. Operator throws cvmexception in case of different sizes of the operands or memory allocation failure.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
double b[] = {10., 20., 30., 40., 50., 60., 70., 80., 90.};
srmatrix m1(a, 3);
srmatrix m2(b, 3);
std::cout << m1 + m2 << std::endl << m1 + m1;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
11 44 77
22 55 88
33 66 99
2 8 14
4 10 16
6 12 18
See Also
sum()
Parameters
[in]msrmatrix to add to calling one.
Returns
Sum of matrices.

Definition at line 15530 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator- ( const basic_srmatrix< TR > &  m) const throw (cvmexception)
inline

Subtraction operator.

Creates object of type srmatrix as difference of calling matrix and matrix m. It throws cvmexception in case of different sizes of the operands or memory allocation failure.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
double b[] = {10., 20., 30., 40., 50., 60., 70., 80., 90.};
srmatrix m1(a, 3);
srmatrix m2(b, 3);
std::cout << m2 - m1 << std::endl << m1 - m1;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
9 36 63
18 45 72
27 54 81
0 0 0
0 0 0
0 0 0
See Also
diff()
Parameters
[in]msrmatrix to subtract from calling one.
Returns
Difference of matrices.

Definition at line 15573 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::sum ( const basic_srmatrix< TR > &  m1,
const basic_srmatrix< TR > &  m2 
) throw (cvmexception)
inline

Sum of matrices.

Assigns sum of square matrices m1 and m2 to calling matrix and returns reference to the matrices changed. It throws cvmexception in case of different sizes of the operands.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
const srmatrix m1(a, 3);
srmatrix m2(3);
srmatrix m(3);
m2.set(1.);
std::cout << m.sum(m1, m2) << std::endl;
std::cout << m.sum(m, m2);
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
2 5 8
3 6 9
4 7 10
3 6 9
4 7 10
5 8 11
See Also
operator +()
Parameters
[in]m1First srmatrix summand.
[in]m2Second srmatrix summand.
Returns
Reference to changed calling matrix.

Definition at line 15618 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::diff ( const basic_srmatrix< TR > &  m1,
const basic_srmatrix< TR > &  m2 
) throw (cvmexception)
inline

Difference of matrices.

Assigns difference of square matrices m1 and m2 to calling matrix and returns reference to the matrix changed. It throws cvmexception in case of different sizes of the operands.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
const srmatrix m1(a, 3);
srmatrix m2(3);
srmatrix m(3);
m2.set(1.);
std::cout << m.diff(m1, m2) << std::endl;
std::cout << m.diff(m, m2);
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
0 3 6
1 4 7
2 5 8
-1 2 5
0 3 6
1 4 7
See Also
operator -()
Parameters
[in]m1First srmatrix subtrahend.
[in]m2Second srmatrix subtrahend.
Returns
Reference to changed calling matrix.

Definition at line 15663 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator+= ( const basic_srmatrix< TR > &  m) throw (cvmexception)
inline

Increment operator.

Adds srmatrix m to calling matrix and returns reference to the matrix changed. It throws cvmexception in case of different sizes of the operands.

Example:
using namespace cvm;
try {
srmatrix m1(3);
srmatrix m2(3);
m1.set(1.);
m2.set(2.);
m1 += m2;
std::cout << m1 << std::endl;
// well, you can do this too, but temporary object would be created
m2 += m2;
std::cout << m2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
3 3 3
3 3 3
3 3 3
4 4 4
4 4 4
4 4 4
See Also
operator +()
sum()
Parameters
[in]msrmatrix to increment by.
Returns
Reference to changed calling matrix.

Definition at line 15710 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator-= ( const basic_srmatrix< TR > &  m) throw (cvmexception)
inline

Decrement operator.

Subtracts srmatrix m from calling matrix and returns reference to the matrix changed. It throws cvmexception in case of different sizes of the operands.

Example:
using namespace cvm;
try {
srmatrix m1(3);
srmatrix m2(3);
m1.set(1.);
m2.set(2.);
m1 -= m2;
std::cout << m1 << std::endl;
// well, you can do this too, but temporary object would be created
m2 -= m2;
std::cout << m2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-1 -1 -1
-1 -1 -1
-1 -1 -1
0 0 0
0 0 0
0 0 0
See Also
operator -()
diff()
Parameters
[in]msrmatrix to decrement by.
Returns
Reference to changed calling matrix.

Definition at line 15758 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator- ( ) const
inline

Unary minus operator.

Creates object of type srmatrix as calling matrix multiplied by -1. It throws cvmexception in case of memory allocation failure.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m(a, 3);
std::cout << -m;
prints
-1 -4 -7
-2 -5 -8
-3 -6 -9
Returns
Result object.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15786 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator++ ( )
inline

Plus identity, prefix.

Adds identity matrix to calling square matrix and returns reference to the matrix changed.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m (a, 3);
m++;
std::cout << m << std::endl;
std::cout << ++m;
prints
2 4 7
2 6 8
3 6 10
3 4 7
2 7 8
3 6 11
Returns
Reference to changed calling matrix.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15819 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator++ ( int  )
inline

Plus identity, postfix.

Adds identity matrix to calling square matrix and returns reference to the matrix changed.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m (a, 3);
m++;
std::cout << m << std::endl;
std::cout << ++m;
prints
2 4 7
2 6 8
3 6 10
3 4 7
2 7 8
3 6 11
Returns
Reference to changed calling matrix.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15850 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator-- ( )
inline

Minus identity, prefix.

Subtracts identity matrix from calling square matrix and returns reference to the matrix changed.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m (a, 3);
m--;
std::cout << m << std::endl;
std::cout << --m;
prints
0 4 7
2 4 8
3 6 8
-1 4 7
2 3 8
3 6 7
Returns
Reference to changed calling matrix.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15881 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator-- ( int  )
inline

Minus identity, postfix.

Subtracts identity matrix from calling square matrix and returns reference to the matrix changed.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m (a, 3);
m--;
std::cout << m << std::endl;
std::cout << --m;
prints
0 4 7
2 4 8
3 6 8
-1 4 7
2 3 8
3 6 7
Returns
Reference to changed calling matrix.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15912 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator* ( TR  dMult) const
inline

Multiply by number operator.

Creates object of type srmatrix as product of calling square matrix and number dMult. It throws cvmexception in case of memory allocation failure.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m(a, 3);
std::cout << m * 5.;
prints
5 20 35
10 25 40
15 30 45
See Also
operator *=()
Parameters
[in]dMultNumber to multiply by.
Returns
Result object.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15941 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator/ ( TR  dDiv) const throw (cvmexception)
inline

Divide by number operator.

Creates object of type srmatrix as quotient of calling square matrix and number dDiv. It throws cvmexception if dDiv has absolute value equal or less than cvmMachMin() (the smallest normalized positive number). It also throws exception in case of memory allocation failure.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m(a, 3);
std::cout << m / 4.;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
0.25 1 1.75
0.5 1.25 2
0.75 1.5 2.25
See Also
operator /=()
Parameters
[in]dDivNumber to divide by.
Returns
Result object.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15978 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator*= ( TR  dMult)
inline

Multiply by number and assign.

Multiplies calling matrix by number dMult and returns reference to the matrix changed.

Example:
using namespace cvm;
double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m(a, 2, 3);
m *= 2.;
std::cout << m;
prints
2 6 10
4 8 12
See Also
operator *(TR) const
Parameters
[in]dMultNumber to multiply by.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15985 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator/= ( TR  dDiv)
inline

Divide by number and assign.

Divides calling matrix by number dDiv and returns reference to the matrix changed. It throws cvmexception if dDiv has absolute value equal or less than cvmMachMin() (the smallest normalized positive number).

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
try {
double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m(a, 2, 3);
m /= 2.;
std::cout << m;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
5.00e-01 1.50e+00 2.50e+00
1.00e+00 2.00e+00 3.00e+00
See Also
operator /(TR) const
Parameters
[in]dDivNumber to divide by.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15991 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::normalize ( )
inline

Matrix normalizer.

Normalizes calling matrix so its Euclidean norm() becomes equal to 1 if it was greater than cvmMachMin() (the smallest normalized positive number) before the call. Does nothing otherwise.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (2);
double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m(a, 2, 3);
m.normalize();
std::cout << m << m.norm() << std::endl;
prints
1.05e-01 3.14e-01 5.24e-01
2.10e-01 4.19e-01 6.29e-01
1
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 15997 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator~ ( ) const throw (cvmexception)
inline

Matrix transposition.

Creates object of type srmatrix as transposed calling matrix.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m(a,3);
srmatrix mt(3);
std::cout << ~m << std::endl ;
mt.transpose(m);
std::cout << mt << std::endl;
mt.transpose();
std::cout << mt;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1 2 3
4 5 6
7 8 9
1 2 3
4 5 6
7 8 9
1 4 7
2 5 8
3 6 9
Returns
Result object.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 16040 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::transpose ( const basic_srmatrix< TR > &  m) throw (cvmexception)
inline

Matrix transposition.

Sets calling matrix to be equal to square matrix m transposed. Function throws cvmexception in case of not appropriate sizes of the operands.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m(a,3);
srmatrix mt(3);
std::cout << ~m << std::endl ;
mt.transpose(m);
std::cout << mt << std::endl;
mt.transpose();
std::cout << mt;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1 2 3
4 5 6
7 8 9
1 2 3
4 5 6
7 8 9
1 4 7
2 5 8
3 6 9
Parameters
[in]msrmatrix to transpose.
Returns
Reference to changed calling matrix.

Definition at line 16085 of file cvm.h.

Here is the caller graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::transpose ( )
inline

Matrix transposition (in-place)

Makes calling matrix to be equal to transposed itself. Function throws cvmexception in case of memory allocation failure.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m(a,2,3);
std::cout << m << std::endl;
std::cout << m.transpose();
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1 3 5
2 4 6
1 2
3 4
5 6
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 16091 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::operator* ( const RVector v) const throw (cvmexception)
inline

Matrix-vector product.

Creates object of type rvector as product of calling matrix and vector v. Function throws cvmexception if the number of columns of the calling matrix differs from size of the vector v. Use basic_rvector::mult(const basic_rmatrix<TR>&,const basic_rvector<TR>&) to avoid new object creation.

Example:
using namespace cvm;
try {
rmatrix m(2, 3);
rvector v(3);
m.set(1.);
v.set(1.);
std::cout << m * v;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
3 3
See Also
basic_rvector::mult(const basic_rmatrix<TR>&,const basic_rvector<TR>&)
Parameters
[in]vrvector to compute product with.
Returns
Result object.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 16097 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
BaseRMatrix basic_srmatrix< TR >::operator* ( const BaseRMatrix m) const throw (cvmexception)
inline

Matrix-matrix product.

Creates object of type rmatrix as product of calling square matrix and matrix m. Operator throws cvmexception if number of columns of calling matrix differs from number of rows of the matrix m. Use mult() to avoid new object creation.

Example:
using namespace cvm;
try {
srmatrix ms(3);
rmatrix m(3,2);
ms.set(1.);
m.set(1.);
std::cout << ms * m;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
3 3
3 3
3 3
Parameters
[in]mrmatrix to compute product with.
Returns
Result object.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 16132 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::operator* ( const basic_srmatrix< TR > &  m) const throw (cvmexception)
inline

Matrix-matrix product.

Creates object of type srmatrix as product of calling matrix and square matrix m. Operator throws cvmexception if number of columns of calling matrix differs from number of rows of the matrix m. Use mult() to avoid new object creation.

Example:
using namespace cvm;
try {
srmatrix m1(3);
srmatrix m2(3);
m1.set(1.);
m2.set(1.);
std::cout << m1 * m2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
3 3 3
3 3 3
3 3 3
Parameters
[in]msrmatrix to compute product with.
Returns
Result object.

Definition at line 16166 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::operator*= ( const basic_srmatrix< TR > &  m) throw (cvmexception)
inline

Matrix-matrix product with assignment.

Sets calling matrix to be equal to product of itself and square matrix m and returns reference to the object it changes. Operator throws cvmexception in case of different dimensions of the operands.

Example:
using namespace cvm;
try {
srmatrix m1(3);
srmatrix m2(3);
m1.set(1.);
m2.set(1.);
m1 *= m2;
std::cout << m1 << std::endl;
m1 *= m1;
std::cout << m1;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
3 3 3
3 3 3
3 3 3
27 27 27
27 27 27
27 27 27
Parameters
[in]msrmatrix to compute product with.
Returns
Reference to changed calling matrix.

Definition at line 16210 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::swap_rows ( tint  n1,
tint  n2 
) throw (cvmexception)
inline

Rows swap.

Swaps two rows of calling matrix and returns reference to the matrix changed. n1 and n2 are indexes of rows to be swapped, both are CVM0 based. Function throws cvmexception if one of the parameters is outside of boundaries.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m (a, 3, 2);
std::cout << m << std::endl;
std::cout << m.swap_rows(2,3);
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1 4
2 5
3 6
1 4
3 6
2 5
Parameters
[in]n1Row index to swap.
[in]n2Row index to swap.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Definition at line 16218 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::swap_cols ( tint  n1,
tint  n2 
) throw (cvmexception)
inline

Columns swap.

Swaps two columnss of calling matrix and returns reference to the matrix changed. n1 and n2 are indexes of columns to be swapped, both are CVM0 based. Function throws cvmexception if one of the parameters is outside of boundaries.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 3., 4., 5., 6.};
rmatrix m (a, 2, 3);
std::cout << m << std::endl;
std::cout << m.swap_cols(2,3);
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1 3 5
2 4 6
1 5 3
2 6 4
Parameters
[in]n1Column index to swap.
[in]n2Column index to swap.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Definition at line 16224 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::solve ( const RVector vB,
TR &  dErr 
) const throw (cvmexception)
inline

Linear solver.

Creates rvector object as solution $x$ of linear equation $A*x=b$ where calling matrix is square matrix $A$ and parameter vB is vector $b$. It also sets output parameter dErr to be equal to the norm of computation error. Function throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
double dErr;
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve (mb, dErr);
std::cout << mx << dErr
<< std::endl << ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve (vb, dErr);
std::cout << vx << dErr << std::endl << ma * vx - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-3.333e-01 -6.667e-01 -1.000e+00 1.000e+00
3.333e-01 6.667e-01 1.000e+00 0.000e+00
6.661e-16 1.332e-15 0.000e+00 0.000e+00
3.301e-14
0.000e+00 0.000e+00 0.000e+00 0.000e+00
-1.110e-16 -2.220e-16 0.000e+00 0.000e+00
2.220e-16 4.441e-16 0.000e+00 0.000e+00
-6.667e-01 6.667e-01 1.332e-15
3.301e-14
0.000e+00 -2.220e-16 4.441e-16
See Also
solve(const RVector&)const
Parameters
[in]vBrvector $b$.
[out]dErrNorm of computation error.
Returns
Result object.

Definition at line 16291 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::solve_tran ( const RVector vB,
TR &  dErr 
) const throw (cvmexception)
inline

Linear solver (transposed)

Creates rvector object as solution $x$ of linear equation $A^T*x=b$ (which is equivalent to $x*A=b$) where calling matrix is square matrix $A$ and parameter vB is vector $b$. It also sets output parameter dErr to be equal to the norm of computation error. Function throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
double dErr;
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve_tran (mb, dErr);
std::cout << mx << dErr
<< std::endl << ~ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve_tran (vb, dErr);
std::cout << vx << dErr << std::endl << vx * ma - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-1.000e+000 -2.000e+000 -3.000e+000 -3.333e-001
+1.000e+000 +2.000e+000 +3.000e+000 +6.667e-001
+0.000e+000 +0.000e+000 -1.332e-015 +0.000e+000
+3.513e-014
+0.000e+000 +0.000e+000 +0.000e+000 +0.000e+000
+0.000e+000 +0.000e+000 +8.882e-016 -2.220e-016
+0.000e+000 +0.000e+000 -2.665e-015 +0.000e+000
-2.000e+000 +2.000e+000 +0.000e+000
+3.168e-014
+0.000e+000 +0.000e+000 +0.000e+000
See Also
solve_tran(const RVector&)const
Parameters
[in]vBrvector $b$.
[out]dErrNorm of computation error.
Returns
Result object.

Definition at line 16356 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::solve ( const RVector vB) const throw (cvmexception)
inline

Linear solver.

Creates rvector object as solution $x$ of linear equation $A*x=b$ where calling matrix is square matrix $A$ and parameter vB is vector $b$. Function throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve (mb);
std::cout << mx << std::endl << ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve (vb);
std::cout << vx << std::endl << ma * vx - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-3.333e-01 -6.667e-01 -1.000e+00 1.000e+00
3.333e-01 6.667e-01 1.000e+00 0.000e+00
6.661e-16 1.332e-15 0.000e+00 0.000e+00
0.000e+00 0.000e+00 0.000e+00 0.000e+00
-1.110e-16 -2.220e-16 0.000e+00 0.000e+00
2.220e-16 4.441e-16 0.000e+00 0.000e+00
-6.667e-01 6.667e-01 1.332e-15
0.000e+00 -2.220e-16 4.441e-16
See Also
solve(const RVector&,TR&)const
Parameters
[in]vBrvector $b$.
Returns
Result object.

Definition at line 16415 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::solve_tran ( const RVector vB) const throw (cvmexception)
inline

Linear solver (transposed)

Creates rvector object as solution $x$ of linear equation $A^T*x=b$ (which is equivalent to $x*A=b$) where calling matrix is square matrix $A$ and parameter vB is vector $b$. Function throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve_tran (mb, dErr);
std::cout << mx << std::endl << ~ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve_tran (vb, dErr);
std::cout << vx << std::endl << vx * ma - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-1.000e+000 -2.000e+000 -3.000e+000 -3.333e-001
+1.000e+000 +2.000e+000 +3.000e+000 +6.667e-001
+0.000e+000 +0.000e+000 -1.332e-015 +0.000e+000
+0.000e+000 +0.000e+000 +0.000e+000 +0.000e+000
+0.000e+000 +0.000e+000 +8.882e-016 -2.220e-016
+0.000e+000 +0.000e+000 -2.665e-015 +0.000e+000
-2.000e+000 +2.000e+000 +0.000e+000
+0.000e+000 +0.000e+000 +0.000e+000
See Also
solve_tran(const RVector&,TR&)const
Parameters
[in]vBrvector $b$.
Returns
Result object.

Definition at line 16477 of file cvm.h.

template<typename TR>
BaseRMatrix basic_srmatrix< TR >::solve ( const BaseRMatrix mB,
TR &  dErr 
) const throw (cvmexception)
inline

Linear solver.

Creates rmatrix object as solution $X$ of matrix linear equation $A*X=B$ where calling square matrix is matrix $A$ and parameter mB is matrix $B$. It also sets output parameter dErr to be equal to the norm of computation error. Function throws cvmexception in case of inappropriate sizes of the operands or when the matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
double dErr;
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve (mb, dErr);
std::cout << mx << dErr
<< std::endl << ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve (vb, dErr);
std::cout << vx << dErr << std::endl << ma * vx - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-3.333e-01 -6.667e-01 -1.000e+00 1.000e+00
3.333e-01 6.667e-01 1.000e+00 0.000e+00
6.661e-16 1.332e-15 0.000e+00 0.000e+00
3.301e-14
0.000e+00 0.000e+00 0.000e+00 0.000e+00
-1.110e-16 -2.220e-16 0.000e+00 0.000e+00
2.220e-16 4.441e-16 0.000e+00 0.000e+00
-6.667e-01 6.667e-01 1.332e-15
3.301e-14
0.000e+00 -2.220e-16 4.441e-16
See Also
solve(const BaseRMatrix&)const
Parameters
[in]mBrmatrix $B$.
[out]dErrNorm of computation error.
Returns
Result object.

Definition at line 16542 of file cvm.h.

template<typename TR>
BaseRMatrix basic_srmatrix< TR >::solve_tran ( const BaseRMatrix mB,
TR &  dErr 
) const throw (cvmexception)
inline

Linear solver (transposed)

Creates rmatrix object as solution $X$ of matrix linear equation $A^T*X=B$ (which is equivalent to $X^T*A=B^T$) where calling matrix is square matrix $A$ and parameter mB is matrix $B$. It also sets output parameter dErr to be equal to the norm of computation error. Function throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve_tran (mb, dErr);
std::cout << mx << std::endl << ~ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve_tran (vb, dErr);
std::cout << vx << std::endl << vx * ma - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-1.000e+000 -2.000e+000 -3.000e+000 -3.333e-001
+1.000e+000 +2.000e+000 +3.000e+000 +6.667e-001
+0.000e+000 +0.000e+000 -1.332e-015 +0.000e+000
+0.000e+000 +0.000e+000 +0.000e+000 +0.000e+000
+0.000e+000 +0.000e+000 +8.882e-016 -2.220e-016
+0.000e+000 +0.000e+000 -2.665e-015 +0.000e+000
-2.000e+000 +2.000e+000 +0.000e+000
+0.000e+000 +0.000e+000 +0.000e+000
See Also
solve_tran(const BaseRMatrix&)const
Parameters
[in]mBrmatrix $B$.
[out]dErrNorm of computation error.
Returns
Result object.

Definition at line 16605 of file cvm.h.

template<typename TR>
BaseRMatrix basic_srmatrix< TR >::solve ( const BaseRMatrix mB) const throw (cvmexception)
inline

Linear solver.

Creates rmatrix object as solution $X$ of matrix linear equation $A*X=B$ where calling square matrix is matrix $A$ and parameter mB is matrix $B$. Function throws cvmexception in case of inappropriate sizes of the operands or when the matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve (mb);
std::cout << mx << std::endl << ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve (vb);
std::cout << vx << std::endl << ma * vx - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-3.333e-01 -6.667e-01 -1.000e+00 1.000e+00
3.333e-01 6.667e-01 1.000e+00 0.000e+00
6.661e-16 1.332e-15 0.000e+00 0.000e+00
0.000e+00 0.000e+00 0.000e+00 0.000e+00
-1.110e-16 -2.220e-16 0.000e+00 0.000e+00
2.220e-16 4.441e-16 0.000e+00 0.000e+00
-6.667e-01 6.667e-01 1.332e-15
0.000e+00 -2.220e-16 4.441e-16
See Also
solve(const BaseRMatrix&,TR&)const
Parameters
[in]mBrmatrix $B$.
Returns
Result object.

Definition at line 16664 of file cvm.h.

template<typename TR>
BaseRMatrix basic_srmatrix< TR >::solve_tran ( const BaseRMatrix mB) const throw (cvmexception)
inline

Linear solver (transposed)

Creates rmatrix object as solution $X$ of matrix linear equation $A^T*X=B$ (which is equivalent to $X^T*A=B^T$) where calling matrix is square matrix $A$ and parameter mB is matrix $B$. Function throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix ma(a, 3);
rmatrix mb(3,4);
rmatrix mx(3,4);
mb(1).set(1.);
mb(2).set(2.);
mb(3).set(3.);
mb(1,4) = 1.; mb(2,4) = 2.; mb(3,4) = 3.;
mx = ma.solve_tran (mb, dErr);
std::cout << mx << std::endl << ~ma * mx - mb << std::endl;
rvector vb(3), vx(3);
vb = mb(2);
vx = ma.solve_tran (vb, dErr);
std::cout << vx << std::endl << vx * ma - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-1.000e+000 -2.000e+000 -3.000e+000 -3.333e-001
+1.000e+000 +2.000e+000 +3.000e+000 +6.667e-001
+0.000e+000 +0.000e+000 -1.332e-015 +0.000e+000
+0.000e+000 +0.000e+000 +0.000e+000 +0.000e+000
+0.000e+000 +0.000e+000 +8.882e-016 -2.220e-016
+0.000e+000 +0.000e+000 -2.665e-015 +0.000e+000
-2.000e+000 +2.000e+000 +0.000e+000
+0.000e+000 +0.000e+000 +0.000e+000
See Also
solve_tran(const BaseRMatrix&,TR&)const
Parameters
[in]mBrmatrix $B$.
Returns
Result object.

Definition at line 16727 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::operator% ( const RVector vB) const throw (cvmexception)
inline

Linear solver operator (transposed)

Creates rvector object as solution $x$ of linear equation $A^T*x=b$ (which is equivalent to $x*A=b$) where calling matrix is square matrix $A$ and parameter vB is vector $b$. Operator throws cvmexception in case of inappropriate sizes of the operands or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (12);
try {
double m[] = {1., -1., 1., 2., -2., 1., 3., -2., 1.};
double b[] = {1., 2., 3.};
srmatrix ma(m, 3);
rvector vb(b, 3);
rvector vx(3);
vx = ma % vb;
std::cout << vx * ma - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
+0.000000000000e+000 +0.000000000000e+000 +0.000000000000e+000
Parameters
[in]vBrvector $b$.
Returns
Result object.

Definition at line 16770 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::operator/ ( const RVector vB) const throw (cvmexception)
inline

Linear solver operator.

Creates rvector object as solution $x$ of linear equation $A*x=b$ where calling matrix is square matrix $A$ and parameter vB is vector $b$. Operator throws cvmexception in case of inappropriate sizes of the operands or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (12);
try {
double m[] = {1., -1., 1., 2., -2., 1., 3., -2., 1.};
double b[] = {1., 2., 3.};
srmatrix ma(m, 3);
rvector vb(b, 3);
rvector vx(3);
vx = ma / vb;
std::cout << ma * vx - vb;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
+0.000000000000e+000 +0.000000000000e+000 +0.000000000000e+000
Parameters
[in]vBrvector $b$.
Returns
Result object.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 16811 of file cvm.h.

template<typename TR>
RVector basic_srmatrix< TR >::solve_lu ( const basic_srmatrix< TR > &  mLU,
const tint pPivots,
const RVector vB,
TR &  dErr 
) const throw (cvmexception)
inline

LU factorization based linear solver.

Creates object of type rvector as solution $x$ of linear equation $A*x=b$ where calling matrix is square matrix $A$, parameter mLU is LU factorization (see low_up() ) of matrix $A$, parameter pPivots is an array of pivot numbers created while factorizing matrix $A$ and parameter vB is vector $b$. Function also sets output parameter dErr to be equal to the norm of computation error. This function is useful when you need to solve few linear equations of kind $A*x=b$ with the same matrix $A$ and different vectors $b$. In such case you save on matrix $A$ factorization since it's needed to be performed just one time. It also sets output parameter dErr to be equal to the norm of computation error and throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., -1., 1., 2., -2., 1., 3., -2., 1.};
srmatrix ma(a,3);
srmatrix mLU(3);
rmatrix mb1(3,2); rvector vb1(3);
rmatrix mb2(3,2); rvector vb2(3);
rmatrix mx1(3,2); rvector vx1(3);
rmatrix mx2(3,2); rvector vx2(3);
iarray nPivots(3);
double dErr = 0.;
mb1.randomize(-1.,3.); vb1.randomize(-2.,4.);
mb2.randomize(-2.,5.); vb2.randomize(-3.,1.);
mLU.low_up(ma, nPivots);
mx1 = ma.solve_lu (mLU, nPivots, mb1, dErr);
std::cout << mx1 << dErr << std::endl;
mx2 = ma.solve_lu (mLU, nPivots, mb2);
std::cout << mx2 << std::endl;;
std::cout << ma * mx1 - mb1 << std::endl << ma * mx2 - mb2;
vx1 = ma.solve_lu (mLU, nPivots, vb1, dErr);
std::cout << vx1 << dErr << std::endl;
vx2 = ma.solve_lu (mLU, nPivots, vb2);
std::cout << vx2 << std::endl;;
std::cout << ma * vx1 - vb1 << std::endl << ma * vx2 - vb2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
2.807e+00 1.107e+00
-3.651e-01 -4.843e+00
-5.412e-01 3.095e+00
6.438e-15
-7.639e-01 1.082e+01
-2.869e-01 -1.110e+01
4.890e-01 3.443e+00
0.000e+00 -4.441e-16
1.110e-16 -4.441e-16
-4.441e-16 4.441e-16
0.000e+00 -4.441e-16
0.000e+00 8.882e-16
0.000e+00 -4.441e-16
-1.651e+00 2.361e-01 -6.384e-02
3.828e-15
-5.886e+00 7.038e+00 -3.125e+00
0.000e+00 0.000e+00 0.000e+00
0.000e+00 0.000e+00 2.220e-16
See Also
low_up()
Parameters
[in]mLULU factorization of matrix $A$.
[in]pPivotspivots vector.
[in]vBrvector $b$.
[out]dErrNorm of computation error.
Returns
Result object.

Definition at line 16899 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
RVector basic_srmatrix< TR >::solve_lu ( const basic_srmatrix< TR > &  mLU,
const tint pPivots,
const RVector vB 
) const throw (cvmexception)
inline

LU factorization based linear solver.

Creates object of type rvector as solution $x$ of linear equation $A*x=b$ where calling matrix is square matrix $A$, parameter mLU is LU factorization (see low_up() ) of matrix $A$, parameter pPivots is an array of pivot numbers created while factorizing matrix $A$ and parameter vB is vector $b$. Function also sets output parameter dErr to be equal to the norm of computation error. This function is useful when you need to solve few linear equations of kind $A*x=b$ with the same matrix $A$ and different vectors $b$. In such case you save on matrix $A$ factorization since it's needed to be performed just one time. It throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., -1., 1., 2., -2., 1., 3., -2., 1.};
srmatrix ma(a,3);
srmatrix mLU(3);
rmatrix mb1(3,2); rvector vb1(3);
rmatrix mb2(3,2); rvector vb2(3);
rmatrix mx1(3,2); rvector vx1(3);
rmatrix mx2(3,2); rvector vx2(3);
iarray nPivots(3);
double dErr = 0.;
mb1.randomize(-1.,3.); vb1.randomize(-2.,4.);
mb2.randomize(-2.,5.); vb2.randomize(-3.,1.);
mLU.low_up(ma, nPivots);
mx1 = ma.solve_lu (mLU, nPivots, mb1, dErr);
std::cout << mx1 << dErr << std::endl;
mx2 = ma.solve_lu (mLU, nPivots, mb2);
std::cout << mx2 << std::endl;;
std::cout << ma * mx1 - mb1 << std::endl << ma * mx2 - mb2;
vx1 = ma.solve_lu (mLU, nPivots, vb1, dErr);
std::cout << vx1 << dErr << std::endl;
vx2 = ma.solve_lu (mLU, nPivots, vb2);
std::cout << vx2 << std::endl;;
std::cout << ma * vx1 - vb1 << std::endl << ma * vx2 - vb2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
2.807e+00 1.107e+00
-3.651e-01 -4.843e+00
-5.412e-01 3.095e+00
6.438e-15
-7.639e-01 1.082e+01
-2.869e-01 -1.110e+01
4.890e-01 3.443e+00
0.000e+00 -4.441e-16
1.110e-16 -4.441e-16
-4.441e-16 4.441e-16
0.000e+00 -4.441e-16
0.000e+00 8.882e-16
0.000e+00 -4.441e-16
-1.651e+00 2.361e-01 -6.384e-02
3.828e-15
-5.886e+00 7.038e+00 -3.125e+00
0.000e+00 0.000e+00 0.000e+00
0.000e+00 0.000e+00 2.220e-16
See Also
low_up()
Parameters
[in]mLULU factorization of matrix $A$.
[in]pPivotspivots vector.
[in]vBrvector $b$.
Returns
Result object.

Definition at line 16990 of file cvm.h.

template<typename TR>
BaseRMatrix basic_srmatrix< TR >::solve_lu ( const basic_srmatrix< TR > &  mLU,
const tint pPivots,
const BaseRMatrix mB,
TR &  dErr 
) const throw (cvmexception)
inline

LU factorization based linear solver.

Creates object of type rmatrix as solution $X$ of matrix linear equation $A*X=B$ where calling matrix is square matrix $A$, parameter mLU is LU factorization (see low_up() ) of matrix $A$, parameter pPivots is an array of pivot numbers created while factorizing matrix $A$ and parameter mB is matrix $B$. This function is useful when you need to solve few linear equations of kind $A*X=B$ with the same matrix $A$ and different matrices $B$. In such case you save on matrix $A$ factorization since it's needed to be performed just one time. Function also sets output parameter dErr to be equal to the norm of computation error and throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., -1., 1., 2., -2., 1., 3., -2., 1.};
srmatrix ma(a,3);
srmatrix mLU(3);
rmatrix mb1(3,2); rvector vb1(3);
rmatrix mb2(3,2); rvector vb2(3);
rmatrix mx1(3,2); rvector vx1(3);
rmatrix mx2(3,2); rvector vx2(3);
iarray nPivots(3);
double dErr = 0.;
mb1.randomize(-1.,3.); vb1.randomize(-2.,4.);
mb2.randomize(-2.,5.); vb2.randomize(-3.,1.);
mLU.low_up(ma, nPivots);
mx1 = ma.solve_lu (mLU, nPivots, mb1, dErr);
std::cout << mx1 << dErr << std::endl;
mx2 = ma.solve_lu (mLU, nPivots, mb2);
std::cout << mx2 << std::endl;;
std::cout << ma * mx1 - mb1 << std::endl << ma * mx2 - mb2;
vx1 = ma.solve_lu (mLU, nPivots, vb1, dErr);
std::cout << vx1 << dErr << std::endl;
vx2 = ma.solve_lu (mLU, nPivots, vb2);
std::cout << vx2 << std::endl;;
std::cout << ma * vx1 - vb1 << std::endl << ma * vx2 - vb2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
2.807e+00 1.107e+00
-3.651e-01 -4.843e+00
-5.412e-01 3.095e+00
6.438e-15
-7.639e-01 1.082e+01
-2.869e-01 -1.110e+01
4.890e-01 3.443e+00
0.000e+00 -4.441e-16
1.110e-16 -4.441e-16
-4.441e-16 4.441e-16
0.000e+00 -4.441e-16
0.000e+00 8.882e-16
0.000e+00 -4.441e-16
-1.651e+00 2.361e-01 -6.384e-02
3.828e-15
-5.886e+00 7.038e+00 -3.125e+00
0.000e+00 0.000e+00 0.000e+00
0.000e+00 0.000e+00 2.220e-16
See Also
low_up()
Parameters
[in]mLULU factorization of matrix $A$.
[in]pPivotspivots vector.
[in]mBrmatrix $B$.
[out]dErrNorm of computation error.
Returns
Result object.

Definition at line 17079 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
BaseRMatrix basic_srmatrix< TR >::solve_lu ( const basic_srmatrix< TR > &  mLU,
const tint pPivots,
const BaseRMatrix mB 
) const throw (cvmexception)
inline

LU factorization based linear solver.

Creates object of type rmatrix as solution $X$ of matrix linear equation $A*X=B$ where calling matrix is square matrix $A$, parameter mLU is LU factorization (see low_up() ) of matrix $A$, parameter pPivots is an array of pivot numbers created while factorizing matrix $A$ and parameter mB is matrix $B$. This function is useful when you need to solve few linear equations of kind $A*X=B$ with the same matrix $A$ and different matrices $B$. In such case you save on matrix $A$ factorization since it's needed to be performed just one time. Function throws cvmexception in case of inappropriate sizes of the objects or when matrix $A$ is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., -1., 1., 2., -2., 1., 3., -2., 1.};
srmatrix ma(a,3);
srmatrix mLU(3);
rmatrix mb1(3,2); rvector vb1(3);
rmatrix mb2(3,2); rvector vb2(3);
rmatrix mx1(3,2); rvector vx1(3);
rmatrix mx2(3,2); rvector vx2(3);
iarray nPivots(3);
double dErr = 0.;
mb1.randomize(-1.,3.); vb1.randomize(-2.,4.);
mb2.randomize(-2.,5.); vb2.randomize(-3.,1.);
mLU.low_up(ma, nPivots);
mx1 = ma.solve_lu (mLU, nPivots, mb1, dErr);
std::cout << mx1 << dErr << std::endl;
mx2 = ma.solve_lu (mLU, nPivots, mb2);
std::cout << mx2 << std::endl;;
std::cout << ma * mx1 - mb1 << std::endl << ma * mx2 - mb2;
vx1 = ma.solve_lu (mLU, nPivots, vb1, dErr);
std::cout << vx1 << dErr << std::endl;
vx2 = ma.solve_lu (mLU, nPivots, vb2);
std::cout << vx2 << std::endl;;
std::cout << ma * vx1 - vb1 << std::endl << ma * vx2 - vb2;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
2.807e+00 1.107e+00
-3.651e-01 -4.843e+00
-5.412e-01 3.095e+00
6.438e-15
-7.639e-01 1.082e+01
-2.869e-01 -1.110e+01
4.890e-01 3.443e+00
0.000e+00 -4.441e-16
1.110e-16 -4.441e-16
-4.441e-16 4.441e-16
0.000e+00 -4.441e-16
0.000e+00 8.882e-16
0.000e+00 -4.441e-16
-1.651e+00 2.361e-01 -6.384e-02
3.828e-15
-5.886e+00 7.038e+00 -3.125e+00
0.000e+00 0.000e+00 0.000e+00
0.000e+00 0.000e+00 2.220e-16
See Also
low_up()
Parameters
[in]mLULU factorization of matrix $A$.
[in]pPivotspivots vector.
[in]mBrmatrix $B$.
Returns
Result object.

Definition at line 17169 of file cvm.h.

template<typename TR>
TR basic_srmatrix< TR >::det ( ) const throw (cvmexception)
inline

Matrix determinant.

Returns determinant of calling matrix. It uses the LU factorization internally and may throw the same exceptions as the factorizer.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix m(a, 3);
std::cout << m << std::endl << m.det() << std::endl;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1.000e+00 4.000e+00 7.000e+00
2.000e+00 5.000e+00 8.000e+00
3.000e+00 6.000e+00 1.000e+01
-3.000e+00
See Also
low_up()
Returns
Determinant value.

Definition at line 17205 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::low_up ( const basic_srmatrix< TR > &  m,
tint nPivots 
) throw (cvmexception)
inline

Low-up (LU) factorization.

Compute LU factorization of square matrix $A$ as

\[ A=PLU \]

where $P$ is permutation matrix, $L$ is lower triangular matrix with unit diagonal elements and $U$ is upper triangular matrix. Function stores result as the matrix $L$ without main diagonal combined with $U$. Function returns pivot indices as array of integers (it should support at least msize() elements) pointed to by nPivots so $i$-th row was interchanged with nPivots[ $i$]-th row. This version sets calling matrix to be equal to matrix m's LU factorization. Function throws cvmexception in case of inappropriate sizes of the operands or when the matrix to be factorized is close to singular. It is recommended to use iarray for pivot values.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix m(a, 3);
srmatrix mLU(3), mLo(3), mUp(3);
iarray naPivots(3);
mLU.low_up (m, naPivots);
mLo.identity ();
mLo.diag(-2) = mLU.diag(-2);
mLo.diag(-1) = mLU.diag(-1);
mUp.diag(0) = mLU.diag(0);
mUp.diag(1) = mLU.diag(1);
mUp.diag(2) = mLU.diag(2);
std::cout << mLo << std::endl << mUp
<< std::endl << naPivots << std::endl;
mLU = mLo * mUp;
for (int i = 3; i >= 1; i--) {
mLU.swap_rows (i, naPivots[i]);
}
std::cout << mLU;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
1.000e+00 0.000e+00 0.000e+00
3.333e-01 1.000e+00 0.000e+00
6.667e-01 5.000e-01 1.000e+00
3.000e+00 6.000e+00 1.000e+01
0.000e+00 2.000e+00 3.667e+00
0.000e+00 0.000e+00 -5.000e-01
3 3 3
1.000e+00 4.000e+00 7.000e+00
2.000e+00 5.000e+00 8.000e+00
3.000e+00 6.000e+00 1.000e+01
Parameters
[in]msrmatrix to compute LU factorization for.
[out]nPivotsArray of pivot indices.
Returns
Reference to changed calling matrix.

Definition at line 17284 of file cvm.h.

Here is the caller graph for this function:

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::low_up ( tint nPivots) const throw (cvmexception)
inline

Low-up (LU) factorization.

Compute LU factorization of square matrix $A$ as

\[ A=PLU \]

where $P$ is permutation matrix, $L$ is lower triangular matrix with unit diagonal elements and $U$ is upper triangular matrix. Function stores result as the matrix $L$ without main diagonal combined with $U$. Function returns pivot indices as array of integers (it should support at least msize() elements) pointed to by nPivots so $i$-th row was interchanged with nPivots[ $i$]-th row. This version creates object of type srmatrix as calling matrix LU factorization. Function throws cvmexception in case of inappropriate sizes of the operands or when the matrix to be factorized is close to singular. It is recommended to use iarray for pivot values.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix m(a, 3);
srmatrix mLU(3), mLo(3), mUp(3);
iarray naPivots(3);
mLU = m.low_up (naPivots);
mLo.identity ();
mLo.diag(-2) = mLU.diag(-2);
mLo.diag(-1) = mLU.diag(-1);
mUp.diag(0) = mLU.diag(0);
mUp.diag(1) = mLU.diag(1);
mUp.diag(2) = mLU.diag(2);
std::cout << mLo << std::endl << mUp
<< std::endl << naPivots << std::endl;
mLU = mLo * mUp;
for (int i = 3; i >= 1; i--) {
mLU.swap_rows (i, naPivots[i]);
}
std::cout << mLU;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
+1.000e+000 +0.000e+000 +0.000e+000
+3.333e-001 +1.000e+000 +0.000e+000
+6.667e-001 +5.000e-001 +1.000e+000
+3.000e+000 +6.000e+000 +1.000e+001
+0.000e+000 +2.000e+000 +3.667e+000
+0.000e+000 +0.000e+000 -5.000e-001
+3 +3 +3
+1.000e+000 +4.000e+000 +7.000e+000
+2.000e+000 +5.000e+000 +8.000e+000
+3.000e+000 +6.000e+000 +1.000e+001
Parameters
[out]nPivotsArray of pivot indices.
Returns
Result object.

Reimplemented in basic_srbmatrix< TR >.

Definition at line 17361 of file cvm.h.

template<typename TR>
TR basic_srmatrix< TR >::cond ( ) const throw (cvmexception)
inline

Condition number reciprocal.

Returns condition number reciprocal of calling matrix $A$ in the infinity-norm

\[ \kappa_\infty=\|A\|_\infty \|A^{-1}\|_\infty. \]

Less value returned means that matrix $A$ is closer to singular. Zero value returned means estimation underflow or that matrix $A$ is singular. The condition number is used for error analysis of systems of linear equations. Function throws cvmexception in case of LAPACK subroutines failure.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
srmatrix m(a, 3);
std::cout << m.cond() << std::endl
<< m.det() << std::endl << std::endl;
m(3,3) = 10.;
std::cout << m.cond() << std::endl << m.det() << std::endl;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
0.000e+00
0.000e+00
7.519e-03
-3.000e+00
See Also
basic_array::norminf()
Returns
Result value.

Definition at line 17409 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::inv ( const basic_srmatrix< TR > &  m) throw (cvmexception)
inline

Matrix inversion.

This version sets calling matrix to be equal to square matrix m inverted. Function throws cvmexception in case of inappropriate sizes of the operands or when the matrix to be inverted is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (10);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix m(a, 3);
srmatrix mi(3);
mi.inv (m);
std::cout << mi << std::endl << mi * m - eye_real(3);
std::cout << std::endl << mi.inv() * mi - eye_real(3);
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-6.6666666667e-01 -6.6666666667e-01 1.0000000000e+00
-1.3333333333e+00 3.6666666667e+00 -2.0000000000e+00
1.0000000000e+00 -2.0000000000e+00 1.0000000000e+00
0.0000000000e+00 0.0000000000e+00 1.7763568394e-15
1.7763568394e-15 3.5527136788e-15 0.0000000000e+00
0.0000000000e+00 0.0000000000e+00 1.7763568394e-15
0.0000000000e+00 1.7763568394e-15 -1.7763568394e-15
-8.8817841970e-16 3.5527136788e-15 -3.5527136788e-15
0.0000000000e+00 0.0000000000e+00 -1.7763568394e-15
Parameters
[in]msrmatrix to invert.
Returns
Reference to changed calling matrix.

Definition at line 17457 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::inv ( ) const throw (cvmexception)
inline

Matrix inversion.

This version creates srmatrix object equal to calling matrix inverted. Function throws cvmexception in case of memory allocation failure or when the matrix to be inverted is close to singular.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (10);
try {
double a[] = {1., 2., 3., 4., 5., 6., 7., 8., 10.};
srmatrix m(a, 3);
srmatrix mi(3);
mi.inv (m);
std::cout << mi << std::endl << mi * m - eye_real(3);
std::cout << std::endl << mi.inv() * mi - eye_real(3);
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-6.6666666667e-01 -6.6666666667e-01 1.0000000000e+00
-1.3333333333e+00 3.6666666667e+00 -2.0000000000e+00
1.0000000000e+00 -2.0000000000e+00 1.0000000000e+00
0.0000000000e+00 0.0000000000e+00 1.7763568394e-15
1.7763568394e-15 3.5527136788e-15 0.0000000000e+00
0.0000000000e+00 0.0000000000e+00 1.7763568394e-15
0.0000000000e+00 1.7763568394e-15 -1.7763568394e-15
-8.8817841970e-16 3.5527136788e-15 -3.5527136788e-15
0.0000000000e+00 0.0000000000e+00 -1.7763568394e-15
Returns
Result object.

Reimplemented in basic_srsmatrix< TR >.

Definition at line 17503 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::exp ( const basic_srmatrix< TR > &  m,
TR  tol = basic_cvmMachSp<TR>() 
) throw (cvmexception)
inline

Matrix exponent.

Computes exponent of square matrix using Pade approximation defined as

\[ R_{pq}(z)=D_{pq}(z)^{-1}N_{pq}(z)=1+z+\dots+z^p/p!\,, \]

where

\[\begin{aligned} N_{pq}(z)&=\sum_{k=0}^p\frac{(p+q-k)!p!}{(p+q)!k!(q-k)!}z^k,\\ D_{pq}(z)&=\sum_{k=0}^q\frac{(p+q-k)!p!}{(p+q)!k!(q-k)!}(-z)^k \end{aligned}\]

along with the matrix normalizing as described in Gene H. Golub, Charles F. Van Loan. Matrix Computations, The Johns Hopkins University Press, 1996, 694~p., ISBN 0-8018-5413-X, p. 572. Function uses DMEXP (or SMEXP for float version) FORTRAN subroutine implementing the algorithm. This version sets calling matrix to be equal to exponent of m and returns reference to the matrix changed. The algorithm uses parameter tol as $\varepsilon(p,q)$ in order to choose constants $p$ and $q$ so that

\[ \varepsilon(p,q)\ge 2^{3-(p+q)}\frac{p!q!}{(p+q)!(p+q+1)!}. \]

This parameter is equal to the largest relative spacing (see cvmMachSp() ) by default. Function throws cvmexception in case of inappropriate sizes of the operands or when LAPACK subroutine fails.

Example
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (15);
try {
srmatrix m(2), me(2);
m(1,1) = -49.;
m(1,2) = 24.;
m(2,1) = -64.;
m(2,2) = 31.;
me.exp(m);
std::cout << m << std::endl << me;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-4.900000000000000e+001 2.400000000000000e+001
-6.400000000000000e+001 3.100000000000000e+001
-7.357587581448284e-001 5.518190996581556e-001
-1.471517599088415e+000 1.103638240715692e+000
Matlab output:
-7.357587581446907e-001 5.518190996580505e-001
-1.471517599088136e+000 1.103638240715478e+000
Parameters
[in]msrmatrix to compute exponent for.
[in]tolComputation tolerance.
Returns
Reference to changed calling matrix.

Definition at line 17571 of file cvm.h.

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::exp ( TR  tol = basic_cvmMachSp<TR>()) const throw (cvmexception)
inline

Matrix exponent.

Computes exponent of square matrix using Pade approximation defined as

\[ R_{pq}(z)=D_{pq}(z)^{-1}N_{pq}(z)=1+z+\dots+z^p/p!\,, \]

where

\[\begin{aligned} N_{pq}(z)&=\sum_{k=0}^p\frac{(p+q-k)!p!}{(p+q)!k!(q-k)!}z^k,\\ D_{pq}(z)&=\sum_{k=0}^q\frac{(p+q-k)!p!}{(p+q)!k!(q-k)!}(-z)^k \end{aligned}\]

along with the matrix normalizing as described in Gene H. Golub, Charles F. Van Loan. Matrix Computations, The Johns Hopkins University Press, 1996, 694~p., ISBN 0-8018-5413-X, p. 572. Function uses DMEXP (or SMEXP for float version) FORTRAN subroutine implementing the algorithm. This version creates object of type srmatrix as exponent of calling matrix. The algorithm uses parameter tol as $\varepsilon(p,q)$ in order to choose constants $p$ and $q$ so that

\[ \varepsilon(p,q)\ge 2^{3-(p+q)}\frac{p!q!}{(p+q)!(p+q+1)!}. \]

This parameter is equal to the largest relative spacing (see cvmMachSp() ) by default. Function throws cvmexception in case of inappropriate sizes of the operands or when LAPACK subroutine fails.

Example
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (15);
try {
srmatrix m(2);
m(1,1) = -49.;
m(1,2) = 24.;
m(2,1) = -64.;
m(2,2) = 31.;
std::cout << m << std::endl << m.exp();
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-4.900000000000000e+01 2.400000000000000e+01
-6.400000000000000e+01 3.100000000000000e+01
-7.357587581448284e-01 5.518190996581556e-01
-1.471517599088415e+00 1.103638240715692e+00
Matlab output:
-7.357587581446907e-001 5.518190996580505e-001
-1.471517599088136e+000 1.103638240715478e+000
Parameters
[in]tolComputation tolerance.
Returns
Result object.

Reimplemented in basic_srsmatrix< TR >.

Definition at line 17635 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::polynom ( const basic_srmatrix< TR > &  m,
const RVector v 
) throw (cvmexception)
inline

Matrix polynomial.

Computes matrix polynomial defined as

\[ p(A)=b_0I+b_1A+\dots+b_qA^q \]

using the Horner's rule:

\[ \newcommand{\floor}{\operatorname{floor}}p(A)=\sum_{k=0}^r B_k(A^s)^k,\quad s=\floor(\!\sqrt{q}\,),\ r=\floor(q/s) \]

where

\[ B_k=\begin{cases} \sum\limits_{i=0}^{s-1}b_{sk+i} A^i, & k=0,1,\dots,r-1\\ \sum\limits_{i=0}^{q-sr}b_{sr+i} A^i, & k=r. \end{cases} \]

See also Gene H. Golub, Charles F. Van Loan. Matrix Computations, The Johns Hopkins University Press, 1996, 694~p., ISBN 0-8018-5413-X, p. 568. The coefficients $b_0,b_1,\dots,b_q$ are passed in parameter v, where $q$ is equal to v.size()-1, so the function computes matrix polynomial equal to

\[ \mathtt{v[1]*I + v[2]*m +\cdots + v[v.size()]*m^{v.size()-1}} \]

This version sets calling matrix to be equal to the polynomial of m. Function uses DPOLY (or SPOLY for float version) FORTRAN subroutine implementing the Horner's algorithm. Function throws cvmexception in case of inappropriate sizes of the operands.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (15);
try {
double a[] = {2.2, 1.3, 1.1, -0.9, 0.2,
-0.45, 45, -30, 10, 3, 3.2};
const rvector v(a, 11);
srmatrix m(2), mp(2);
m(1,1) = 1.;
m(1,2) = 0.5;
m(2,1) = -1.;
m(2,2) = 0.3;
mp.polynom (m, v);
std::cout << mp;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-7.963641665999998e+00 -7.551532476200001e+00
1.510306495240000e+01 2.608503800680002e+00
Matlab output:
-7.963641665999999e+000 -7.551532476200002e+000
1.510306495240000e+001 2.608503800680002e+000
Parameters
[in]msrmatrix to compute polynomial for.
[in]vVector of coefficients.
Returns
Reference to changed calling matrix.

Definition at line 17706 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix basic_srmatrix< TR >::polynom ( const RVector v) const
inline

Matrix polynomial.

Computes matrix polynomial defined as

\[ p(A)=b_0I+b_1A+\dots+b_qA^q \]

using the Horner's rule:

\[ \newcommand{\floor}{\operatorname{floor}}p(A)=\sum_{k=0}^r B_k(A^s)^k,\quad s=\floor(\!\sqrt{q}\,),\ r=\floor(q/s) \]

where

\[ B_k=\begin{cases} \sum\limits_{i=0}^{s-1}b_{sk+i} A^i, & k=0,1,\dots,r-1\\ \sum\limits_{i=0}^{q-sr}b_{sr+i} A^i, & k=r. \end{cases} \]

See also Gene H. Golub, Charles F. Van Loan. Matrix Computations, The Johns Hopkins University Press, 1996, 694~p., ISBN 0-8018-5413-X, p. 568. The coefficients $b_0,b_1,\dots,b_q$ are passed in parameter v, where $q$ is equal to v.size()-1, so the function computes matrix polynomial equal to

\[ \mathtt{v[1]*I + v[2]*m +\cdots + v[v.size()]*m^{v.size()-1}} \]

This version creates object of type srmatrix as the polynomial of calling matrix. Function uses DPOLY (or SPOLY for float version) FORTRAN subroutine implementing the Horner's algorithm. Function throws cvmexception in case of memory allocation failure.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (15);
try {
double a[] = {2.2, 1.3, 1.1, -0.9, 0.2,
-0.45, 45, -30, 10, 3, 3.2};
const rvector v(a, 11);
srmatrix m(2);
m(1,1) = 1.;
m(1,2) = 0.5;
m(2,1) = -1.;
m(2,2) = 0.3;
srmatrix mp = m.polynom (v);
std::cout << mp;
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
-7.963641665999998e+000 -7.551532476200001e+000
1.510306495240000e+001 2.608503800680002e+000
Matlab output:
-7.963641665999999e+000 -7.551532476200002e+000
1.510306495240000e+001 2.608503800680002e+000
Parameters
[in]vVector of coefficients.
Returns
Result object.

Reimplemented in basic_srsmatrix< TR >.

Definition at line 17778 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
CVector basic_srmatrix< TR >::eig ( basic_scmatrix< TR, TC > &  mEigVect,
bool  bRightVect = true 
) const throw (cvmexception)
inline

Eigenvalues and eigenvectors.

Solves eigenvalue problem and creates cvector object equal to eigenvalues of calling matrix. The eigenvalue problem is defined as follows: given square matrix $A$, find the eigenvalues $\lambda$ and the corresponding eigenvectors $z$ that satisfy the equation

\[ Az = \lambda z. \]

Some eigenvalues may be complex even for real matrix $A$. Moreover, if real nonsymmetric matrix has a complex eigenvalue $a+bi$ corresponding to an eigenvector $z$, then $a-bi$ is also an eigenvalue. The eigenvalue $a-bi$ corresponds to the eigenvector whose elements are complex conjugate to the elements of $z$. Function sets output parameter mEigVect to be equal to square matrix containing eigenvectors as columns. Function also computes "left" eigenvectors if parameter bRightVect is set to false. Left eigencectors satisfy

\[ z^HA = \lambda z^H. \]

Function throws cvmexception in case of inappropriate calling object sizes or in case of convergence error.

Example:
using namespace cvm;
try {
scmatrix m(3), me(3);
cvector vl(3);
m(1,1) = 0.1; m(1,2) = 0.2; m(1,3) = 0.1;
m(2,1) = 0.11; m(2,2) = -2.9; m(2,3) = -8.4;
m(3,1) = 0.; m(3,2) = 2.91; m(3,3) = 8.2;
vl = m.eig (me);
std::cout << vl;
m(2,2) = 2.9;
vl = m.eig (me);
std::cout << vl << std::endl;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (1);
std::cout << m * me(1) - me(1) * vl(1);
std::cout << m * me(2) - me(2) * vl(2);
std::cout << m * me(3) - me(3) * vl(3);
}
catch (std::exception& e) {
std::cout << "Exception: " << e.what () << std::endl;
}
prints
(-0.0555784,0) (0.285327,0) (5.17025,0)
(0.0968985,-1.38778e-017) (5.55155,4.1733) (5.55155,-4.1733)
(+1.4e-017,+2.8e-017) (-1.0e-016,-5.4e-017) (+3.9e-017,+4.3e-017)
(-5.0e-016,+1.7e-016) (-7.1e-015,+2.9e-015) (+3.6e-015,+2.7e-015)
(-3.1e-016,+2.8e-017) (-2.7e-015,-8.9e-016) (+8.9e-016,+3.6e-015)
See Also
cvector::eig()
Parameters
[out]mEigVectEigenvectors of calling matrix.
[in]bRightVecttrue (default) to compute right eigenvectors.
Returns
Result object.

Definition at line 17853 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
CVector basic_srmatrix< TR >::eig ( ) const throw (cvmexception)
inline

Eigenvalues.

Solves eigenvalue problem and creates cvector object equal to eigenvalues of calling matrix. The eigenvalue problem is defined as follows: given square matrix $A$, find the eigenvalues $\lambda$ and the corresponding eigenvectors $z$ that satisfy the equation

\[ Az = \lambda z. \]

Some eigenvalues may be complex even for real matrix $A$. Moreover, if real nonsymmetric matrix has a complex eigenvalue $a+bi$ corresponding to an eigenvector $z$, then $a-bi$ is also an eigenvalue. The eigenvalue $a-bi$ corresponds to the eigenvector whose elements are complex conjugate to the elements of $z$. Function throws cvmexception in case of in caso of memory allocation failure or convergence error.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::showpos);
std::cout.precision (2);
try {
srmatrix m(3);
m(1,1) = 0.1; m(1,2) = 0.2; m(1,3) = 0.1;
m(2,1) = 0.11; m(2,2) = 2.9; m(2,3) = -8.4;
m(3,1) = 0.; m(3,2) = 2.91; m(3,3) = 8.2;
std::cout << m.eig();
}
catch (std::exception& e) {
std::cout << "Exception " << e.what () << std::endl;
}
prints
(+9.69e-002,+0.00e+000) (+5.55e+000,+4.17e+000) (+5.55e+000,-4.17e+000)
See Also
cvector::eig()
Returns
Result object.

Reimplemented in basic_srsmatrix< TR >.

Definition at line 17899 of file cvm.h.

Here is the call graph for this function:

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::cholesky ( const basic_srsmatrix< TR > &  m) throw (cvmexception)
inline

Cholesky factorization.

Forms the Cholesky factorization of symmetric positive-definite matrix $A$ defined as

\[ A=U^T U, \]

where $U$ is upper triangular matrix. It utilizes one of ZPOTRF routines of the LAPACK library. Function sets calling matrix to be equal to the factorization of symmetric positive-definite matrix m. Function throws cvmexception in case of inappropriate sizes of the operands or in case of convergence error.

Example:
using namespace cvm;
try {
double a[] = {1., 2., 1., 2., 5., -1., 1., -1., 20.};
const srsmatrix m(a, 3);
srmatrix h(3);
h.cholesky(m);
std::cout << h << std::endl;
std::cout << ~h * h - m;
}
catch (std::exception& e) {
std::cout << "Exception: " << e.what () << std::endl;
}
prints
1 2 1
0 1 -3
0 0 3.16228
0 0 0
0 0 0
0 0 0
See Also
http://www.netlib.org/lapack
Parameters
[in]msrsmatrix to factorize.
Returns
Reference to changed calling matrix.

Definition at line 17947 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::bunch_kaufman ( const basic_srsmatrix< TR > &  m,
tint nPivots 
) throw (cvmexception)
inline

Bunch-Kaufman factorization.

Forms the Bunch-Kaufman factorization of symmetric matrix (cited from the MKL library documentation):

\[ A=PUDU^TP^T, \]

where $A$ is the input matrix passed in parameter m, $P$ is a permutation matrix, $U$ and $L$ are upper and lower triangular matrices with unit diagonal, and $D$ is a symmetric block-diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks. $U$ and $L$ have 2-by-2 unit diagonal blocks corresponding to the 2-by-2 blocks of $D$. It utilizes one of DSYTRF routines of the LAPACK library. Function sets calling matrix to be equal to the factorization of symmetric positive-definite matrix m. Function throws cvmexception in case of inappropriate sizes of the operands or in case of convergence error. Function is mostly designed to be used for subsequent calls of DSYTRS, DSYCON and DSYTRI routines of the LAPACK library. Currently it's used internally in det() flow when argument is symmetric but not positive-definite.

See Also
http://www.netlib.org/lapack
Parameters
[in]msrsmatrix to factorize.
[out]nPivotsPivot indices array.
Returns
Reference to changed calling matrix.

Definition at line 17987 of file cvm.h.

template<typename TR>
void basic_srmatrix< TR >::qr ( basic_srmatrix< TR > &  mQ,
basic_srmatrix< TR > &  mR 
) const throw (cvmexception)
inline

QR factorization.

Computes QR factorization as

\[ A=QR \]

where $A$ is calling square matrix, orthogonal matrix $Q$ and upper triangular matrix $R$ are mQ and mR respectively. Function throws cvmexception in case if inappropriate sizes of the operands passed.

Example:
using namespace cvm;
treal a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
const cvm::srmatrix m(a, 3);
cvm::srmatrix q(3), r(3);
m.qr(q,r);
std::cout << (eye_real(3) - ~q * q).norm()
<< " " << (m - q * r).norm() << std::endl;
prints
+5.2889959e-16 +7.0854500e-15
Parameters
[out]mQsrmatrix $Q$.
[out]mRsrmatrix $R$.

Definition at line 18022 of file cvm.h.

template<typename TR>
void basic_srmatrix< TR >::lq ( basic_srmatrix< TR > &  mL,
basic_srmatrix< TR > &  mQ 
) const throw (cvmexception)
inline

LQ factorization.

Computes LQ factorization as

\[ A=LQ \]

where $A$ is calling square matrix, lower triangular matrix $L$ and orthogonal matrix $Q$ are mL and mQ respectively. Function throws cvmexception in case if inappropriate sizes of the operands passed.

Example:
using namespace cvm;
treal a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
const cvm::srmatrix m(a, 3);
cvm::srmatrix l(3), q(3);
m.lq(l,q);
std::cout << (eye_real(3) - ~q * q).norm()
<< " " << (m - l * q).norm() << std::endl;
prints
+7.3329369e-016 +8.1523942e-015
Parameters
[out]mLsrmatrix $L$.
[out]mQsrmatrix $Q$.

Definition at line 18053 of file cvm.h.

template<typename TR>
void basic_srmatrix< TR >::ql ( basic_srmatrix< TR > &  mQ,
basic_srmatrix< TR > &  mL 
) const throw (cvmexception)
inline

QL factorization.

Computes QL factorization as

\[ A=QL \]

where $A$ is calling square matrix, orthogonal matrix $Q$ and lower triangular matrix $L$ are mQ and mL respectively. Function throws cvmexception in case if inappropriate sizes of the operands passed.

Example:
using namespace cvm;
treal a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
const cvm::srmatrix m(a, 3);
cvm::srmatrix q(3), l(3);
m.ql(q,l);
std::cout << (eye_real(3) - ~q * q).norm()
<< " " << (m - q * l).norm() << std::endl;
prints
+1.6146017e-015 +4.3341378e-015
Parameters
[out]mQsrmatrix $Q$.
[out]mLsrmatrix $L$.

Definition at line 18084 of file cvm.h.

template<typename TR>
void basic_srmatrix< TR >::rq ( basic_srmatrix< TR > &  mR,
basic_srmatrix< TR > &  mQ 
) const throw (cvmexception)
inline

RQ factorization.

Computes RQ factorization as

\[ A=RQ \]

where $A$ is calling square matrix, upper triangular matrix $R$ and orthogonal matrix $Q$ are mR and mQ respectively. Function throws cvmexception in case if inappropriate sizes of the operands passed.

Example:
using namespace cvm;
treal a[] = {1., 2., 3., 4., 5., 6., 7., 8., 9.};
const cvm::srmatrix m(a, 3);
cvm::srmatrix r(3), q(3);
m.rq(r,q);
std::cout << (eye_real(3) - ~q * q).norm()
<< " " << (m - r * q).norm() << std::endl;
prints
+3.7030557e-016 +2.2752801e-015
Parameters
[out]mRsrmatrix $R$.
[out]mQsrmatrix $Q$.

Definition at line 18115 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::identity ( )
inline

Identity matrix.

Sets calling matrix to be equal to identity matrix and returns reference to the matrix changed.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
srmatrix m(3);
m.randomize(0.,1.);
std::cout << m << std::endl;
std::cout << m.identity();
prints
9.423e-01 2.950e-01 8.429e-01
2.013e-01 3.250e-01 2.904e-01
7.920e-01 2.405e-02 7.801e-01
1.000e+00 0.000e+00 0.000e+00
0.000e+00 1.000e+00 0.000e+00
0.000e+00 0.000e+00 1.000e+00
Returns
Reference to changed calling matrix.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 18145 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::vanish ( )
inline

Set matrix to zero.

Sets every element of calling square matrix to be equal to zero and returns reference to the matrix changed. This function is faster than, for example, set(TR) with zero parameter passed.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (3);
srmatrix m(3);
m.randomize(0.,1.);
std::cout << m << std::endl;
std::cout << m.vanish ();
prints
1.747e-01 7.563e-01 5.163e-01
9.657e-01 6.619e-01 8.036e-01
6.392e-01 6.658e-01 6.495e-01
0.000e+00 0.000e+00 0.000e+00
0.000e+00 0.000e+00 0.000e+00
0.000e+00 0.000e+00 0.000e+00
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 18180 of file cvm.h.

template<typename TR>
basic_srmatrix& basic_srmatrix< TR >::randomize ( TR  dFrom,
TR  dTo 
)
inline

Randomizer.

Fills calling matrix with pseudo-random numbers distributed between dFrom and dTo. Function returns reference to the matrix changed.

Example:
using namespace cvm;
std::cout.setf (std::ios::scientific | std::ios::left);
std::cout.precision (7);
rmatrix m(3,4);
m.randomize(-2.,3.);
std::cout << m;
prints
9.6853542e-01 2.7761467e+00 2.3791009e+00 -3.4452345e-01
2.9029511e+00 -9.5519883e-01 -4.9131748e-01 -1.2561113e+00
1.5219886e+00 -1.4494461e+00 2.8193304e+00 4.8817408e-01
Parameters
[in]dFromFirst limit.
[in]dToSecond limit.
Returns
Reference to changed calling matrix.

Reimplemented from basic_rmatrix< TR >.

Reimplemented in basic_srsmatrix< TR >, and basic_srbmatrix< TR >.

Definition at line 18186 of file cvm.h.


The documentation for this class was generated from the following file: