Teuchos - Trilinos Tools Package  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Related Functions | List of all members
Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType > Class Template Reference

This class creates and provides basic support for symmetric, positive-definite dense matrices of templated type. More...

#include <Teuchos_SerialSymDenseMatrix.hpp>

Inheritance diagram for Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >:
Teuchos::CompObject Teuchos::BLAS< OrdinalType, ScalarType > Teuchos::DefaultBLASImpl< OrdinalType, ScalarType >

Public Types

typedef OrdinalType ordinalType
 Typedef for ordinal type. More...
 
typedef ScalarType scalarType
 Typedef for scalar type. More...
 
- Public Types inherited from Teuchos::DefaultBLASImpl< OrdinalType, ScalarType >
typedef details::GivensRotator
< ScalarType >::c_type 
rotg_c_type
 The type used for c in ROTG. More...
 

Public Member Functions

Constructor/Destructor Methods
 SerialSymDenseMatrix ()=default
 Default constructor; defines a zero size object. More...
 
 SerialSymDenseMatrix (OrdinalType numRowsCols, bool zeroOut=true)
 Basic constructor; defines a matrix of numRowsCols size and (optionally) initializes it. More...
 
 SerialSymDenseMatrix (DataAccess CV, bool upper, ScalarType *values, OrdinalType stride, OrdinalType numRowsCols)
 Set object values from two-dimensional array. More...
 
 SerialSymDenseMatrix (const SerialSymDenseMatrix< OrdinalType, ScalarType > &Source)
 Teuchos::SerialSymDenseMatrix copy constructor. More...
 
 SerialSymDenseMatrix (DataAccess CV, const SerialSymDenseMatrix< OrdinalType, ScalarType > &Source, OrdinalType numRowCols, OrdinalType startRowCol=0)
 Submatrix Copy Constructor. More...
 
virtual ~SerialSymDenseMatrix ()
 Teuchos::SerialSymDenseMatrix destructor. More...
 
Shaping Methods
int shape (OrdinalType numRowsCols)
 Set dimensions of a Teuchos::SerialSymDenseMatrix object; init values to zero. More...
 
int shapeUninitialized (OrdinalType numRowsCols)
 Set dimensions of a Teuchos::SerialSymDenseMatrix object; don't initialize values. More...
 
int reshape (OrdinalType numRowsCols)
 Reshape a Teuchos::SerialSymDenseMatrix object. More...
 
void setLower ()
 Specify that the lower triangle of the this matrix should be used. More...
 
void setUpper ()
 Specify that the upper triangle of the this matrix should be used. More...
 
Set methods.
SerialSymDenseMatrix
< OrdinalType, ScalarType > & 
operator= (const SerialSymDenseMatrix< OrdinalType, ScalarType > &Source)
 Copies values from one matrix to another. More...
 
SerialSymDenseMatrix
< OrdinalType, ScalarType > & 
assign (const SerialSymDenseMatrix< OrdinalType, ScalarType > &Source)
 Copies values from one matrix to another. More...
 
SerialSymDenseMatrix
< OrdinalType, ScalarType > & 
operator= (const ScalarType value)
 Set all values in the matrix to a constant value. More...
 
int putScalar (const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero(), bool fullMatrix=false)
 Set all values in the matrix to a constant value. More...
 
void swap (SerialSymDenseMatrix< OrdinalType, ScalarType > &B)
 Swap values between this matrix and incoming matrix. More...
 
int random (const ScalarType bias=0.1 *Teuchos::ScalarTraits< ScalarType >::one())
 Set all values in the active area (upper/lower triangle) of this matrix to be random numbers. More...
 
Accessor methods.
ScalarType & operator() (OrdinalType rowIndex, OrdinalType colIndex)
 Element access method (non-const). More...
 
const ScalarType & operator() (OrdinalType rowIndex, OrdinalType colIndex) const
 Element access method (const). More...
 
ScalarType * values () const
 Returns the pointer to the ScalarType data array contained in the object. More...
 
Query methods
bool upper () const
 Returns true if upper triangular part of this matrix has and will be used. More...
 
char UPLO () const
 Returns character value of UPLO used by LAPACK routines. More...
 
Mathematical Methods
SerialSymDenseMatrix
< OrdinalType, ScalarType > & 
operator*= (const ScalarType alpha)
 Inplace scalar-matrix product A = alpha*A. More...
 
