ConicBundle
Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
CH_Matrix_Classes::Symmatrix Class Reference

Matrix class of symmetric matrices with real values of type Real More...

#include <symmat.hxx>

Inheritance diagram for CH_Matrix_Classes::Symmatrix:
CH_Matrix_Classes::Memarrayuser ConicBundle::DensePSCPrimal

Public Member Functions

Constructors, Destructor, and Initialization (Members)
 Symmatrix ()
 empty matrix
 
 Symmatrix (const Symmatrix &A, double d=1.)
 copy constructor, *this=d*A
 
 Symmatrix (Integer nr)
 generate a matrix of size nr x nr but WITHOUT initializing the memory More...
 
 Symmatrix (Integer nr, Real d)
 generate a matrix of size nr x nr initializing all elements to the value d
 
 Symmatrix (Integer nr, const Real *dp)
 generate a matrix of size nr x nr initializing the elements from the (one dimensional) array dp, which must have the elements arranged consecutively in internal order
 
 ~Symmatrix ()
 
void set_init (bool)
 after external initialization, call matrix.set_init(true) (not needed if CONICBUNDLE_DEBUG is undefined)
 
bool get_init () const
 returns true if the matrix has been declared initialized (not needed if CONICBUNDLE_DEBUG is undefined)
 
Symmatrixinit (const Symmatrix &A, double d=1.)
 initialize to *this=A*d
 
Symmatrixinit (const Matrix &A, double d=1.)
 initialize to this=d(A+transpose(A))/2.
 
Symmatrixinit (const Indexmatrix &A, double d=1.)
 initialize to this=d(A+transpose(A))/2.
 
Symmatrixinit (const Sparsesym &A, Real d=1.)
 initialize to *this=A*d
 
Symmatrixinit (Integer nr, Real d)
 intialize *this to a matrix of size nr x nr initializing all elements to the value d
 
Symmatrixinit (Integer nr, const Real *dp)
 generate a matrix of size nr x nc initializing the elements from the (one dimensional) array dp which must have the elements arranged consecutively in internal order
 
void newsize (Integer n)
 resize the matrix to nr x nr elements but WITHOUT initializing the memory More...
 
Conversions from other Matrix Classes (Members)
 Symmatrix (const Matrix &, double d=1.)
 (this)=d(A+transpose(A))/2.
 
 Symmatrix (const Indexmatrix &, double d=1.)
 (this)=d(A+transpose(A))/2.
 
 Symmatrix (const Sparsesym &A, Real d=1.)
 (*this)=d*A
 
Size and Type Information (Members)
void dim (Integer &_nr, Integer &_nc) const
 returns the number of rows in _nr and _nc
 
Integer dim () const
 returns the dimension rows * columns when the matrix is regarded as a vector
 
Integer rowdim () const
 returns the row dimension
 
Integer coldim () const
 returns the column dimension
 
Mtype get_mtype () const
 returns the type of the matrix, MTsymmetric
 
Indexing and Submatrices (Members)
Realoperator() (Integer i, Integer j)
 returns reference to element (i,j) of the matrix (rowindex i, columnindex j)
 
Realoperator() (Integer i)
 returns reference to element (i) of the matrix if regarded as vector of stacked columns [element (irowdim, i/rowdim)]
 
Real operator() (Integer i, Integer j) const
 returns value of element (i,j) of the matrix (rowindex i, columnindex j)
 
Real operator() (Integer i) const
 returns value of element (i) of the matrix if regarded as vector of stacked columns [element (irowdim, i/rowdim)]
 
Matrix col (Integer i) const
 returns column i copied to a new Matrix
 
Matrix row (Integer i) const
 returns row i copied to a new Matrix
 
Matrix cols (const Indexmatrix &vec) const
 returns a matrix of size this->rowdim() x vec.dim(), with column i a copy of column vec(i) of *this
 
Matrix rows (const Indexmatrix &vec) const
 returns a matrix of size vec.dim() x this->coldim(), with row i a copy of row vec(i) of *this
 
Symmatrixswapij (Integer i, Integer j)
 swaps rows (and columns) i and j
 
Symmatrixpivot_permute (const Indexmatrix &piv, bool inverse=false)
 for i=0 to rowdim row (and column) i of this matrix is swapped with row piv(j); for inverse=true the inverse permutation is generated
 
