AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects
Version of the Day
|
Abstract interface for mutable coordinate vectors {abstract}. More...
#include <AbstractLinAlgPack_VectorMutable.hpp>
Virtual methods with default implementations | |
virtual VectorMutable & | operator= (value_type alpha) |
Assign the elements of this vector to a scalar. More... | |
virtual VectorMutable & | operator= (const Vector &v) |
Assign the elements of a vector to this . More... | |
virtual VectorMutable & | operator= (const VectorMutable &v) |
Default implementation calls operator=((const &Vector)v) . More... | |
virtual void | set_ele (index_type i, value_type val) |
Set a specific element of a vector. More... | |
virtual vec_mut_ptr_t | sub_view (const Range1D &rng) |
Create a mutable abstract view of a vector object. More... | |
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... | |
virtual void | get_sub_vector (const Range1D &rng, RTOpPack::MutableSubVector *sub_vec) |
Get a mutable explicit view of a sub-vector. More... | |
virtual void | commit_sub_vector (RTOpPack::MutableSubVector *sub_vec) |
Free a mutable explicit view of a sub-vector. More... | |
virtual void | set_sub_vector (const RTOpPack::SparseSubVector &sub_vec) |
Set a specific sub-vector. More... | |
virtual void | Vp_StMtV (value_type alpha, const GenPermMatrixSlice &P, BLAS_Cpp::Transp P_trans, const Vector &x, value_type beta) |
Perform a gather or scatter operation with a vector. More... | |
Overridden from Vector | |
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... | |
Additional Inherited Members | |
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 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 const VectorSpace & | space () const =0 |
Return the vector space that this vector belongs to. More... | |
virtual index_type | dim () const |
Return the dimension of this vector. 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 | get_ele (index_type i) const |
Fetch an element in the vector. 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... | |
virtual void | get_sub_vector (const Range1D &rng, RTOpPack::SubVector *sub_vec) const |
Get a non-mutable explicit view of a sub-vector. More... | |
virtual void | free_sub_vector (RTOpPack::SubVector *sub_vec) const |
Free an explicit view of a sub-vector. More... | |
Protected Member Functions inherited from AbstractLinAlgPack::Vector | |
virtual 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 =0 |
Apply a reduction/transformation,operation over a set of vectors: op(op(v[0]...v[nv-1],z[0]...z[nz-1]),(*reduct_obj)) -> z[0]...z[nz-1],(*reduct_obj) . More... | |
virtual void | finalize_apply_op (const size_t num_targ_vecs, VectorMutable **targ_vecs) const |
This method usually needs to be called by subclasses at the end of the apply_op() method implementation to insure that has_changed() is called on the transformed vector objects. More... | |
Abstract interface for mutable coordinate vectors {abstract}.
Objects of this type can act as a target vector of a transformation operation. Similarly to Vector
this interface contains very few (only one extra) pure virtual methods that must be overridden. However, more efficient and more general implementations will choose to override more methods.
In addition to being able to create non-mutable (const
) abstract sub-views of a vector object thorugh the Vector
interface, this interface allows the creation of mutable (non-const
) sub-views using sub_view()
. Also, in addition to being able to extract an explicit non-mutable view of some (small?) sub-set of elements, this interface allows a client to either extract a explicit mutable sub-views using get_sub_vector()
or to set sub-vectors using set_sub_vector()
. As much as possible, abstract views should be preferred (i.e. sub_view()
) over explict views (i.e. get_sub_vector() and set_sub_vector()).
There are only two pure virtual methods that a concreate VectorMutable
subclass must override. The space()
and apply_op()
methods from the Vector
base class inteface must be defined.
The non-mutable (const
) sub_view(...)
method from the Vector
interface has a default implementation defined here that will be adequate for most subclasses.
Definition at line 72 of file AbstractLinAlgPack_VectorMutable.hpp.
|
virtual |
Assign the elements of this
vector to a scalar.
The default implementation of this function uses a transforamtion operator class (see RTOp_TOp_assign_scalar.h) and calls this->apply_op()
.
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorMutableDense.
Definition at line 92 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Assign the elements of a vector to this
.
The default implementation of this function uses a transforamtion operator class (see RTOp_TOp_assign_vectors.h) and calls this->apply_op()
.
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorMutableDense.
Definition at line 100 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Default implementation calls operator=((const &Vector)v)
.
Reimplemented in AbstractLinAlgPack::VectorMutableDense.
Definition at line 110 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Set a specific element of a vector.
Preconditions:
1 <= i <= this->dim()
(throw std::out_of_range
) Postconditions:
this->get(i) == val
The default implementation uses a transforamtion operator class (see RTOp_TOp_set_ele.h) and calls this->apply_op()
.
i | [in] Index of the element value to set. |
val | [in] Value of the element to set. |
Reimplemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorMutableSubView.
Definition at line 115 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Create a mutable abstract view of a vector object.
This is only a transient view of a sub-vector that is to be immediately used and then released by RCP<>
. This function is declared as non-constant because the object returned has the capacity to alter this
object.
The compatibility of sub-views goes along with the compatibility of sub-spaces (see VectorSpace). For example, given the vector objects where x.space().is_compatible(y.space()) == true
then if x.space().sub_space(rng1)->is_compatible(*y.space().sub_space(rng2)) == true
then the sub-vector views *x.sub_view(rng1)
and *y.sub_view(rng2)
should be compatible and can be combined in vector operations.
Preconditions:
rng.in_range(this->dim()) == true
(throw std::out_of_range
) rng | [in] The range of the elements to extract the sub-vector view. |
return.get() == NULL
. In most applications, only specific views are every required. The default implementation uses the subclass VectorSubView to represent any arbitrary sub-view but this can be inefficient if the sub-view is very small compared this this full vector space but not necessarily. Note that the underlying vector this
is not guarrenteed to show the changes made the sub-view *return .get()
until the smart reference counted pointer return
is destroyed. Reimplemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorMutableSubView.
Definition at line 126 of file AbstractLinAlgPack_VectorMutable.cpp.
|
inline |
Inline member function that simply calls this->sub_view(Range1D(l,u))
.
Definition at line 307 of file AbstractLinAlgPack_VectorMutable.hpp.
|
virtual |
Zeros this vector.
Calls operator=(0.0)
.
Definition at line 148 of file AbstractLinAlgPack_VectorMutable.cpp.
Adds a linear combination of another vector to this vector object.
Calls this->apply_op()
with an operator class (see RTOp_TOp_axpy.h).
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked.
Definition at line 153 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Get a mutable explicit view of a sub-vector.
This is only a transient view of a sub-vector that is to be immediately used and then released with a call to release_sub_vector()
.
Note that calling this operation might require some internal allocations and temporary memory. Therefore, it is critical that this->release_sub_vector(sub_vec)
is called to clean up memory and avoid memory leaks after the sub-vector is used.
If this->get_sub_vector(...,sub_vec)
was previously called on sub_vec
then it may be possible to reuse this memory if it is sufficiently sized. The user is encouraged to make multiple calls to this->get_sub_vector(...,sub_vec)
before this->release_sub_vector(sub_vec)
to finally clean up all of the memory. Of course the same sub_vec
object must be passed to the same vector object for this to work correctly.
Changes to the underlying sub-vector are not guarrenteed to become permanent until this->get_sub_vector(...,sub_vec)
or this->commit_sub_vector(sub_vec)
is called.
Preconditions:
!rng.full_range()
] (rng.ubound() <= this->dim()) == true
(throw std::out_of_range
) This method has a default implementation based on a vector reduction operator class (see RTOp_ROp_get_sub_vector.h) and calls apply_op()
. Note that the footprint of the reduction object (both internal and external state) will be O(rng.size()
). For serial applications this is faily adequate and will not be a major performance penalty. For parallel applications, this will be a terrible implementation and must be overridden if rng.size()
is large at all. If a subclass does override this method, it must also override release_sub_vector()
which has a default implementation which is a companion to this method's default implementation.
rng | [in] The range of the elements to extract the sub-vector view. |
sub_vec | [in/out] Mutable view of the sub-vector. Prior to the first call RTOp_mutable_sub_vector_null(sub_vec) must have been called for the correct behavior. Technically *sub_vec owns the memory but this memory can be freed only by calling this->commit_sub_vector(sub_vec) . |
Reimplemented in AbstractLinAlgPack::VectorMutableDense, and AbstractLinAlgPack::VectorMutableSubView.
Definition at line 161 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Free a mutable explicit view of a sub-vector.
The sub-vector view must have been allocated by this->get_sub_vector()
first.
This method has a default implementation which is a companion to the default implementation for get_sub_vector(...)
. If get_sub_vector(...)
is overridden by a subclass then this method must be overridden also!
sub_vec | [in/out] The memory refered to by sub_vec->values and sub_vec->indices will be released if it was allocated and *sub_vec will be zeroed out using RTOp_mutable_sub_vector_null(sub_vec) . |
Reimplemented in AbstractLinAlgPack::VectorMutableDense, and AbstractLinAlgPack::VectorMutableSubView.
Definition at line 178 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Set a specific sub-vector.
After this function returns, the corresponding elements in this
vector object will be set equal to those in the input vector (the post conditions are obvious).
Preconditions:
sub_vec.global_offset + sub_dim <= this->dim()
(throw std::out_of_range
) The default implementation of this operation uses a transformation operator class (see RTOp_TOp_set_sub_vector.h) and calls apply_op()
. Be forewarned however, that the operator objects state data (both internal and external) will be O(sub_vec.sub_nz
). For serial applications, this is entirely adequate. For parallel applications this will be very bad!
sub_vec | [in] Represents the elements in the subvector to be set. |
Reimplemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorMutableSubView.
Definition at line 191 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Perform a gather or scatter operation with a vector.
this = alpha * op(P) * x + beta * this
The default implementation is based on a transformation or reduction operator (depending if a gather or scatter is being performed).
Reimplemented in AbstractLinAlgPack::VectorMutableDense.
Definition at line 221 of file AbstractLinAlgPack_VectorMutable.cpp.
|
virtual |
Default implementation calls this->sub_view()
(non-const
) and then performs an cast to vec_ptr_t
.
This function override is actually needed here for another reason. Without, the override, the non-const version defined in this interface hides the const version defined in Vector.
Reimplemented from AbstractLinAlgPack::Vector.
Reimplemented in AbstractLinAlgPack::VectorMutableSubView.
Definition at line 242 of file AbstractLinAlgPack_VectorMutable.cpp.