MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Private Member Functions | Private Attributes | List of all members
ConstrainedOptPack::MatrixIdentConcatStd Class Reference

Concrete implementation class for a matrix vertically concatonated with an identity matrix. More...

#include <ConstrainedOptPack_MatrixIdentConcatStd.hpp>

Inheritance diagram for ConstrainedOptPack::MatrixIdentConcatStd:
Inheritance graph
[legend]

Private Member Functions

void assert_initialized () const
 

Private Attributes

AbstractLinAlgPack::VectorSpacespace_cols
 
AbstractLinAlgPack::VectorSpacespace_rows
 
AbstractLinAlgPack::MatrixOpD
 
RangePack::Range1D D_rng
 
RangePack::Range1D I_rng
 
value_type alpha_
 
BLAS_Cpp::Transp D_trans_
 

Public types

enum  ETopBottom { TOP, BOTTOM }
 
typedef Teuchos::RCP< const
MatrixOp > 
D_ptr_t
 

Constructors/initializers.

 MatrixIdentConcatStd ()
 Constructs to uninitialized. More...
 
virtual void initialize (const VectorSpace::space_ptr_t &space_cols, const VectorSpace::space_ptr_t &space_rows, ETopBottom top_or_bottom, value_type alpha, const D_ptr_t &D_ptr, BLAS_Cpp::Transp D_trans)
 Setup with a matrix object. More...
 
virtual void set_uninitialized ()
 Set the matrix to uninitialized. More...
 
virtual const D_ptr_tD_ptr () const
 Return the smart reference counted point to the D matrix. More...
 

Overridden form MatrixIdentConcat

Range1D D_rng () const
 
Range1D I_rng () const
 
value_type alpha () const
 
const MatrixOp & D () const
 
BLAS_Cpp::Transp D_trans () const
 

Overridden from MatrixOp

const VectorSpace & space_cols () const
 
const VectorSpace & space_rows () const
 
MatrixOp & operator= (const MatrixOp &m)
 The default just performs a shallow copy and just copies the underlying smart reference counted pointer. If other behavior is desired then this method must be overridden. More...
 

Additional Inherited Members

- Public Types inherited from AbstractLinAlgPack::MatrixOp
enum  EMatNormType { MAT_NORM_INF, MAT_NORM_2, MAT_NORM_1, MAT_NORM_FORB }
 Type of matrix norm. More...
 
- Public Member Functions inherited from ConstrainedOptPack::MatrixIdentConcat
size_type rows () const
 
size_type cols () const
 
size_type nz () const
 
std::ostream & output (std::ostream &out) const
 
void Vp_StMtV (VectorMutable *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &vs_rhs2, value_type beta) const
 
void Vp_StMtV (VectorMutable *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const
 
- Public Member Functions inherited from AbstractLinAlgPack::MatrixOp
virtual void zero_out ()
 M_lhs = 0 : Zero out the matrix. More...
 
virtual void Mt_S (value_type alpha)
 M_lhs *= alpha : Multiply a matrix by a scalar. More...
 
virtual MatrixOpoperator= (const MatrixOp &mwo_rhs)
 M_lhs = mwo_rhs : Virtual assignment operator. More...
 
virtual mat_mut_ptr_t clone ()
 Clone the non-const matrix object (if supported). More...
 
virtual mat_ptr_t clone () const
 Clone the const matrix object (if supported). More...
 
const MatNorm calc_norm (EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const
 Compute a norm of this matrix. More...
 
virtual mat_ptr_t sub_view (const Range1D &row_rng, const Range1D &col_rng) const
 Create a transient constant sub-matrix view of this matrix (if supported). More...
 
mat_ptr_t sub_view (const index_type &rl, const index_type &ru, const index_type &cl, const index_type &cu) const
 Inlined implementation calls this->sub_view(Range1D(rl,ru),Range1D(cl,cu)). More...
 
virtual mat_ptr_t perm_view (const Permutation *P_row, const index_type row_part[], int num_row_part, const Permutation *P_col, const index_type col_part[], int num_col_part) const
 Create a permuted view: M_perm = P_row' * M * P_col. More...
 
virtual mat_ptr_t perm_view_update (const Permutation *P_row, const index_type row_part[], int num_row_part, const Permutation *P_col, const index_type col_part[], int num_col_part, const mat_ptr_t &perm_view) const
 Reinitialize a permuted view: M_perm = P_row' * M * P_col. More...
 
- Public Member Functions inherited from AbstractLinAlgPack::MatrixBase
virtual ~MatrixBase ()
 Virtual destructor. More...
 
- Protected Member Functions inherited from AbstractLinAlgPack::MatrixOp
virtual bool Mp_StM (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs) const
 mwo_lhs += alpha * op(M_rhs) (BLAS xAXPY). More...
 
virtual bool Mp_StM (value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
 M_lhs += alpha * op(mwo_rhs) (BLAS xAXPY). More...
 
virtual bool Mp_StMtP (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans) const
 mwo_lhs += alpha * op(M_rhs) * op(P_rhs). More...
 
virtual bool Mp_StMtP (value_type alpha, const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans)
 M_lhs += alpha * op(mwo_rhs) * op(P_rhs). More...
 
virtual bool Mp_StPtM (MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, BLAS_Cpp::Transp M_trans) const
 mwo_lhs += alpha * op(P_rhs) * op(M_rhs). More...
 
virtual bool Mp_StPtM (value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans)
 M_lhs += alpha * op(P_rhs) * op(mwo_rhs). More...
 
virtual bool Mp_StPtMtP (MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans) const
 mwo_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2). More...
 
virtual bool Mp_StPtMtP (value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans)
 M_lhs += alpha * op(P_rhs1) * op(mwo_rhs) * op(P_rhs2). More...
 
virtual void Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const =0
 v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs (BLAS xGEMV) More...
 
virtual void Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const
 v_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * v_lhs (BLAS xGEMV) More...
 
virtual void Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const Vector &v_rhs3, value_type beta) const
 v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * v_rhs3 + beta * v_rhs More...
 
virtual void Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const SpVectorSlice &sv_rhs3, value_type beta) const
 v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * v_rhs More...
 