Symmatrixprincipal_submatrix (const Indexmatrix &ind, Symmatrix &S) const
 returns S and in S the principal submatrix indexed by ind (multiple indices are allowed)
 
Symmatrix principal_submatrix (const Indexmatrix &ind) const
 returns the principal submatrix indexed by ind (multiple indices are allowed)
 
Symmatrixdelete_principal_submatrix (const Indexmatrix &ind, bool sorted_increasingly=false)
 returns this afte deleting the principal submatrix indexed by ind (no repetitions!);
 
Symmatrixenlarge_below (Integer addn)
 increases the order of the matrix by appending storage for further addn rows and columns (marked as not initiliazed if addn>0, no changes if addn<=0)
 
Symmatrixenlarge_below (Integer addn, Real d)
 increases the order of the matrix by appending storage for further addn rows and columns initialized to d (no changes if addn<=0);
 
Realget_store ()
 returns the current address of the internal value array; use cautiously, do not use delete!
 
const Realget_store () const
 returns the current address of the internal value array; use cautiously!
 
BLAS-like Routines (Members)
Symmatrixxeya (const Symmatrix &A, Real d=1.)
 sets *this=d*A and returns *this
 
Symmatrixxpeya (const Symmatrix &A, Real d=1.)
 sets *this+=d*A and returns *this
 
Usual Arithmetic Operators (Members)
Symmatrixoperator= (const Symmatrix &A)
 
Symmatrixoperator+= (const Symmatrix &A)
 
Symmatrixoperator-= (const Symmatrix &A)
 
Symmatrixoperator%= (const Symmatrix &A)
 ATTENTION: this is redefined as the Hadamard product, (*this)(i,j)=(*this)(i,j)*A(i,j) for all i<=j.
 
Symmatrix operator- () const
 
Symmatrixoperator*= (Real d)
 
Symmatrixoperator/= (Real d)
 ATTENTION: d is NOT checked for 0.
 
Symmatrixoperator+= (Real d)
 sets (*this)(i,j)+=d for all i<=j
 
Symmatrixoperator-= (Real d)
 sets (*this)(i,j)-=d for all i<=j
 
Symmatrixtranspose ()
 transposes itself (at almost no cost)
 
Connections to other Classes (Members)
Symmatrixxeya (const Matrix &A, Real d=1.)
 sets this=d(A+transpose(A))/2. and returns *this
 
Symmatrixxpeya (const Matrix &A, Real d=1.)
 sets this+=d(A+transpose(A))/2. and returns *this
 
Symmatrixxeya (const Indexmatrix &A, Real d=1.)
 sets this=d(A+transpose(A))/2. and returns *this
 
Symmatrixxpeya (const Indexmatrix &A, Real d=1.)
 sets this+=d(A+transpose(A))/2. and returns *this
 
Symmatrixxeya (const Sparsesym &A, Real d=1.)
 sets *this=d*A and returns *this
 
Symmatrixxpeya (const Sparsesym &A, Real d=1.)
 sets *this+=d*A and returns *this
 
Symmatrixxetriu_yza (const Matrix &A, const Matrix &B, Real d=1.)
 sets *this(i,j), i<=j to the upper triangle of the matrix product d*transpose(A)*B
 
Symmatrixxpetriu_yza (const Matrix &A, const Matrix &B, Real d=1.)
 adds to *this(i,j), i<=j the upper triangle of the matrix product d*transpose(A)*B
 
Symmatrixxetriu_yza (const Sparsemat &A, const Matrix &B, Real d=1.)
 sets *this(i,j), i<=j to the upper triangle of the matrix product d*transpose(A)*B
 
Symmatrixxpetriu_yza (const Sparsemat &A, const Matrix &B, Real d=1.)
 adds to *this(i,j), i<=j the upper triangle of the matrix product d*transpose(A)*B
 
Symmatrixoperator= (const Sparsesym &A)
 
Symmatrixoperator+= (const Sparsesym &A)
 
Symmatrixoperator-= (const Sparsesym &A)
 
Numerical Methods (Members)
Symmatrixshift_diag (Real s)
 shifts the diagonal by s, i.e., (*this)(i,i)+=s for all i
 
int LDLfactor (Real tol=1e-10)
 computes LDLfactorization (implemented only for positive definite matrices so far, no pivoting), (*this) is overwritten by the factorization; returns 1 if diagonal elements go below tol
 
