AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects
Version of the Day
|
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}. More...
#include <AbstractLinAlgPack_Vector.hpp>
Public Types | |
typedef Teuchos::RCP< const Vector > | vec_ptr_t |
typedef Teuchos::RCP < VectorMutable > | vec_mut_ptr_t |
Public Member Functions | |
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... | |
Friends | |
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) |
Pure virtual methods (must be overridden by subclass) | |
virtual const VectorSpace & | space () const =0 |
Return the vector space that this vector belongs to. More... | |
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... | |
Miscellaneous virtual methods with default implementations | |
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 vec_ptr_t | sub_view (const Range1D &rng) const |
Create an abstract view of a vector object . More... | |
Vector norms | |
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... | |
Inner product | |
virtual value_type | inner_product (const Vector &v) const |
Return the inner product of *this with v . More... | |
Explicit sub-vector access | |
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 helper functions | |
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 immutable, finite dimensional, coordinate vectors {abstract}.
This interface contains a mimimal set of operations. The main feature of this interface is the operation apply_op()
. Almost every standard (i.e. BLAS) and non-standard element-wise operation that can be performed on a set of coordinate vectors without changing (mutating) the vectors can be performed through reduction operators. More standard vector operations could be included in this interface and allow for specialized implementations but in general, assuming the sub-vectors are large enough, such implementations would not be significantly faster than those implemented through reduction/transformation operators. There are some operations however that can not always be efficiently with reduction/transforamtion operators and a few of these important methods are included in this interface. The apply_op()
method allows to client to specify a sub-set of the vector elements to include in reduction/transformation operation. This greatly increases the generality of this vector interface as vector objects can be used as sub objects in larger composite vectors and sub-views of a vector can be created.
This interface allows clients to create sub-views of a vector using sub_view()
that in turn are fully functional Vector
objects. This functionality is supported by default by using a default vector subclass VectorSubView
which in turn calls apply_op()
but the client need not ever worry about how this is done.
This interface also allows a client to extract a sub-set of elements in an explicit form as an RTOpPack::SubVector object using the method
get_sub_vector()
. In general, this is very bad thing to do and should be avoided at all costs. However, there are some applications where this is needed and therefore it is supported. The default implementation of this method uses a reduction/transformation operator with apply_op()
in order to extract the needed elements.
In order to create a concreate subclass of this interface, only two methods must be overridden. The
space()
method must be overridden which in turn requires defining a concreate VectorSpace
class (which has only two pure virtual methods). And, as mentioned above, the apply_op()
method must be overridden as well.
The fact that this interface defines
space()
which returns a VectorSpace
object (which in turn can create mutable vectors) implies that for every possible vector object, it is possible to associate with it a mutable vector object that can be the target of transformation operations. This is not a serious limitation. For any application area, mutable vectors should be able to defined and should be usable with the non-mutable vectors.
This interface includes methods for the common vector norms:
norm_1()
, norm_2()
, norm_inf()
. The default implementation of this class uses reduction operator classes (See RTOp_ROp_norms.h) and caches away the values of the norms that are computed since it is common that the norms will be accessed many times before a vector is changed. The operations in any subclass that modifies the underlying vector must call the method this->
has_changed() in order to alert this implementation that the norms are no longer valid.
ToDo: Add example code!
Definition at line 187 of file AbstractLinAlgPack_Vector.hpp.
typedef Teuchos::RCP<const Vector> AbstractLinAlgPack::Vector::vec_ptr_t |
Definition at line 191 of file AbstractLinAlgPack_Vector.hpp.
Definition at line 193 of file AbstractLinAlgPack_Vector.hpp.
AbstractLinAlgPack::Vector::Vector | ( | ) |
Definition at line 123 of file AbstractLinAlgPack_Vector.cpp.
|
inlinevirtual |
Definition at line 212 of file AbstractLinAlgPack_Vector.hpp.
|
pure virtual |
Return the vector space that this vector belongs to.
Note that the vectors space object returned is specifically bound to this vector object. The vector space object returned should only be considered to be transient and may become invalid if this
is modified in some way (though some other interface).
Implemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorSubView.
|
protectedpure virtual |
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)
.
The vector this
that this method is called on is assumed to be one of the vectors in v[0]...v[nv-1],z[0]...z[nz-1]
. This method is not called directly by the client but is instead TSFCore::applyOp()
.
See the documentation for the method AbstractLinAlgPack::apply_op()
for a description of what this method does.
Implemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorSubView.
|
virtual |
Return the dimension of this vector.
It is allowed for a vector to return a dimension of 0
in which case the vector should be considered uninitialized in which the client should not call any of the member functions (except space()). The default implementation returns this->space().dim()
.
Reimplemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorSubView.
Definition at line 128 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Return the number of nonzero elements in the vector.
The default implementation just uses a reduction operator with the apply_op()
method (See RTOp_ROp_num_nonzeros.h).
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked.
Definition at line 134 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Virtual output function.
The default implementation just uses get_sub_vector(...)
to convert to a dense vector and then prints this.
ToDo: Finish documentation!
out | [in/out] Receives the output. |
print_dim | [in] (default = true) If true, then the dimension is printed first followed by a newline. |
newline | [in] (default = true) If true, then a newline is printed after the last element is printed. |
global_offset | [in] (default = 0) The offset added to the vector element indexes when they are printed. |
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked.
Definition at line 150 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Create a clone of this vector objet.
The vector object returned in a smart reference counted pointer to a functional copy of the current vector object. The vector object this
and the vector returned by this method can be modified independently.
The default implementation of this function calls on this->space().create_member()
and then copies over the elements from this
using operator=()
.
Definition at line 166 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Fetch an element in the vector.
Preconditions:
1 <= i <= this->dim()
(throw std::out_of_range
) The default implementation uses a C reduction operator class (See RTOp_ROp_get_ele.h C).
i | [in] Index of the element value to get. |
Reimplemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorSubView.
Definition at line 175 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Create an 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<>
.
It is important to understand what the minimum postconditions are for the sub vector objects returned from this method. If two vector objects x
and y
are compatible (possibly of different types) it is assumed that *x.sub_view(rng)
and *y.sub_view(rng)
will also be compatible vector objects no mater what range rng
represents. However, if i1 < i2 < i3 < i4
with i2-i1 == i4-i3
, then in general, one can not expect that the vector objects *x.sub_view(Range1D(i2,i1))
and *x.sub_view(Range1D(i4,i5))
will be compatible objects. For some vector implementaions they may be (i.e. serial vectors) but for others they most certainly will not be (i.e. parallel vectors). This limitation must be kept in mind by all vector subclass implementors and vector interface clients.
Preconditions:
!rng.full_range()
] (rng.ubound() <= this->dim()) == true
(throw std::out_of_range
) Postconditions:
[return.get() != NULL] return->get_ele(i) == this->get_ele(i+rng.lbound()-1)
, for i = 1...rng.size()
. rng | [in] The range of the elements to extract the sub-vector view. It is allowed for rng.full_range() == true in which case it implicitly treated as rng = [1,this->dim()] . |
return->get() == NULL
for some selections of rng
. Only some rng
ranges may be allowed but they will be appropriate for the application at hand. However, a very good implementation should be able to accommodate any valid rng
that meets the basic preconditions. 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. Reimplemented in AbstractLinAlgPack::VectorMutable, AbstractLinAlgPack::VectorSubView, and AbstractLinAlgPack::VectorMutableSubView.
Definition at line 241 of file AbstractLinAlgPack_Vector.cpp.
|
inline |
Inline member function that simply calls this->sub_view(Range1D(l,u))
.
Definition at line 510 of file AbstractLinAlgPack_Vector.hpp.
|
virtual |
One norm. ||v||_1 = sum( |v(i)|, i = 1,,,this->dim() )
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked.
Definition at line 187 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Two norm. ||v||_2 = sqrt( sum( v(i)^2, i = 1,,,this->dim() ) )
Definition at line 203 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Infinity norm. ||v||_inf = max( |v(i)|, i = 1,,,this->dim() )
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked.
Definition at line 219 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Return the inner product of *this
with v
.
this->space().inner_prod()->inner_prod(*this,v)
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked.
Definition at line 235 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Get a non-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.
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] View of the sub-vector. Prior to the first call RTOp_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->free_sub_vector(sub_vec) . |
Reimplemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorSubView.
Definition at line 261 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Free an 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_sub_vector_null(sub_vec) . |
Reimplemented in AbstractLinAlgPack::VectorMutableDense, AbstractLinAlgPack::VectorMutableBlocked, and AbstractLinAlgPack::VectorSubView.
Definition at line 297 of file AbstractLinAlgPack_Vector.cpp.
|
virtual |
Must be called by any vector subclass that modifies this vector object!
The way to use this method by subclasses is to call it when ever there is a chance that the vector may have changed. Therefore, this method should be called in every non-const member function in every subclass. This is a little bit of a pain but overall this is a good design in that it allows for efficient cacheing of information for multiple retreval. For example, if the subclass SomeVector
has cashed data and has a method SomeVector::foo()
may modify the vector then SomeVector
should override the method has_changed()
and its implementation should look someting likde like this!
void SomeVector::has_changed() { BaseClass::has_changed(); // Called on most direct subclass that // has overridden this method as well. ... // Reinitialize your own cached data to uninitialized! }
Reimplemented in AbstractLinAlgPack::VectorMutableBlocked.
Definition at line 306 of file AbstractLinAlgPack_Vector.cpp.
|
protectedvirtual |
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.
Definition at line 316 of file AbstractLinAlgPack_Vector.cpp.
|
friend |
The logical vector passed to the op.apply_op(...)
method is:
v(k + global_offset) = this->get_ele(first_ele + k - 1) , for k = 1 ... sub_dim
where v
represents any one of the input or input/output vectors. The situation where first_ele == 1
and global_offset > 1
corresponds to the case where the vectors represent consitituent vectors in a larger aggregrate vector. The situation where first_ele > 1
and global_offset == 0
is for when a sub-view of the vectors are being treated as full vectors. Other combinations of these arguments are also possible.
Preconditions:
num_vecs > 0
] vecs[k]->space().isCompatible(this->space()) == true
, for k = 0...num_vecs-1
(throw VectorSpace::IncompatibleVectorSpaces
) num_targ_vecs > 0
] targ_vecs[k]->space().isCompatible(this->space()) == true
, for k = 0...num_targ_vecs-1
(throw VectorSpace::IncompatibleVectorSpaces
) 1 <= first_ele <= this->dim()
(throw std::out_of_range
) global_offset >= 0
(throw std::invalid_argument
) sub_dim - (first_ele - 1) <= this->dim()
(throw std::length_error
). op | [in] Reduction/transformation operator to apply over each sub-vector and assemble the intermediate targets into reduct_obj (if reduct_obj != NULL ). |
num_vecs | [in] Number of nonmutable vectors in vecs[] . If vecs==NULL then this argument is ignored but should be set to zero. |
vecs | [in] Array (length num_vecs ) of a set of pointers to nonmutable vectors to include in the operation. The order of these vectors is significant to op . |
num_targ_vecs | [in] Number of mutable vectors in targ_vecs[] . If targ_vecs==NULL then this argument is ignored but should be set to zero. |
targ_vecs | [in] Array (length num_targ_vecs ) of a set of pointers to mutable vectors to include in the operation. The order of these vectors is significant to op . If targ_vecs==NULL then op is called with no mutable vectors. |
reduct_obj | [in/out] Target object of the reduction operation. This object must have been created by the op.reduct_obj_create_raw(&reduct_obj) function first. The reduction operation will be added to (*reduct_obj) if (*reduct_obj) has already been through a reduction. By allowing the info in (*reduct_obj) to be added to the reduction over all of these vectors, the reduction operation can be accumulated over a set of abstract vectors which can be useful for implementing composite vectors for instance. If op.get_reduct_type_num_entries(...) returns num_values == 0 , num_indexes == 0 and num_chars == 0 then reduct_obj should be set to NULL and no reduction will be performed. |
first_ele | [in] (default = 1) The index of the first element in this to be included. |
sub_dim | [in] (default = 0) The number of elements in these vectors to include in the reduction/transformation operation. The value of sub_dim == 0 means to include all available elements. |
global_offset | [in] (default = 0) The offset applied to the included vector elements. |