virtual value_type transVtMtV (const Vector &v_rhs1, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3) const
 result = v_rhs1' * op(M_rhs2) * v_rhs3 More...
 
virtual value_type transVtMtV (const SpVectorSlice &sv_rhs1, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3) const
 result = sv_rhs1' * op(M_rhs2) * sv_rhs3 More...
 
virtual void syr2k (BLAS_Cpp::Transp M_trans, value_type alpha, const GenPermMatrixSlice &P1, BLAS_Cpp::Transp P1_trans, const GenPermMatrixSlice &P2, BLAS_Cpp::Transp P2_trans, value_type beta, MatrixSymOp *symwo_lhs) const
 Perform a specialized rank-2k update of a dense symmetric matrix of the form: More...
 
virtual bool Mp_StMtM (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) const
 mwo_lhs = alpha * op(M_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (left) (xGEMM). More...
 
virtual bool Mp_StMtM (MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const
 mwo_lhs = alpha * op(mwo_rhs1) * op(M_rhs2) + beta * mwo_lhs (right) (xGEMM) More...
 
virtual bool Mp_StMtM (value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta)
 M_lhs = alpha * op(mwo_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (left) (xGEMM) More...
 
virtual bool syrk (BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs) const
 Perform a rank-k update of a symmetric matrix of the form: More...
 
virtual bool syrk (const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta)
 Perform a rank-k update of a symmetric matrix of the form: More...
 

Detailed Description

Concrete implementation class for a matrix vertically concatonated with an identity matrix.

Represents an interface for a matrix that represents:

M = [ alpha*op(D) ]    (TOP)
    [     I       ]

or

M = [     I       ]
    [ alpha*op(D) ]   (BOTTOM)

This subclass allows a client to set the representation matrix D.

Definition at line 67 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

Member Typedef Documentation

Definition at line 75 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

Member Enumeration Documentation

Enumerator
TOP 
BOTTOM 

Definition at line 73 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

Constructor & Destructor Documentation

ConstrainedOptPack::MatrixIdentConcatStd::MatrixIdentConcatStd ( )

Constructs to uninitialized.

Definition at line 51 of file ConstrainedOptPack_MatrixIdentConcatStd.cpp.

Member Function Documentation

void ConstrainedOptPack::MatrixIdentConcatStd::initialize ( const VectorSpace::space_ptr_t &  space_cols,
const VectorSpace::space_ptr_t &  space_rows,
ETopBottom  top_or_bottom,
value_type  alpha,
const D_ptr_t D_ptr,
BLAS_Cpp::Transp  D_trans 
)
virtual

Setup with a matrix object.

Parameters
top_or_bottom[in] If TOP then M = [ alpha*op(D); I ] and if BOTTOM then M = [ I; alpha*op(D) ]
alpha[in] Scalar that multiplies op(D)
D_ptr[in] Smart pointer to a matrix object that represents D. The matrix object pointed to must not be altered until this object is no longer in use or this->set_uninitialized() has been called.
D_trans[in] Determines if op(D) = D (no_trans#) or op(D) = D' (trans).

Preconditions:

  • D.get() != NULL (throw std::invalid_argument)
  • space_cols->dim() == rows(op(D)) + cols(op(D)) (throw std::invalid_argument)
  • space_rows->dim() == cols(op(D)) (throw std::invalid_argument)
  • space_cols->sub_space(D_rng)->is_compatible(op(D).space_cols()) (throw std::invalid_argument) See D_rng defined below
  • space_rows->is_compatible(op(D).space_rows()) (throw std::invalid_argument)

Postconditions:

  • this->D_ptr().get() == D_ptr.get()
  • &this->D() == this->D_ptr().get()
  • this->D_trans() == D_trans
  • this->alpha() == alpha
  • this->rows() == rows(op(D)) + cols(op(D))
  • this->cols() == cols(op(D))
  • &this->space_cols() == space_cols.get()
  • &this->space_rows() == space_rows.get()
  • [top_or_bottom == TOP] this->D_rng() = [1,rows(op(D))]
  • [top_or_bottom == TOP] this->I_rng() = [rows(op(D))+1,rows(op(D))+cols(op(D))]
  • [top_or_bottom == BOTTOM] this->D_rng() = [cols(op(D))+1,rows(op(D))+cols(op(D))]
  • [top_or_bottom == BOTTOM] this->I_rng() = [1,cols(op(D))]

Definition at line 56 of file ConstrainedOptPack_MatrixIdentConcatStd.cpp.

void ConstrainedOptPack::MatrixIdentConcatStd::set_uninitialized ( )
virtual

Set the matrix to uninitialized.

Postconditions:

  • this->space_cols() throws an exception
  • this->space_rows() throws an exception
  • this->D_ptr().get() == NULL
  • &this->D() throws an exception
  • this->D_trans() == no_trans
  • this->alpha() == 0.0
  • this->rows() == 0
  • this->cols() == 0
  • [top_or_bottom == TOP] this->D_rng() = [1,rows(op(D))]
  • [top_or_bottom == TOP] this->I_rng() = [rows(op(D))+1,rows(op(D))+cols(op(D))]
  • [top_or_bottom == BOTTOM] this->D_rng() = [cols(op(D))+1,rows(op(D))+cols(op(D))]
  • [top_or_bottom == BOTTOM] this->I_rng() = [1,cols(op(D))]

Definition at line 94 of file ConstrainedOptPack_MatrixIdentConcatStd.cpp.

const MatrixIdentConcatStd::D_ptr_t & ConstrainedOptPack::MatrixIdentConcatStd::D_ptr ( ) const
virtual

Return the smart reference counted point to the D matrix.

If the matrix object D is owned exclusively by this matrix object then this->D_ptr().count() == 1 on return.

Definition at line 106 of file ConstrainedOptPack_MatrixIdentConcatStd.cpp.

Range1D ConstrainedOptPack::MatrixIdentConcatStd::D_rng ( ) const
virtual
Range1D ConstrainedOptPack::MatrixIdentConcatStd::I_rng ( ) const
virtual
value_type ConstrainedOptPack::MatrixIdentConcatStd::alpha ( ) const
virtual
const MatrixOp& ConstrainedOptPack::MatrixIdentConcatStd::D ( ) const
virtual
BLAS_Cpp::Transp ConstrainedOptPack::MatrixIdentConcatStd::D_trans ( ) const
virtual
const VectorSpace& ConstrainedOptPack::MatrixIdentConcatStd::space_cols ( ) const
virtual
const VectorSpace& ConstrainedOptPack::MatrixIdentConcatStd::space_rows ( ) const
virtual
MatrixOp & ConstrainedOptPack::MatrixIdentConcatStd::operator= ( const MatrixOp &  m)

The default just performs a shallow copy and just copies the underlying smart reference counted pointer. If other behavior is desired then this method must be overridden.

Definition at line 150 of file ConstrainedOptPack_MatrixIdentConcatStd.cpp.

void ConstrainedOptPack::MatrixIdentConcatStd::assert_initialized ( ) const
private

Definition at line 158 of file ConstrainedOptPack_MatrixIdentConcatStd.cpp.

Member Data Documentation

const VectorSpace & ConstrainedOptPack::MatrixIdentConcatStd::space_cols
private

Definition at line 182 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

const VectorSpace & ConstrainedOptPack::MatrixIdentConcatStd::space_rows
private

Definition at line 183 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

const MatrixOp & ConstrainedOptPack::MatrixIdentConcatStd::D
private

Definition at line 184 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

Range1D ConstrainedOptPack::MatrixIdentConcatStd::D_rng
private

Definition at line 185 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

Range1D ConstrainedOptPack::MatrixIdentConcatStd::I_rng
private

Definition at line 186 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

value_type ConstrainedOptPack::MatrixIdentConcatStd::alpha_
private

Definition at line 194 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.

BLAS_Cpp::Transp ConstrainedOptPack::MatrixIdentConcatStd::D_trans_
private

Definition at line 195 of file ConstrainedOptPack_MatrixIdentConcatStd.hpp.


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