MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_MultiVector.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef ALAP_MULTI_VECTOR_H
43 #define ALAP_MULTI_VECTOR_H
44 
46 #include "RTOpPack_RTOpT.hpp"
47 #include "Teuchos_RCP.hpp"
48 
49 namespace AbstractLinAlgPack {
50 
52 enum EApplyBy {
55 };
56 
62 void apply_op(
63  EApplyBy apply_by
64  ,const RTOpPack::RTOp &primary_op
65  ,const size_t num_multi_vecs
66  ,const MultiVector* multi_vecs[]
67  ,const size_t num_targ_multi_vecs
68  ,MultiVectorMutable* targ_multi_vecs[]
69  ,RTOpPack::ReductTarget* reduct_objs[] = NULL
70  ,const index_type primary_first_ele = 1
71  ,const index_type primary_sub_dim = 0
72  ,const index_type primary_global_offset = 0
73  ,const index_type secondary_first_ele = 1
74  ,const index_type secondary_sub_dim = 0
75  );
76 
82 void apply_op(
83  EApplyBy apply_by
84  ,const RTOpPack::RTOp &primary_op
85  ,const RTOpPack::RTOp &secondary_op
86  ,const size_t num_multi_vecs
87  ,const MultiVector* multi_vecs[]
88  ,const size_t num_targ_multi_vecs
89  ,MultiVectorMutable* targ_multi_vecs[]
90  ,RTOpPack::ReductTarget *reduct_obj
91  ,const index_type primary_first_ele = 1
92  ,const index_type primary_sub_dim = 0
93  ,const index_type primary_global_offset = 0
94  ,const index_type secondary_first_ele = 1
95  ,const index_type secondary_sub_dim = 0
96  );
97 
170 class MultiVector : virtual public MatrixOp {
171 public:
172 
174  using MatrixOp::clone;
176  using MatrixOp::Mp_StMtM;
177 
179  typedef int access_by_t;
181  enum {
182  ROW_ACCESS = 0x1
183  ,COL_ACCESS = 0x2
184  ,DIAG_ACCESS = 0x4
185  };
190 
193 
195  friend void apply_op(
196  EApplyBy apply_by
197  ,const RTOpPack::RTOp &primary_op
198  ,const size_t num_multi_vecs
199  ,const MultiVector* multi_vecs[]
200  ,const size_t num_targ_multi_vecs
201  ,MultiVectorMutable* targ_multi_vecs[]
202  ,RTOpPack::ReductTarget* reduct_objs[]
203  ,const index_type primary_first_ele
204  ,const index_type primary_sub_dim
205  ,const index_type primary_global_offset
206  ,const index_type secondary_first_ele
207  ,const index_type secondary_sub_dim
208  );
210  friend void apply_op(
211  EApplyBy apply_by
212  ,const RTOpPack::RTOp &primary_op
213  ,const RTOpPack::RTOp &secondary_op
214  ,const size_t num_multi_vecs
215  ,const MultiVector* multi_vecs[]
216  ,const size_t num_targ_multi_vecs
217  ,MultiVectorMutable* targ_multi_vecs[]
218  ,RTOpPack::ReductTarget *reduct_obj
219  ,const index_type primary_first_ele
220  ,const index_type primary_sub_dim
221  ,const index_type primary_global_offset
222  ,const index_type secondary_first_ele
223  ,const index_type secondary_sub_dim
224  );
225 
227 
230 
235  virtual multi_vec_ptr_t mv_clone() const;
236 
238 
241 
248  virtual access_by_t access_by() const = 0;
249 
257  virtual vec_ptr_t col(index_type j) const = 0;
265  virtual vec_ptr_t row(index_type i) const = 0;
272  virtual vec_ptr_t diag(int k) const = 0;
273 
275 
278 
286  virtual multi_vec_ptr_t mv_sub_view(const Range1D& row_rng, const Range1D& col_rng) const;
287 
291  const index_type& rl, const index_type& ru
292  ,const index_type& cl, const index_type& cu
293  ) const;
294 
296 
297 protected:
298 
301 
315  virtual void apply_op(
316  EApplyBy apply_by, const RTOpPack::RTOp& primary_op
317  ,const size_t num_multi_vecs, const MultiVector* multi_vecs[]
318  ,const size_t num_targ_multi_vecs, MultiVectorMutable* targ_multi_vecs[]
319  ,RTOpPack::ReductTarget* reduct_objs[]
320  ,const index_type primary_first_ele, const index_type primary_sub_dim, const index_type primary_global_offset
321  ,const index_type secondary_first_ele, const index_type secondary_sub_dim
322  ) const;
323 
337  virtual void apply_op(
338  EApplyBy apply_by, const RTOpPack::RTOp& primary_op, const RTOpPack::RTOp& secondary_op
339  ,const size_t num_multi_vecs, const MultiVector* multi_vecs[]
340  ,const size_t num_targ_multi_vecs, MultiVectorMutable* targ_multi_vecs[]
341  ,RTOpPack::ReductTarget* reduct_obj
342  ,const index_type primary_first_ele, const index_type primary_sub_dim, const index_type primary_global_offset
343  ,const index_type secondary_first_ele, const index_type secondary_sub_dim
344  ) const;
345 
347 
348 public:
349 
352 
355  mat_ptr_t clone() const;
356 
359  mat_ptr_t sub_view(const Range1D& row_rng, const Range1D& col_rng) const;
360 
374  bool Mp_StMtM(
375  MatrixOp* mwo_lhs, value_type alpha
376  ,const MatrixOp& mwo_rhs1, BLAS_Cpp::Transp trans_rhs1
377  ,BLAS_Cpp::Transp trans_rhs2
378  ,value_type beta
379  ) const;
380 
394  bool Mp_StMtM(
395  MatrixOp* mwo_lhs, value_type alpha
396  ,BLAS_Cpp::Transp trans_rhs1
397  ,const MatrixOp& mwo_rhs2, BLAS_Cpp::Transp trans_rhs2
398  ,value_type beta
399  ) const;
400 
402 
403 private:
404 
405 #ifdef DOXYGEN_COMPILE
409 #endif
410 
411 }; // end class MultiVector
412 
413 // //////////////////////////////////////////////////
414 // Inlined functions
415 
416 inline
419  const index_type& rl, const index_type& ru
420  ,const index_type& cl, const index_type& cu
421  ) const
422 {
423  return this->mv_sub_view(Range1D(rl,ru),Range1D(cl,cu));
424 }
425 
426 } // end namespace AbstractLinAlgPack
427 
428 #endif // ALAP_MULTI_VECTOR_H
mat_ptr_t clone() const
Returns this->mv_clone().
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
virtual access_by_t access_by() const =0
Return a bit field for the types of access that are the most convenient.
mat_ptr_t sub_view(const Range1D &row_rng, const Range1D &col_rng) const
Returns this->mv_sub_view(row_rng,col_rng) casted to a MatrixOp.
virtual mat_mut_ptr_t clone()
Clone the non-const matrix object (if supported).
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
Provides a specialized implementation for mwo_rhs1 of type MatrixSymDiag.
. One-based subregion index range class.
RTOpT< RTOp_value_type > RTOp
Teuchos::RCP< const MultiVector > multi_vec_ptr_t
friend void Mp_StMtM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta)
virtual vec_ptr_t diag(int k) const =0
Get a non-mutable diagonal vector.
Base class for all matrices that support basic matrix operations.
Interface for a collection of non-mutable vectors (multi-vector, matrix).
virtual multi_vec_ptr_t mv_sub_view(const Range1D &row_rng, const Range1D &col_rng) const
Returns a sub-view of the multi vector.
Interface for a collection of mutable vectors (multi-vector, matrix).
void apply_op(EApplyBy apply_by, const RTOpPack::RTOp &primary_op, const size_t num_multi_vecs, const MultiVector *multi_vecs[], const size_t num_targ_multi_vecs, MultiVectorMutable *targ_multi_vecs[], RTOpPack::ReductTarget *reduct_objs[]=NULL, const index_type primary_first_ele=1, const index_type primary_sub_dim=0, const index_type primary_global_offset=0, const index_type secondary_first_ele=1, const index_type secondary_sub_dim=0)
Apply a reduction/transformation operator column by column and return an array of the reduction objec...
friend void apply_op(EApplyBy apply_by, const RTOpPack::RTOp &primary_op, const size_t num_multi_vecs, const MultiVector *multi_vecs[], const size_t num_targ_multi_vecs, MultiVectorMutable *targ_multi_vecs[], RTOpPack::ReductTarget *reduct_objs[], const index_type primary_first_ele, const index_type primary_sub_dim, const index_type primary_global_offset, const index_type secondary_first_ele, const index_type secondary_sub_dim)
virtual multi_vec_ptr_t mv_clone() const
Clone the non-const multi-vector object.
virtual vec_ptr_t col(index_type j) const =0
Get a non-mutable column vector.
virtual vec_ptr_t row(index_type i) const =0
Get a non-mutable row vector.
Transp
TRANS.