int LDLsolve (Matrix &x) const
 computes, after LDLfactor was executed succesfully, the solution to (*old_this)x=rhs; rhs is overwritten by the solution; always returns 0; NOTE: there is NO check against division by zero
 
int LDLinverse (Symmatrix &S) const
 computes, after LDLfactor was executed succesfully, the inverse to (*old_this) and stores it in S (numerically not too wise); always returns 0; NOTE: there is NO check against division by zero
 
int Chol_factor (Real tol=1e-10)
 computes the Cholesky factorization, for positive definite matrices only, (*this) is overwritten by the factorization; there is no pivoting; returns 1 if diagonal elements go below tol
 
int Chol_solve (Matrix &x) const
 computes, after Chol_factor was executed succesfully, the solution to (*old_this)x=rhs; rhs is overwritten by the solution; always returns 0; NOTE: there is NO check against division by zero
 
int Chol_inverse (Symmatrix &S) const
 computes, after Chol_factor was executed succesfully, the inverse to (*old_this) and stores it in S (numerically not too wise); always returns 0; NOTE: there is NO check against division by zero
 
int Chol_Lsolve (Matrix &rhs) const
 computes, after Chol_factor into LL^T was executed succesfully, the solution to Lx=rhs; rhs is overwritten by the solution; always returns 0; NOTE: there is NO check against division by zero
 
int Chol_Ltsolve (Matrix &rhs) const
 computes, after Chol_factor into LL^T was executed succesfully, the solution to L^Tx=rhs; rhs is overwritten by the solution; always returns 0; NOTE: there is NO check against division by zero
 
int Chol_scaleLi (Symmatrix &S) const
 computes, after Chol_factor into LL^T was executed succesfully, L^{-1}SL^{-T} overwriting S
 
int Chol_scaleLt (Symmatrix &S) const
 computes, after Chol_factor into LL^T was executed succesfully, L^TSL overwriting S
 
int Chol_Lmult (Matrix &rhs) const
 computes, after Chol_factor into LL^T was executed succesfully, L*rhs, overwriting rhs by the result; always returns 0;
 
int Chol_Ltmult (Matrix &rhs) const
 computes, after Chol_factor into LL^T was executed succesfully, L^Trhs, overwriting rhs by the result; always returns 0;
 
int Chol_factor (Indexmatrix &piv, Real tol=1e-10)
 computes the Cholesky factorization with pivoting, for positive semidefinite matrices only, (*this) is overwritten by the factorization; on termination piv.dim() is the number of positive pivots>=tol; returns 1 if negative diagonal element is encountered during computations, 0 otherwise.
 
int Chol_solve (Matrix &x, const Indexmatrix &piv) const
 computes, after Chol_factor(Indexmatrix&,Real) with pivoting was executed succesfully, the solution to (*old_this)*x=rhs(piv); rhs is overwritten by the solution arranged in original unpermuted order; always returns 0; NOTE: there is NO check against division by zero
 
int Chol_inverse (Symmatrix &S, const Indexmatrix &piv) const
 computes, after Chol_factor(Indexmatrix&,Real) with pivoting was executed succesfully, the inverse to (*old_this) and stores it in S (the pivoting permutation is undone in S); NOTE: there is NO check against division by zero
 
int Aasen_factor (Indexmatrix &piv)
 computes Aasen factorization LTL^T with pivoting, where L is unit lower triangular with first colum e_1 and T is tridiagonal; (*this) is overwritten by the factorization, with column i of L being stored in column i-1 of (*this); always returns 0;
 
int Aasen_Lsolve (Matrix &x) const
 computes, after Aasen_factor into LTL^T was executed, the solution to Lx=rhs; rhs is overwritten by the solution; always returns 0;
 
int Aasen_Ltsolve (Matrix &x) const
 computes, after Aasen_factor into LTL^T was executed, the solution to L^Tx=rhs; rhs is overwritten by the solution; always returns 0;
 
int Aasen_tridiagsolve (Matrix &x) const
 computes, after Aasen_factor into LTL^T was executed, the solution to Tx=rhs; rhs is overwritten by the solution;if the solution fails due to division by zero (=system not solvable) the return value is -(rowindex+1) where this occured in the backsolve
 
int Aasen_solve (Matrix &x, const Indexmatrix &piv) const
 computes, after Aasen_factor into LTL^T was executed, the solution to (*old_this)x=rhs; rhs is overwritten by the solution; if the solution fails due to division by zero (=system not solvable) the return value is -(rowindex+1) where this occured in the backsolve
 
