MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_VectorMutableThyra.cpp
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 #include <assert.h>
43 
48 #include "RTOpPack_RTOpSubRangeDecorator.hpp"
49 #include "Teuchos_Assert.hpp"
50 #include "Teuchos_Workspace.hpp"
51 #include "Teuchos_dyn_cast.hpp"
52 #include "Teuchos_as.hpp"
53 
54 namespace AbstractLinAlgPack {
55 
56 // Constructors / Initializers
57 
59 {}
60 
62  const Teuchos::RCP<Thyra::VectorBase<value_type> >& thyra_vec
63  )
64 {
65  this->initialize(thyra_vec);
66 }
67 
69  const Teuchos::RCP<Thyra::VectorBase<value_type> >& thyra_vec
70  )
71 {
72  namespace mmp = MemMngPack;
74  thyra_vec.get()==NULL, std::invalid_argument
75  ,"VectorMutableThyra::initialize(thyra_vec): Error!"
76  );
77  thyra_vec_ = thyra_vec;
78  space_.initialize(thyra_vec->space());
79  this->has_changed();
80 }
81 
84 {
85  Teuchos::RCP<Thyra::VectorBase<value_type> > tmp_thyra_vec = thyra_vec_;
86  thyra_vec_ = Teuchos::null;
88  this->has_changed();
89  return tmp_thyra_vec;
90 }
91 
92 // Methods overridden from Vector
93 
94 const VectorSpace&
96 {
97  return space_;
98 }
99 
101  const RTOpPack::RTOp &op
102  ,const size_t num_vecs
103  ,const Vector* vecs[]
104  ,const size_t num_targ_vecs
105  ,VectorMutable* targ_vecs[]
106  ,RTOpPack::ReductTarget *reduct_obj
107  ,const index_type first_ele
108  ,const index_type sub_dim
109  ,const index_type global_offset
110  ) const
111 {
112  using Teuchos::as;
113  using Teuchos::dyn_cast;
114  using Teuchos::Workspace;
115  using Teuchos::rcpFromRef;
116  using Teuchos::RCP;
117  using Teuchos::Ptr;
119 
120  // If these are in-core vectors then just let a default implementation
121  // take care of this.
122  if( space_.is_in_core() ) {
123  this->apply_op_serial(
124  op,num_vecs,vecs,num_targ_vecs,targ_vecs,reduct_obj
125  ,first_ele,sub_dim,global_offset
126  );
127  return;
128  }
129 
130  // Convert the non-mutable vectors into non-mutable Thyra vectors
131  Workspace< Teuchos::RCP<const Thyra::VectorBase<value_type> > > thyra_vecs_sptr(wss,num_vecs);
132  Workspace<Ptr<const Thyra::VectorBase<value_type> > > thyra_vecs(wss,num_vecs);
133  for(int k = 0; k < as<int>(num_vecs); ++k ) {
134  get_thyra_vector( space_, *vecs[k], &thyra_vecs_sptr[k] );
135  thyra_vecs[k] = thyra_vecs_sptr[k].ptr();
136  }
137 
138  // Convert the mutable vetors into mutable Thyra vectors
139  Workspace< Teuchos::RCP<Thyra::VectorBase<value_type> > > targ_thyra_vecs_sptr(wss,num_targ_vecs);
140  Workspace<Ptr<Thyra::VectorBase<value_type> > > targ_thyra_vecs(wss,num_targ_vecs);
141  for(int k = 0; k < as<int>(num_targ_vecs); ++k ) {
142  get_thyra_vector( space_, targ_vecs[k], &targ_thyra_vecs_sptr[k] );
143  targ_thyra_vecs[k] = targ_thyra_vecs_sptr[k].ptr();
144  }
145 
146  // Call the Thyra::apply_op(...)
147  RTOpPack::RTOpSubRangeDecorator<value_type>
148  subRangeOp(rcpFromRef(op), first_ele-1, sub_dim==0 ? -1 : sub_dim);
149  Thyra::applyOp<value_type>(subRangeOp, thyra_vecs(), targ_thyra_vecs(),
150  Teuchos::ptr(reduct_obj), global_offset);
151 
152  // Free/commit the Thyra vector views
153  for(int k = 0; k < as<int>(num_vecs); ++k ) {
154  free_thyra_vector( space_, *vecs[k], &thyra_vecs_sptr[k] );
155  }
156  for(int k = 0; k < as<int>(num_targ_vecs); ++k ) {
157  commit_thyra_vector( space_, targ_vecs[k], &targ_thyra_vecs_sptr[k] );
158  }
159 
160 }
161 
162 
164 {
165  return space_.dim();
166 }
167 
169  const Range1D& rng, RTOpPack::SubVector* sub_vec
170  ) const
171 {
172  RTOpPack::ConstSubVectorView<RTOp_value_type> _sub_vec = *sub_vec;
173  thyra_vec_->acquireDetachedView(convert(rng),&_sub_vec);
174  *sub_vec = _sub_vec;
175 }
176 
178  RTOpPack::SubVector* sub_vec
179  ) const
180 {
181  RTOpPack::ConstSubVectorView<RTOp_value_type> _sub_vec = *sub_vec;
182  thyra_vec_->releaseDetachedView(&_sub_vec);
183  *sub_vec = _sub_vec;
184 }
185 
186 // Methods overridden from VectorMutable
187 
189 {
190  RTOpPack::SubVectorView<RTOp_value_type> _sub_vec = *sub_vec;
191  thyra_vec_->acquireDetachedView(convert(rng),&_sub_vec);
192  *sub_vec = _sub_vec;
193 }
194 
196 {
197  RTOpPack::SubVectorView<RTOp_value_type> _sub_vec = *sub_vec;
198  thyra_vec_->commitDetachedView(&_sub_vec);
199  *sub_vec = _sub_vec;
200  this->has_changed();
201 }
202 
204 {
205  thyra_vec_->setSubVector(sub_vec);
206  this->has_changed();
207 }
208 
210  value_type alpha
211  ,const GenPermMatrixSlice &P
212  ,BLAS_Cpp::Transp P_trans
213  ,const Vector &x
214  ,value_type beta
215  )
216 {
217 #ifdef TEUCHOS_DEBUG
218  TEUCHOS_TEST_FOR_EXCEPT(!(space().is_in_core()));
219 #endif
220  VectorDenseMutableEncap y_de(*this);
221  VectorDenseEncap x_de(x);
222  AbstractLinAlgPack::Vp_StMtV( &y_de(), alpha, P, P_trans, x_de(), beta );
223 }
224 
225 } // end namespace AbstractLinAlgPack
void set_sub_vector(const RTOpPack::SparseSubVector &sub_vec)
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const Thyra::VectorSpaceBase< value_type > > set_uninitialized()
Set to uninitialized and return smart pointer to the internal Thyra::VectorSpaceBase<value_type> obj...
void get_thyra_vector(const VectorSpaceThyra &thyra_vec_spc, const Vector &vec, Teuchos::RCP< const Thyra::VectorBase< value_type > > *thyra_vec)
void commit_thyra_vector(const VectorSpaceThyra &thyra_vec_spc, VectorMutable *vec, Teuchos::RCP< Thyra::VectorBase< value_type > > *thyra_vec)
. One-based subregion index range class.
virtual void has_changed() const
Must be called by any vector subclass that modifies this vector object!
RTOpT< RTOp_value_type > RTOp
Abstract interface for objects that represent a space for mutable coordinate vectors.
T_To & dyn_cast(T_From &from)
Extract a constant DenseLinAlgPack::DVectorSlice view of a Vector object.
void get_sub_vector(const Range1D &rng, RTOpPack::SubVector *sub_vec) const
void commit_sub_vector(RTOpPack::MutableSubVector *sub_vec)
void Vp_StMtV(VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta=1.0)
v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs (BLAS xGEMV)
TypeTo as(const TypeFrom &t)
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
void Vp_StMtV(value_type alpha, const GenPermMatrixSlice &P, BLAS_Cpp::Transp P_trans, const Vector &x, value_type beta)
void initialize(const Teuchos::RCP< const Thyra::VectorSpaceBase< value_type > > &thyra_vec_spc, const inner_prod_ptr_t &inner_prod=Teuchos::null)
Initalize given a smart pointer to a Thyra::VetorSpace object.
void free_sub_vector(RTOpPack::SubVector *sub_vec) const
Teuchos::RCP< Thyra::VectorBase< value_type > > set_uninitialized()
Set to uninitialized and return smart pointer to the internal Thyra::VectorBase object.
Abstract interface for mutable coordinate vectors {abstract}.
SparseSubVectorT< RTOp_value_type > SparseSubVector
void free_thyra_vector(const VectorSpaceThyra &thyra_vec_spc, const Vector &vec, Teuchos::RCP< const Thyra::VectorBase< value_type > > *thyra_vec)
Extract a non-const DenseLinAlgPack::DVectorSlice view of a VectorMutable object. ...
void apply_op_serial(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
Class for a mutable sub-vector.
Teuchos::RCP< const Thyra::VectorBase< value_type > > thyra_vec() const
Return a (converted) smart pointer to the internal smart pointer to the Thyra::VectorBase object...
Transp
TRANS.
void initialize(const Teuchos::RCP< Thyra::VectorBase< value_type > > &thyra_vec)
Initalize given a smart pointer to a Thyra::Vetor object.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
Class for a non-mutable sub-vector.
TEUCHOSCORE_LIB_DLL_EXPORT Teuchos::RCP< WorkspaceStore > get_default_workspace_store()
Concrete matrix type to represent general permutation (mapping) matrices.