AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | List of all members
AbstractLinAlgPack::VectorMutableDense Class Reference

DVector "Adaptor" subclass for DenseLinAlgPack::DVectorSlice or DenseLinAlgPack::DVector objects. More...

#include <AbstractLinAlgPack_VectorMutableDense.hpp>

Inheritance diagram for AbstractLinAlgPack::VectorMutableDense:
Inheritance graph
[legend]

Public Types

typedef Teuchos::RCP
< MemMngPack::ReleaseResource
release_resource_ptr_t
 
- Public Types inherited from AbstractLinAlgPack::Vector
typedef Teuchos::RCP< const
Vector
vec_ptr_t
 
typedef Teuchos::RCP
< VectorMutable
vec_mut_ptr_t
 

Public Member Functions

VectorMutableDenseoperator& ()
 Hack. More...
 
- Public Member Functions inherited from AbstractLinAlgPack::VectorMutable
vec_mut_ptr_t sub_view (const index_type &l, const index_type &u)
 Inline member function that simply calls this->sub_view(Range1D(l,u)). More...
 
virtual void zero ()
 Zeros this vector. More...
 
virtual void axpy (value_type alpha, const Vector &x)
 Adds a linear combination of another vector to this vector object. More...
 
vec_ptr_t sub_view (const Range1D &rng) const
 Default implementation calls this->sub_view() (non-const) and then performs an cast to vec_ptr_t. More...
 
- Public Member Functions inherited from AbstractLinAlgPack::Vector
 Vector ()
 
virtual ~Vector ()
 
vec_ptr_t sub_view (const index_type &l, const index_type &u) const
 Inline member function that simply calls this->sub_view(Range1D(l,u)). More...
 
virtual void has_changed () const
 Must be called by any vector subclass that modifies this vector object! More...
 
virtual index_type nz () const
 Return the number of nonzero elements in the vector. More...
 
virtual std::ostream & output (std::ostream &out, bool print_dim=true, bool newline=true, index_type global_offset=0) const
 Virtual output function. More...
 
virtual vec_mut_ptr_t clone () const
 Create a clone of this vector objet. More...
 
virtual value_type norm_1 () const
 One norm. ||v||_1 = sum( |v(i)|, i = 1,,,this->dim() ) More...
 
virtual value_type norm_2 () const
 Two norm. ||v||_2 = sqrt( sum( v(i)^2, i = 1,,,this->dim() ) ) More...
 
virtual value_type norm_inf () const
 Infinity norm. ||v||_inf = max( |v(i)|, i = 1,,,this->dim() ) More...
 
virtual value_type inner_product (const Vector &v) const
 Return the inner product of *this with v. More...
 

Constructors/initializers

 VectorMutableDense (const size_type dim=0)
 Calls this->initialize(dim). More...
 
 VectorMutableDense (DVectorSlice v, const release_resource_ptr_t &v_release)
 Calls this->initialize(v,v_release). More...
 
void initialize (const size_type dim)
 Call this->initialize(v,v_release) with an allocated DenseLinAlgPack::DVector object. More...
 
void initialize (DVectorSlice v, const release_resource_ptr_t &v_release)
 Initialize with a dense vector slice. More...
 

Access

DVectorSlice set_vec ()
 Return the non-const dense vector. More...
 
const DVectorSlice get_vec () const
 Return a const dense vector. More...
 
const release_resource_ptr_tvec_release () const
 Return a RCP<> pointer to the object that will release the associated resource. More...
 

Overriddenn from Vector

const VectorSpacespace () const
 
void apply_op (const RTOpPack::RTOp &op, const size_t num_vecs, const Vector *vecs[], const size_t num_targ_vecs, VectorMutable *targ_vecs[], RTOpPack::ReductTarget *reduct_obj, const index_type first_ele, const index_type sub_dim, const index_type global_offset) const
 
index_type dim () const
 
value_type get_ele (index_type i) const
 
void get_sub_vector (const Range1D &rng, RTOpPack::SubVector *sub_vec) const
 
void free_sub_vector (RTOpPack::SubVector *sub_vec) const
 

Overriddenn from VectorMutable

VectorMutableoperator= (value_type alpha)
 
VectorMutableoperator= (const Vector &v)
 
VectorMutableoperator= (const VectorMutable &v)
 
void set_ele (index_type i, value_type val)
 
vec_mut_ptr_t sub_view (const Range1D &rng)
 
void get_sub_vector (const Range1D &rng, RTOpPack::MutableSubVector *sub_vec)
 
void commit_sub_vector (RTOpPack::MutableSubVector *sub_vec)
 
void set_sub_vector (const RTOpPack::SparseSubVector &sub_vec)
 
void Vp_StMtV (value_type alpha, const GenPermMatrixSlice &P, BLAS_Cpp::Transp P_trans, const Vector &x, value_type beta)
 

Detailed Description

DVector "Adaptor" subclass for DenseLinAlgPack::DVectorSlice or DenseLinAlgPack::DVector objects.

This class can be used either as a view of a DenseLinAlgPack::DVectorSlice object or as a storage type for a DenseLinAlgPack::DVector object.

To create a storage type with the dimension of dim just call the constructor VectorMutableDense(dim) or after construction you can call this->initialize(dim).

To simply create a view of a vector, say v, without ownership just call VectorMutableDense(v(),NULL) or after construction call this->initialize(v(),NULL).

Alternately, this can be given a vector with the responsibility to delete any associated memory by calling this->initialize() with a ReleaseResource object to perform the deallocation.

If this has been initialized by this->initialize(dim) and if the client really needs to get at the DenseLinAlgPack::DVector object itself, then it can be obtained as:

void f( VectorMutableDense* v )
namespace rmp = MemMngPack;
DVector &_v = *dynamic_cast<rmp::ReleaseResource_ref_count_ptr<DVector>&>(*v.vec_release()).ptr;

This is not pretty but it is not supposed to be. Of course the above function will throw an exception if the dynamic_cast<> fails.

Definition at line 84 of file AbstractLinAlgPack_VectorMutableDense.hpp.

Member Typedef Documentation

Definition at line 92 of file AbstractLinAlgPack_VectorMutableDense.hpp.

Constructor & Destructor Documentation

AbstractLinAlgPack::VectorMutableDense::VectorMutableDense ( const size_type  dim = 0)

Calls this->initialize(dim).

Definition at line 65 of file AbstractLinAlgPack_VectorMutableDense.cpp.

AbstractLinAlgPack::VectorMutableDense::VectorMutableDense ( DVectorSlice  v,
const release_resource_ptr_t v_release 
)

Calls this->initialize(v,v_release).

Definition at line 74 of file AbstractLinAlgPack_VectorMutableDense.cpp.

Member Function Documentation

void AbstractLinAlgPack::VectorMutableDense::initialize ( const size_type  dim)

Call this->initialize(v,v_release) with an allocated DenseLinAlgPack::DVector object.

Definition at line 84 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::initialize ( DVectorSlice  v,
const release_resource_ptr_t v_release 
)

Initialize with a dense vector slice.

Definition at line 103 of file AbstractLinAlgPack_VectorMutableDense.cpp.

DVectorSlice AbstractLinAlgPack::VectorMutableDense::set_vec ( )
inline

Return the non-const dense vector.

Note that calling this method will result in the vector implementation being modified. Therefore, no other methods on this object should be called until the DVectorSlice returned from this method is discarded.

Note that the underlying implementation calls this->has_changed() before this method returns.

Definition at line 226 of file AbstractLinAlgPack_VectorMutableDense.hpp.

const DVectorSlice AbstractLinAlgPack::VectorMutableDense::get_vec ( ) const
inline

Return a const dense vector.

Definition at line 234 of file AbstractLinAlgPack_VectorMutableDense.hpp.

const VectorMutableDense::release_resource_ptr_t & AbstractLinAlgPack::VectorMutableDense::vec_release ( ) const
inline

Return a RCP<> pointer to the object that will release the associated resource.

Definition at line 241 of file AbstractLinAlgPack_VectorMutableDense.hpp.

const VectorSpace & AbstractLinAlgPack::VectorMutableDense::space ( ) const
virtual
void AbstractLinAlgPack::VectorMutableDense::apply_op ( const RTOpPack::RTOp &  op,
const size_t  num_vecs,
const Vector vecs[],
const size_t  num_targ_vecs,
VectorMutable targ_vecs[],
RTOpPack::ReductTarget reduct_obj,
const index_type  first_ele,
const index_type  sub_dim,
const index_type  global_offset 
) const
virtual
index_type AbstractLinAlgPack::VectorMutableDense::dim ( ) const
virtual

Reimplemented from AbstractLinAlgPack::Vector.

Definition at line 144 of file AbstractLinAlgPack_VectorMutableDense.cpp.

value_type AbstractLinAlgPack::VectorMutableDense::get_ele ( index_type  i) const
virtual

Reimplemented from AbstractLinAlgPack::Vector.

Definition at line 149 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::get_sub_vector ( const Range1D rng,
RTOpPack::SubVector *  sub_vec 
) const
virtual

Reimplemented from AbstractLinAlgPack::Vector.

Definition at line 154 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::free_sub_vector ( RTOpPack::SubVector *  sub_vec) const
virtual

Reimplemented from AbstractLinAlgPack::Vector.

Definition at line 175 of file AbstractLinAlgPack_VectorMutableDense.cpp.

VectorMutable & AbstractLinAlgPack::VectorMutableDense::operator= ( value_type  alpha)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 183 of file AbstractLinAlgPack_VectorMutableDense.cpp.

VectorMutable & AbstractLinAlgPack::VectorMutableDense::operator= ( const Vector v)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 192 of file AbstractLinAlgPack_VectorMutableDense.cpp.

VectorMutable & AbstractLinAlgPack::VectorMutableDense::operator= ( const VectorMutable v)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 204 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::set_ele ( index_type  i,
value_type  val 
)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 215 of file AbstractLinAlgPack_VectorMutableDense.cpp.

VectorMutableDense::vec_mut_ptr_t AbstractLinAlgPack::VectorMutableDense::sub_view ( const Range1D rng)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 223 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::get_sub_vector ( const Range1D rng,
RTOpPack::MutableSubVector *  sub_vec 
)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 242 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::commit_sub_vector ( RTOpPack::MutableSubVector *  sub_vec)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 263 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::set_sub_vector ( const RTOpPack::SparseSubVector &  sub_vec)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 270 of file AbstractLinAlgPack_VectorMutableDense.cpp.

void AbstractLinAlgPack::VectorMutableDense::Vp_StMtV ( value_type  alpha,
const GenPermMatrixSlice P,
BLAS_Cpp::Transp  P_trans,
const Vector x,
value_type  beta 
)
virtual

Reimplemented from AbstractLinAlgPack::VectorMutable.

Definition at line 276 of file AbstractLinAlgPack_VectorMutableDense.cpp.

VectorMutableDense* AbstractLinAlgPack::VectorMutableDense::operator& ( )
inline

Hack.

Definition at line 202 of file AbstractLinAlgPack_VectorMutableDense.hpp.


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