Integer eig (Matrix &P, Matrix &d, bool sort_non_decreasingly=true) const
 computes an eigenvalue decomposition P*Diag(d)*tranpose(P)=(*this) by symmetric QR; returns 0 on success,
 
Input, Output (Members)
void display (std::ostream &out, int precision=0, int width=0, int screenwidth=0) const
 displays a matrix in a pretty way for bounded screen widths; for variables of value zero default values are used. More...
 
void mfile_output (std::ostream &out, int precision=16, int width=0) const
 outputs a matrix A in the format "[ A(0,1) ... A(0,nc-1)\n ... A(nr-1,nc-1)];\n" so that it can be read e.g. by octave as an m-file More...
 

Private Member Functions

void init_to_zero ()
 initialize the matrix to a 0x0 matrix without storage
 
Integer tred2 (Integer nm, Integer n, Real *a, Real *d, Real *e, Real *z) const
 a subroutine needed internally for eigenvalue computations (eigval.cxx)
 
Integer imtql2 (Integer nm, Integer n, Real *d, Real *e, Real *z) const
 a subroutine needed internally for eigenvalue computations (eigval.cxx)
 

Private Attributes

Integer mem_dim
 amount of memory currently allocated
 
Integer nr
 number of rows = number of columns
 
Realm
 lower triangle stored columnwise (a11,a21,...,anr1,a22,.....)
 
bool is_init
 flag whether memory is initialized, it is only used if CONICBUNDLE_DEBUG is defined
 

Static Private Attributes

static const Mtype mtype
 used for MatrixError templates (runtime type information was not yet existing)
 

Friends

class Matrix
 
class Sparsesym
 
class Sparsemat
 
Indexing and Submatrices (Friends)
Matrix diag (const Symmatrix &A)
 returns a column vector v consisting of the elements v(i)=(*this)(i,i), 0<=i<row dimension
 
Symmatrix Diag (const Matrix &A)
 returns a symmetric diagonal matrix S of order A.dim() with vec(A) on the diagonal, i.e., S(i,i)=A(i) for all i and S(i,j)=0 for i!=j
 
void swap (Symmatrix &A, Symmatrix &B)
 swap the content of the two matrices A and B (involves no copying)
 
BLAS-like Routines (Friends)
Symmatrixrankadd (const Matrix &A, Symmatrix &C, Real alpha, Real beta, int trans)
 returns C=beta*C+alpha* A*A^T, where A may be transposed. If beta==0. then C is initiliazed to the correct size.
 
Symmatrixscaledrankadd (const Matrix &A, const Matrix &D, Symmatrix &C, Real alpha, Real beta, int trans)
 returns C=beta*C+alpha* A*D*A^T, where D is a vector representing a diagonal matrix and A may be transposed; if beta==0. then C is initialized to the correct size
 
Symmatrixrank2add (const Matrix &A, const Matrix &B, Symmatrix &C, Real alpha, Real beta, int trans)
 returns C=beta*C+alpha*(A*B^T+B*A^T)/2 [or for transposed (A^T*B+B^T*A)/2]. If beta==0. then C is initiliazed to the correct size.
 
Symmatrixxbpeya (Symmatrix &x, const Symmatrix &y, Real alpha, Real beta)
 returns x= alpha*y+beta*x; if beta==0. then x is initialized to the correct size
 
Symmatrixxeyapzb (Symmatrix &x, const Symmatrix &y, const Symmatrix &z, Real alpha, Real beta)
 returns x= alpha*y+beta*z; x is initialized to the correct size
 
Matrixgenmult (const Symmatrix &A, const Matrix &B, Matrix &C, Real alpha, Real beta, int btrans)
 returns C=beta*C+alpha*A*B, where B may be transposed; C must not be equal to B; if beta==0. then C is initialized to the correct size More...
 
Matrixgenmult (const Matrix &A, const Symmatrix &B, Matrix &C, Real alpha, Real beta, int atrans)
 returns C=beta*C+alpha*A*B, where A may be transposed; C must not be equal to A; if beta==0. then C is initialized to the correct size More...
 
Usual Arithmetic Operators (Friends)
Matrix operator* (const Symmatrix &A, const Symmatrix &B)
 returns a Matrix that equals A*B
 