SerialSymDenseMatrix
< OrdinalType, ScalarType > & 
operator+= (const SerialSymDenseMatrix< OrdinalType, ScalarType > &Source)
 Add another matrix to this matrix. More...
 
SerialSymDenseMatrix
< OrdinalType, ScalarType > & 
operator-= (const SerialSymDenseMatrix< OrdinalType, ScalarType > &Source)
 Subtract another matrix from this matrix. More...
 
Comparison methods.
bool operator== (const SerialSymDenseMatrix< OrdinalType, ScalarType > &Operand) const
 Equality of two matrices. More...
 
bool operator!= (const SerialSymDenseMatrix< OrdinalType, ScalarType > &Operand) const
 Inequality of two matrices. More...
 
Attribute methods.
OrdinalType numRows () const
 Returns the row dimension of this matrix. More...
 
OrdinalType numCols () const
 Returns the column dimension of this matrix. More...
 
OrdinalType stride () const
 Returns the stride between the columns of this matrix in memory. More...
 
bool empty () const
 Returns whether this matrix is empty. More...
 
Norm methods.
ScalarTraits< ScalarType >
::magnitudeType 
normOne () const
 Returns the 1-norm of the matrix. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
normInf () const
 Returns the Infinity-norm of the matrix. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
normFrobenius () const
 Returns the Frobenius-norm of the matrix. More...
 
I/O methods.
virtual std::ostream & print (std::ostream &os) const
 Print method. Defines the behavior of the std::ostream << operator. More...
 
- Public Member Functions inherited from Teuchos::CompObject
 CompObject ()
 Default constructor. More...
 
 CompObject (const CompObject &source)
 Copy Constructor. More...
 
virtual ~CompObject ()
 Destructor. More...
 
void setFlopCounter (const Flops &FlopCounter)
 Set the internal Teuchos::Flops() pointer. More...
 
void setFlopCounter (const CompObject &compObject)
 Set the internal Teuchos::Flops() pointer to the flop counter of another Teuchos::CompObject. More...
 
void unsetFlopCounter ()
 Set the internal Teuchos::Flops() pointer to 0 (no flops counted). More...
 
FlopsgetFlopCounter () const
 Get the pointer to the Teuchos::Flops() object associated with this object, returns 0 if none. More...
 
void resetFlops () const
 Resets the number of floating point operations to zero for this multi-std::vector. More...
 
double getFlops () const
 Returns the number of floating point operations with this multi-std::vector. More...
 
void updateFlops (int addflops) const
 Increment Flop count for this object. More...
 
void updateFlops (long int addflops) const
 Increment Flop count for this object. More...
 
void updateFlops (double addflops) const
 Increment Flop count for this object. More...
 
void updateFlops (float addflops) const
 Increment Flop count for this object. More...
 
- Public Member Functions inherited from Teuchos::BLAS< OrdinalType, ScalarType >
 BLAS (void)
 Default constructor. More...
 
 BLAS (const BLAS< OrdinalType, ScalarType > &)
 Copy constructor. More...
 
virtual ~BLAS (void)
 Destructor. More...
 
- Public Member Functions inherited from Teuchos::DefaultBLASImpl< OrdinalType, ScalarType >
 DefaultBLASImpl (void)
 Default constructor. More...
 
 DefaultBLASImpl (const DefaultBLASImpl< OrdinalType, ScalarType > &)
 Copy constructor. More...
 
virtual ~DefaultBLASImpl (void)
 Destructor. More...
 
template<typename alpha_type , typename A_type , typename x_type , typename beta_type >
void GEMV (ETransp trans, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const x_type *x, const OrdinalType &incx, const beta_type beta, ScalarType *y, const OrdinalType &incy) const
 Performs the matrix-vector operation: y <- alpha*A*x+beta*y or y <- alpha*A'*x+beta*y where A is a general m by n matrix. More...
 
template<typename A_type >
void TRMV (EUplo uplo, ETransp trans, EDiag diag, const OrdinalType &n, const A_type *A, const OrdinalType &lda, ScalarType *x, const OrdinalType &incx) const
 Performs the matrix-vector operation: x <- A*x or x <- A'*x where A is a unit/non-unit n by n upper/lower triangular matrix. More...
 
template<typename alpha_type , typename x_type , typename y_type >
void GER (const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy, ScalarType *A, const OrdinalType &lda) const
 Performs the rank 1 operation: A <- alpha*x*y'+A. More...
 