Symmatrix operator% (const Symmatrix &A, const Symmatrix &B)
 ATTENTION: this is redefined as the Hadamard product and sets (i,j)=A(i,j)*B(i,j) for all i<=j.
 
Symmatrix operator+ (const Symmatrix &A, const Symmatrix &B)
 returns a Matrix that equals A+B
 
Symmatrix operator- (const Symmatrix &A, const Symmatrix &B)
 returns a Matrix that equals A-B
 
Matrix operator* (const Symmatrix &A, const Matrix &B)
 returns a Matrix that equals A*B
 
Matrix operator* (const Matrix &A, const Symmatrix &B)
 returns a Matrix that equals A*B
 
Matrix operator+ (const Symmatrix &A, const Matrix &B)
 returns a Matrix that equals A+B
 
Matrix operator+ (const Matrix &A, const Symmatrix &B)
 returns a Matrix that equals A+B
 
Matrix operator- (const Symmatrix &A, const Matrix &B)
 returns a Matrix that equals A-B
 
Matrix operator- (const Matrix &A, const Symmatrix &B)
 returns a Matrix that equals A-B
 
Symmatrix operator* (const Symmatrix &A, Real d)
 returns a Symmatrix that equals A*d
 
Symmatrix operator* (Real d, const Symmatrix &A)
 returns a Symmatrix that equals A*d
 
Symmatrix operator/ (const Symmatrix &A, Real d)
 returns a Symmatrix that equals A/d; ATTENTION: no check against division by zero
 
Symmatrix operator+ (const Symmatrix &A, Real d)
 returns (i,j)=A(i,j)+d for all i<=j
 
Symmatrix operator+ (Real d, const Symmatrix &A)
 returns (i,j)=A(i,j)+d for all i<=j
 
Symmatrix operator- (const Symmatrix &A, Real d)
 returns (i,j)=A(i,j)-d for all i<=j
 
Symmatrix operator- (Real d, const Symmatrix &A)
 returns (i,j)=d-A(i,j) for all i<=j
 
Symmatrix transpose (const Symmatrix &A)
 (drop it or use a constructor instead)
 
Connections to other Classes (Friends)
Matrixgenmult (const Symmatrix &A, const Sparsemat &B, Matrix &C, Real alpha, Real beta, int btrans)
 returns C=beta*C+alpha*A*B, where B may be transposed; if beta==0. then C is initialized to the correct size
 
Matrixgenmult (const Sparsemat &A, const Symmatrix &B, Matrix &C, Real alpha, Real beta, int atrans)
 returns C=beta*C+alpha*A*B, where A may be transposed; if beta==0. then C is initialized to the correct size
 
Symmatrixrankadd (const Sparsemat &A, Symmatrix &C, Real alpha, Real beta, int trans)
 returns C=beta*C+alpha* A*A^T, where A may be transposed; if beta==0. then C is initialized to the correct size
 
Symmatrixscaledrankadd (const Sparsemat &A, const Matrix &D, Symmatrix &C, Real alpha, Real beta, int trans)
 returns C=beta*C+alpha* A*D*A^T, where D is a vector representing a diagonal matrix and A may be transposed; if beta==0. then C is initialized to the correct size
 
Symmatrixrank2add (const Sparsemat &A, const Matrix &B, Symmatrix &C, Real alpha, Real beta, int trans)
 returns C=beta*C+alpha*(A*B^T+B*A^T)/2 [or for transposed (A^T*B+B^T*A)/2]. If beta==0. then C is initiliazed to the correct size.
 
Elementwise Operations (Friends)
Symmatrix abs (const Symmatrix &A)
 returns a Symmatrix with elements abs(A(i,j))
 
Comparisons, Max, Min, Sort, Find (Friends)
Matrix minrows (const Symmatrix &A)
 returns a row vector holding in each column the minimum over all rows in this column
 
Matrix mincols (const Symmatrix &A)
 returns a column vector holding in each row the minimum over all columns in this row
 
Real min (const Symmatrix &A)
 returns the minimum value over all elements of the matrix
 
Matrix maxrows (const Symmatrix &A)
 returns a row vector holding in each column the maximum over all rows in this column
 
Matrix maxcols (const Symmatrix &A)
 returns a column vector holding in each row the maximum over all columns in this row
 
Real max (const Symmatrix &A)
 returns the maximum value over all elements of the matrix
 
Input, Output (Friends)
std::ostream & operator<< (std::ostream &o, const Symmatrix &A)
 output format (nr and nc are Integer values, all others Real values):
nr nc \n A(1,1) A(1,2) ... A(1,nc) \n A(2,1) ... A(nr,nc) \n
 
std::istream & operator>> (std::istream &i, Symmatrix &A)
 input format (nr and nc are Integer values, all others Real values):
nr nc \n A(1,1) A(1,2) ... A(1,nc) \n A(2,1) ... A(nr,nc) \n
 

Numerical Methods (Friends)

Real trace (const Symmatrix &A)
 returns the sum of the diagonal elements A(i,i) over all i
 
Real ip (const Symmatrix &A, const Symmatrix &B)
 returns the usual inner product of A and B, i.e., the sum of A(i,j)*B(i,j) over all i,j
 
Real ip (const Matrix &A, const Symmatrix &B)
 returns the usual inner product of A and B, i.e., the sum of A(i,j)*B(i,j) over all i,j
 
Real ip (const Symmatrix &A, const Matrix &B)
 returns the usual inner product of A and B, i.e., the sum of A(i,j)*B(i,j) over all i,j
 
Real norm2 (const Symmatrix &A)
 returns the Frobenius norm of A, i.e., the square root of the sum of A(i,j)*A(i,j) over all i,j
 
Matrix sumrows (const Symmatrix &A)
 returns a row vector holding the sum over all rows, i.e., (1 1 ... 1)*A
 
Matrix sumcols (const Symmatrix &A)
 returns a column vector holding the sum over all columns, i.e., A*(1 1 ... 1)^T
 
Real sum (const Symmatrix &A)
 returns the sum over all elements of A, i.e., (1 1 ... 1)*A*(1 1 ... 1)^T
 
void svec (const Symmatrix &A, Matrix &v, Real a, bool add, Integer startindex_vec, Integer startindex_A, Integer blockdim)
 the symmetric vec operator stacks the lower triangle of A to a n*(n+1)/2 vector with the same norm2 as A; here it sets svec(A)=[a11,sqrt(2)a12,...,sqrt(2)a1n,a22,...,sqrt(2)a(n-1,n),ann]', multiplies it by a and sets or adds (if add==true) it to v starting from startindex_vec possibly restricted to the subblock of order blockdim (whenever >=0, else blockdim is set to A.rowdim()-startindex_A) starting from startindex_A (must be >=0); if add==false and startindex_vec<0 then vec is also reinitialzed to the appropriate size
 
Matrix svec (const Symmatrix &A)
 the symmetric vec operator, stacks the lower triangle of A to a n*(n+1)/2 vector with the same norm2 as A; i.e., it returns svec(A)=[a11,sqrt(2)a12,...,sqrt(2)a1n,a22,...,sqrt(2)a(n-1,n),ann]'
 
void sveci (const Matrix &v, Symmatrix &A, Real a, bool add, Integer startindex_vec, Integer startindex_A, Integer blockdim)
 the inverse operator to svec, extracts from v at startindex_vec (>=0) the symmetric matrix of blockdim adding its mutliple by a into A starting at startindex_A; if add==false and startindex_A<0 A is initialized to the size of blockdim; if the latter is also negative then v.dim()-startindex_vec must match an exact order and matrix A is initialized to this size. In all other cases the size of the symmetric matrix determines the missing parameters and vec.dim-startindex_vec
 
Symmatrix skron (const Symmatrix &A, const Symmatrix &B, Real alpha, bool add, Integer startindex_S)
 the symmetric Kronecker product, defined via (A skron B)svec(C)=(BCA'+ACB')/2; sets or adds (if add==true) the symmetric matrix a*(A skron B) into S starting at startindex_S; if add==false and startindex_S<0, S is initialzed to the correct size
 
void skron (const Symmatrix &A, const Symmatrix &B, Symmatrix &S, Real a, bool add, Integer startindex_S)
 def symmetric Kronecker product (A skron B)svec(C)=(BCA'+ACB')/2; sets S=alpha*(A skron B) or S*=... (if add==true) possibly shifted to the block starting at startindex_S; if add==false and startindex_S<0, S is initialzed to the correct size
 
void symscale (const Symmatrix &A, const Matrix &B, Symmatrix &S, Real a, Real b, int btrans)
 sets S=beta*S+alpha*B'*A*B for symmatrix A and matrix B
 