template<typename alpha_type , typename A_type , typename B_type , typename beta_type >
void GEMM (ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 General matrix-matrix multiply. More...
 
void SWAP (const OrdinalType &n, ScalarType *const x, const OrdinalType &incx, ScalarType *const y, const OrdinalType &incy) const
 Swap the entries of x and y. More...
 
template<typename alpha_type , typename A_type , typename B_type , typename beta_type >
void SYMM (ESide side, EUplo uplo, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 Performs the matrix-matrix operation: C <- alpha*A*B+beta*C or C <- alpha*B*A+beta*C where A is an m by m or n by n symmetric matrix and B is a general matrix. More...
 
template<typename alpha_type , typename A_type , typename beta_type >
void SYRK (EUplo uplo, ETransp trans, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 Performs the symmetric rank k operation: C <- alpha*A*A'+beta*C or C <- alpha*A'*A+beta*C, where alpha and beta are scalars, C is an n by n symmetric matrix and A is an n by k matrix in the first case or k by n matrix in the second case. More...
 
template<typename alpha_type , typename A_type >
void TRMM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
 Performs the matrix-matrix operation: B <- alpha*op(A)*B or B <- alpha*B*op(A) where op(A) is an unit/non-unit, upper/lower triangular matrix and B is a general matrix. More...
 
template<typename alpha_type , typename A_type >
void TRSM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
 Solves the matrix equations: op(A)*X=alpha*B or X*op(A)=alpha*B where X and B are m by n matrices, A is a unit/non-unit, upper/lower triangular matrix and op(A) is A or A'. The matrix X is overwritten on B. More...
 
void ROTG (ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
 Computes a Givens plane rotation. More...
 
void ROT (const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
 Applies a Givens plane rotation. More...
 
void SCAL (const OrdinalType &n, const ScalarType &alpha, ScalarType *x, const OrdinalType &incx) const
 Scale the vector x by the constant alpha. More...
 
void COPY (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
 Copy the vector x to the vector y. More...
 
template<typename alpha_type , typename x_type >
void AXPY (const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
 Perform the operation: y <- y+alpha*x. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
ASUM (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Sum the absolute values of the entries of x. More...
 
template<typename x_type , typename y_type >
ScalarType DOT (const OrdinalType &n, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy) const
 Form the dot product of the vectors x and y. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
NRM2 (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Compute the 2-norm of the vector x. More...
 
OrdinalType IAMAX (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Return the index of the element of x with the maximum magnitude. More...
 

Related Functions

(Note that these are not member functions.)

template<typename OrdinalType , typename ScalarType >
void symMatTripleProduct (ETransp transw, const ScalarType alpha, const SerialSymDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &W, SerialSymDenseMatrix< OrdinalType, ScalarType > &B)
 A templated, non-member, helper function for computing the matrix triple-product: B = alpha*W^T*A*W or B = alpha*W*A*W^T. More...
 

Detailed Description

template<typename OrdinalType, typename ScalarType>
class Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >

This class creates and provides basic support for symmetric, positive-definite dense matrices of templated type.

The Teuchos_SerialSymDenseMatrix class enables the construction and use of
symmetric, positive-definite, dense matrices of templated type.

The Teuchos::SerialSymDenseMatrix class is intended to provide full-featured support for solving
linear and eigen system problems for symmetric positive-definite (SPD) matrices.  It is written on
top of BLAS and LAPACK and thus has excellent performance and numerical capabilities.  Using this
class, one can either perform simple factorizations and solves or apply all the tricks available in
LAPACK to get the best possible solution for very ill-conditioned problems.

<b>Teuchos::SerialSymDenseMatrix vs. Teuchos::LAPACK</b>

The Teuchos::LAPACK class provides access to most of the same functionality as Teuchos::SerialSymDenseMatrix.
The primary difference is that Teuchos::LAPACK is a "thin" layer on top of LAPACK and
Teuchos::SerialSymDenseMatrix attempts to provide easy access to the more sophisticated aspects of
solving dense linear and eigensystems.

Constructing Teuchos::SerialSymDenseMatrix Objects

There are three Teuchos::SerialSymDenseMatrix constructors. The first constructs a zero-sized object which should be made to appropriate length using the Shape() or Reshape() functions and then filled with the [] or () operators. The second is a constructor that accepts user data as a 2D array, the third is a copy constructor. The second constructor has two data access modes (specified by the Teuchos::DataAccess argument):

  1. Copy mode - Allocates memory and makes a copy of the user-provided data. In this case, the user data is not needed after construction.
  2. View mode - Creates a "view" of the user data. In this case, the user data is required to remain intact for the life of the object.
Warning
View mode is extremely dangerous from a data hiding perspective. Therefore, we strongly encourage users to develop code using Copy mode first and only use the View mode in a secondary optimization phase.

Extracting Data from Teuchos::SerialSymDenseMatrix Objects

Once a Teuchos::SerialSymDenseMatrix is constructed, it is possible to view the data via access functions.

Warning
Use of these access functions cam be extremely dangerous from a data hiding perspective.

Vector and Utility Functions

Once a Teuchos::SerialSymDenseMatrix is constructed, several mathematical functions can be applied to the object. Specifically:

Examples:
DenseMatrix/cxx_main_sym.cpp.

Definition at line 121 of file Teuchos_SerialSymDenseMatrix.hpp.

Member Typedef Documentation

template<typename OrdinalType, typename ScalarType>
typedef OrdinalType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::ordinalType

Typedef for ordinal type.

Definition at line 126 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
typedef ScalarType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::scalarType

Typedef for scalar type.

Definition at line 128 of file Teuchos_SerialSymDenseMatrix.hpp.

Constructor & Destructor Documentation

template<typename OrdinalType, typename ScalarType>
Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::SerialSymDenseMatrix ( )
default

Default constructor; defines a zero size object.

Teuchos::SerialSymDenseMatrix objects defined by the default constructor should be sized with the Shape() or Reshape() functions. Values should be defined by using the [] or ()operators.

Note: By default the active part of the matrix is assumed to be the lower triangular part. To set the upper part as active, call SetUpper(). See Detailed Description section for further discussion.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::SerialSymDenseMatrix ( OrdinalType  numRowsCols,
bool  zeroOut = true 
)

Basic constructor; defines a matrix of numRowsCols size and (optionally) initializes it.

Parameters
numRowsCols- Number of rows and columns in the matrix.
zeroOut- Initializes values to 0 if true (default)

Creates a shaped matrix with numRowsCols rows and cols. All values are initialized to 0 when zeroOut is true. Values of this matrix should be set using the [] or the () operators.

Note
By default the active part of the matrix is assumed to be the lower triangular part. To set the upper part as active, call SetUpper(). See Detailed Description section for further discussion.

Definition at line 449 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::SerialSymDenseMatrix ( DataAccess  CV,
bool  upper,
ScalarType *  values,
OrdinalType  stride,
OrdinalType  numRowsCols 
)

Set object values from two-dimensional array.

Parameters
CV- Enumerated type set to Teuchos::Copy or Teuchos::View.
values- Pointer to an array of ScalarType. The first column starts at values, the second at values+stride, etc.
stride- The stride between the columns of the matrix in memory.
numRowsCols- Number of rows and columns in the matrix.
Note
By default the active part of the matrix is assumed to be the lower triangular part. To set the upper part as active, call SetUpper(). See Detailed Description section for further discussion.

Definition at line 461 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::SerialSymDenseMatrix ( const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Source)

Teuchos::SerialSymDenseMatrix copy constructor.

Definition at line 482 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::SerialSymDenseMatrix ( DataAccess  CV,
const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Source,
OrdinalType  numRowCols,
OrdinalType  startRowCol = 0 
)

Submatrix Copy Constructor.

Parameters
CV- Enumerated type set to Teuchos::Copy or Teuchos::View.
Source- Reference to another dense matrix from which values are to be copied.
numRowCols- The number of rows and columns in this matrix.
startRowCol- The row and column of Source from which the submatrix copy should start.

Creates a shaped matrix with numRowCols rows and columns, which is a submatrix of Source. If startRowCol are not given, then the submatrix is the leading submatrix of Source.

Definition at line 508 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::~SerialSymDenseMatrix ( )
virtual

Teuchos::SerialSymDenseMatrix destructor.

Definition at line 528 of file Teuchos_SerialSymDenseMatrix.hpp.

Member Function Documentation

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::shape ( OrdinalType  numRowsCols)

Set dimensions of a Teuchos::SerialSymDenseMatrix object; init values to zero.

Parameters
numRowsCols- Number of rows and columns in object.

Allows user to define the dimensions of a Teuchos::SerialSymDenseMatrix at any point. This function can be called at any point after construction. Any values that were previously in this object are destroyed and the resized matrix starts off with all zero values.

Returns
Integer error code, set to 0 if successful.

Definition at line 538 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::shapeUninitialized ( OrdinalType  numRowsCols)

Set dimensions of a Teuchos::SerialSymDenseMatrix object; don't initialize values.

Parameters
numRowsCols- Number of rows and columns in object.

Allows user to define the dimensions of a Teuchos::SerialSymDenseMatrix at any point. This function can be called at any point after construction. Any values that were previously in this object are destroyed. The resized matrix has uninitialized values.

Returns
Integer error code, set to 0 if successful.

Definition at line 550 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::reshape ( OrdinalType  numRowsCols)

Reshape a Teuchos::SerialSymDenseMatrix object.

Parameters
numRowsCols- Number of rows and columns in object.

Allows user to define the dimensions of a Teuchos::SerialSymDenseMatrix at any point. This function can be called at any point after construction. Any values that were previously in this object are copied into the new shape. If the new shape is smaller than the original, the upper left portion of the original matrix (the principal submatrix) is copied to the new matrix.

Returns
Integer error code, set to 0 if successful.

Definition at line 561 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
void Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::setLower ( )

Specify that the lower triangle of the this matrix should be used.

Warning
This may necessitate the matrix values be copied from the upper to lower portion of the matrix.

Definition at line 588 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
void Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::setUpper ( )

Specify that the upper triangle of the this matrix should be used.

Warning
This may necessitate the matrix values be copied from the lower to upper portion of the matrix.

Definition at line 599 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
SerialSymDenseMatrix< OrdinalType, ScalarType > & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator= ( const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Source)

Copies values from one matrix to another.

The operator= copies the values from one existing SerialSymDenseMatrix to another. If Source is a view (i.e. CV = Teuchos::View), then this method will return a view. Otherwise, it will return a copy of Source. this object will be resized if it is not large enough to copy Source into.

Definition at line 689 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
SerialSymDenseMatrix< OrdinalType, ScalarType > & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::assign ( const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Source)

Copies values from one matrix to another.

Copies the values from one existing SerialSymDenseMatrix to another if the dimension of both matrices are the same. If not, this matrix will be returned unchanged.

Definition at line 781 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
SerialSymDenseMatrix<OrdinalType, ScalarType>& Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator= ( const ScalarType  value)
inline

Set all values in the matrix to a constant value.

Parameters
value- Value to use;

Definition at line 264 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::putScalar ( const ScalarType  value = Teuchos::ScalarTraits<ScalarType>::zero(),
bool  fullMatrix = false 
)

Set all values in the matrix to a constant value.

Parameters
value- Value to use; zero if none specified.
fullMatrix- set full matrix entries to value, not just active portion of symmetric matrix.
Returns
Integer error code, set to 0 if successful.

Definition at line 610 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
void Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::swap ( SerialSymDenseMatrix< OrdinalType, ScalarType > &  B)

Swap values between this matrix and incoming matrix.

Swaps pointers and associated state without copying the matrix data.

Definition at line 643 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::random ( const ScalarType  bias = 0.1*Teuchos::ScalarTraits<ScalarType>::one())

Set all values in the active area (upper/lower triangle) of this matrix to be random numbers.

Note
The diagonal will be the sum of the off diagonal elements, plus a bias, so the matrix is SPD.

Definition at line 655 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
ScalarType & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator() ( OrdinalType  rowIndex,
OrdinalType  colIndex 
)
inline

Element access method (non-const).

Returns the element in the ith row and jth column if A(i,j) is specified.

Returns
Element from the specified rowIndex row and colIndex column.
Note
If the requested element lies in the inactive part of the matrix, then A(j,i) will be returned.
Warning
The validity of rowIndex and colIndex will only be checked if Teuchos is configured with –enable-teuchos-abc.

Definition at line 803 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
const ScalarType & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator() ( OrdinalType  rowIndex,
OrdinalType  colIndex 
) const
inline

Element access method (const).

Returns the element in the ith row and jth column if A(i,j) is specified.

Returns
Element from the specified rowIndex row and colIndex column.
Note
If the requested element lies in the inactive part of the matrix, then A(j,i) will be returned.
Warning
The validity of rowIndex and colIndex will only be checked if Teuchos is configured with –enable-teuchos-abc.

Definition at line 825 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
ScalarType* Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::values ( ) const
inline

Returns the pointer to the ScalarType data array contained in the object.

Note
The matrix values are only guaranteed to be stored in the active area of the matrix (upper/lower).

Definition at line 313 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::upper ( ) const
inline

Returns true if upper triangular part of this matrix has and will be used.

Definition at line 321 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
char Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::UPLO ( ) const
inline

Returns character value of UPLO used by LAPACK routines.

Definition at line 324 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
SerialSymDenseMatrix< OrdinalType, ScalarType > & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator*= ( const ScalarType  alpha)

Inplace scalar-matrix product A = alpha*A.

Scale a matrix, entry-by-entry using the value alpha. This method is sensitive to the UPLO() parameter.

Parameters
alpha- Scalar to multiply with A.

Definition at line 750 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
SerialSymDenseMatrix< OrdinalType, ScalarType > & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator+= ( const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Source)

Add another matrix to this matrix.

Add Source to this if the dimension of both matrices are the same. If not, this matrix will be returned unchanged.

Definition at line 757 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
SerialSymDenseMatrix< OrdinalType, ScalarType > & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator-= ( const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Source)

Subtract another matrix from this matrix.

Subtract Source from this if the dimension of both matrices are the same. If not, this matrix will be returned unchanged.

Definition at line 769 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
bool Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator== ( const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Operand) const

Equality of two matrices.

Returns
True if this matrix and Operand are of the same shape (rows / columns) and have the same entries in the active (upper / lower triangular) area of the matrix, else False will be returned.

Definition at line 933 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
bool Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::operator!= ( const SerialSymDenseMatrix< OrdinalType, ScalarType > &  Operand) const

Inequality of two matrices.

Returns
True if this matrix and Operand of not of the same shape (rows / columns) or don't have the same entries in the active (upper / lower triangular), else False will be returned.

Definition at line 953 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::numRows ( ) const
inline

Returns the row dimension of this matrix.

Definition at line 374 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::numCols ( ) const
inline

Returns the column dimension of this matrix.

Definition at line 377 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::stride ( ) const
inline

Returns the stride between the columns of this matrix in memory.

Definition at line 380 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::empty ( ) const
inline

Returns whether this matrix is empty.

Definition at line 383 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
ScalarTraits< ScalarType >::magnitudeType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::normOne ( ) const

Returns the 1-norm of the matrix.

Definition at line 851 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
ScalarTraits< ScalarType >::magnitudeType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::normInf ( ) const

Returns the Infinity-norm of the matrix.

Definition at line 857 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
ScalarTraits< ScalarType >::magnitudeType Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::normFrobenius ( ) const

Returns the Frobenius-norm of the matrix.

Definition at line 900 of file Teuchos_SerialSymDenseMatrix.hpp.

template<typename OrdinalType , typename ScalarType >
std::ostream & Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType >::print ( std::ostream &  os) const
virtual

Print method. Defines the behavior of the std::ostream << operator.

Definition at line 1012 of file Teuchos_SerialSymDenseMatrix.hpp.

Friends And Related Function Documentation

template<typename OrdinalType , typename ScalarType >
void symMatTripleProduct ( ETransp  transw,
const ScalarType  alpha,
const SerialSymDenseMatrix< OrdinalType, ScalarType > &  A,
const SerialDenseMatrix< OrdinalType, ScalarType > &  W,
SerialSymDenseMatrix< OrdinalType, ScalarType > &  B 
)
related

A templated, non-member, helper function for computing the matrix triple-product: B = alpha*W^T*A*W or B = alpha*W*A*W^T.

Parameters
transw- [in] Compute B = alpha*W^T*A*W if transw = Teuchos::TRANS, else compute B = alpha*W*A*W^T if transw = Teuchos::NO_TRANS.
alpha- [in] The scaling factor.
A- [in] SerialSymDenseMatrix
W- [in] SerialDenseMatrix
B- [out] SerialSymDenseMatrix
Note
The syntax for calling this function is: Teuchos::symMatTripleProduct<int,double>( Teuchos::TRANS, alpha, A, W, B )

Definition at line 76 of file Teuchos_SerialDenseHelpers.hpp.


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