void init_svec (Integer nr, const Real *dp, Integer incr=1, Real d=1.)
 
void store_svec (Real *dp, Integer incr=1, Real d=1.) const
 

Additional Inherited Members

- Protected Member Functions inherited from CH_Matrix_Classes::Memarrayuser
 Memarrayuser ()
 if memarray is NULL, then a new Memarray is generated. In any case the number of users of the Memarray is incremented
 
virtual ~Memarrayuser ()
 the number of users is decremented and the Memarray memory manager is destructed, if the number is zero.
 
- Static Protected Attributes inherited from CH_Matrix_Classes::Memarrayuser
static Memarraymemarray
 pointer to common memory manager for all Memarrayusers, instantiated in memarray.cxx
 

Detailed Description

Matrix class of symmetric matrices with real values of type Real

Internally a symmetric matrix of size nr x nr is stored in a one dimensional array of Real variables, the elements of the lower triangle are arranged in columnwise order (a11,a21,...,anr1,a22,a32,...).

Any matrix element can be indexed by (i,j) or directly by the one dimensional index (i+j*nr). The latter view directly corresponds to the vec() operator often used in the linear algebra literature, i.e., the matrix is transformed to a vector by stacking the columns on top of each other.

NOTE: Any change of A(i,j) also changes A(j,i) as both variables are identical!

Constructor & Destructor Documentation

◆ Symmatrix()

CH_Matrix_Classes::Symmatrix::Symmatrix ( Integer  nr)
inline

generate a matrix of size nr x nr but WITHOUT initializing the memory

If initializing the memory externally and CONICBUNDLE_DEBUG is defined, please use set_init() via matrix.set_init(true) in order to avoid warnings concerning improper initialization

References init_to_zero(), and newsize().

Member Function Documentation

◆ display()

void CH_Matrix_Classes::Symmatrix::display ( std::ostream &  out,
int  precision = 0,
int  width = 0,
int  screenwidth = 0 
) const

displays a matrix in a pretty way for bounded screen widths; for variables of value zero default values are used.

Parameters
outoutput stream
precisionnumber of most significant digits, default=4
widthfield width, default = precision+6
screenwidthmaximum number of characters in one output line, default = 80

Referenced by ConicBundle::CMsymdense::display(), and transpose().

◆ mfile_output()

void CH_Matrix_Classes::Symmatrix::mfile_output ( std::ostream &  out,
int  precision = 16,
int  width = 0 
) const

outputs a matrix A in the format "[ A(0,1) ... A(0,nc-1)\n ... A(nr-1,nc-1)];\n" so that it can be read e.g. by octave as an m-file

Parameters
outoutput stream
precisionnumber of most significant digits, default=16
widthfield width, default = precision+6

Referenced by transpose().

◆ newsize()

void CH_Matrix_Classes::Symmatrix::newsize ( Integer  n)

resize the matrix to nr x nr elements but WITHOUT initializing the memory

If initializing the memory externally and CONICBUNDLE_DEBUG is defined, please use set_init() via matrix.set_init(true) in order to avoid warnings concerning improper initialization

Referenced by get_init(), init(), ConicBundle::UQPSolver::init_size(), ConicBundle::CMsingleton::project(), Symmatrix(), and CH_Matrix_Classes::xeyapzb().

Friends And Related Function Documentation

◆ genmult [1/2]

Matrix& genmult ( const Symmatrix A,
const Matrix B,
Matrix C,
Real  alpha,
Real  beta,
int  btrans 
)
friend

returns C=beta*C+alpha*A*B, where B may be transposed; C must not be equal to B; if beta==0. then C is initialized to the correct size

returns C=beta*C+alpha*A*B, where B may be transposed; C must not be equal to B; if beta==0. then C is initialized to the correct size

Referenced by get_store(), CH_Matrix_Classes::Matrix::init(), CH_Matrix_Classes::operator*(), and transpose().

◆ genmult [2/2]

Matrix& genmult ( const Matrix A,
const Symmatrix B,
Matrix C,
Real  alpha,
Real  beta,
int  atrans 
)
friend

returns C=beta*C+alpha*A*B, where A may be transposed; C must not be equal to A; if beta==0. then C is initialized to the correct size

returns C=beta*C+alpha*A*B, where A may be transposed; C must not be equal to A; if beta==0. then C is initialized to the correct size


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