PlayaBelosAdapter.hpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                 Playa: Programmable Linear Algebra
00005 //                 Copyright 2012 Sandia Corporation
00006 // 
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Kevin Long (kevin.long@ttu.edu)
00038 // 
00039 
00040 /* @HEADER@ */
00041 
00042 #ifndef BELOS_PLAYA_ADAPTER_HPP
00043 #define BELOS_PLAYA_ADAPTER_HPP
00044 
00045 
00046 #include "PlayaDefs.hpp"
00047 #include "PlayaAnasaziAdapter.hpp"
00048 #include "BelosMultiVecTraits.hpp"
00049 #include "BelosOperatorTraits.hpp"
00050 
00051 namespace Belos
00052 {
00053 
00054 
00055 using Playa::Vector;
00056 using Teuchos::RCP;
00057 using Teuchos::Array;
00058 using Anasazi::SimpleMV;
00059 
00060 /** */
00061 template<> 
00062 class MultiVecTraits<double, SimpleMV>
00063 {
00064 public:
00065   typedef SimpleMV _MV;
00066   typedef Teuchos::ScalarTraits<double> SCT;
00067   typedef Anasazi::MultiVecTraits<double, _MV> AMVT;
00068 
00069   static double one() {static double rtn = SCT::one(); return rtn;}
00070   static double zero() {static double rtn = SCT::zero(); return rtn;}
00071   
00072   /** */
00073   static RCP<_MV> Clone( const  _MV & mv, const int numvecs )
00074     {return AMVT::Clone(mv, numvecs);}
00075 
00076   /** */
00077   static RCP< _MV > CloneCopy( const  _MV & mv )
00078     {return AMVT::CloneCopy(mv);}
00079 
00080   /** */
00081   static RCP< _MV > CloneCopy( const  _MV & mv, const std::vector<int>& index )
00082     {return AMVT::CloneCopy(mv, index);}
00083 
00084   /** */     
00085   static RCP< _MV > CloneViewNonConst(  _MV & mv, 
00086     const std::vector<int>& index )
00087     {return AMVT::CloneViewNonConst(mv, index);}   
00088 
00089   /** */     
00090   static RCP<_MV> 
00091   CloneViewNonConst (_MV & mv, 
00092          const Teuchos::Range1D& index)
00093   {
00094     return AMVT::CloneViewNonConst (mv, index);
00095   }
00096 
00097   /** */
00098   static RCP<const _MV > CloneView( const _MV & mv, const std::vector<int>& index )
00099     {return AMVT::CloneView(mv, index);}
00100 
00101   /** */
00102   static RCP<const _MV > CloneView( const _MV & mv, const Teuchos::Range1D& index )
00103     {return AMVT::CloneView(mv, index);}
00104   
00105   /** Obtain the vector length of \c mv. */
00106   static int GetVecLength( const  _MV & mv )
00107     {return AMVT::GetVecLength(mv);}
00108 
00109   /** Obtain the number of vectors in \c mv */
00110   static int GetNumberVecs( const  _MV & mv )
00111     {return AMVT::GetNumberVecs(mv);}
00112 
00113 
00114   /*! \brief Update \c mv with \f$ \alpha A B + \beta mv \f$ */
00115   static void MvTimesMatAddMv( const double alpha, const  _MV & A, 
00116     const Teuchos::SerialDenseMatrix<int,double>& B, 
00117     const double beta,  _MV & mv )
00118     {AMVT::MvTimesMatAddMv(alpha, A, B, beta, mv);}
00119 
00120   
00121   /*! \brief Replace \c mv with \f$\alpha A + \beta B\f$. */
00122   static void MvAddMv( const double alpha, const  _MV & A, 
00123     const double beta,  const  _MV & B,  _MV & mv )
00124     {AMVT::MvAddMv(alpha, A, beta, B, mv);}
00125 
00126   /*! \brief Compute a dense matrix \c B through 
00127    * the matrix-matrix multiply \f$ \alpha A^Tmv \f$.
00128    */
00129   static void MvTransMv( const double alpha, const  _MV & A, const  _MV & mv, 
00130     Teuchos::SerialDenseMatrix<int,double>& B )
00131     {AMVT::MvTransMv(alpha, A, mv, B);}
00132 
00133   
00134   /**
00135    * Dot product
00136   */
00137   static void MvDot( const  _MV & mv, const  _MV & A, std::vector<double> &b )
00138     {AMVT::MvDot(mv, A, b);}
00139 
00140   
00141   /** Scale each element of the vectors in \c *this with \c alpha.
00142    */
00143   static void MvScale (  _MV & mv, const double alpha )
00144     {AMVT::MvScale(mv, alpha);}
00145       
00146 
00147     
00148   /** Scale each element of the \c i-th vector in \c *this with \c alpha[i].
00149    */
00150   static void MvScale (  _MV & mv, const std::vector<double>& alpha ) 
00151     {AMVT::MvScale(mv, alpha);}
00152     
00153   /** Compute the 2-norm of each individual vector of \c mv. */
00154   static void MvNorm( const  _MV & mv, 
00155     std::vector<Teuchos::ScalarTraits<double>::magnitudeType> &normvec,
00156     NormType type = TwoNorm)
00157     {
00158       normvec.resize(mv.size());
00159       for (int i=0; i<mv.size(); i++)
00160       {
00161         if (type==OneNorm) 
00162           normvec[i] = mv[i].norm1();
00163         else if (type==TwoNorm) 
00164           normvec[i] = mv[i].norm2();
00165         else 
00166           normvec[i] = mv[i].normInf();
00167       }
00168     }
00169 
00170   /*! \brief Copy the vectors in \c A to a set of vectors in \c mv indicated by the indices given in \c index.
00171    */
00172   static void SetBlock( const  _MV & A, const std::vector<int>& index,  
00173     _MV & mv )
00174     {AMVT::SetBlock(A, index, mv);}
00175 
00176   //! Assign (deep copy) A into mv.
00177   static void Assign( const _MV& A, _MV& mv ) {
00178     AMVT::Assign (A, mv);
00179   }
00180 
00181   /*! \brief Replace the vectors in \c mv with random vectors.
00182    */
00183   static void MvRandom(  _MV & mv )
00184     {AMVT::MvRandom(mv);}
00185 
00186   /*! \brief Replace each element of the vectors in \c mv with \c alpha.
00187    */
00188   static void MvInit(  _MV & mv, 
00189     double alpha = Teuchos::ScalarTraits<double>::zero() )
00190     { AMVT::MvInit(mv, alpha);}
00191     
00192   /** Print the \c mv multi-vector to the \c os output stream. */
00193   static void MvPrint( const  _MV & mv, std::ostream& os )
00194     { AMVT::MvPrint(mv, os);}
00195 
00196 #ifdef HAVE_BELOS_TSQR
00197   /// \typedef tsqr_adaptor_type
00198   /// \brief TSQR adapter for the multivector type SimpleMV.
00199   ///
00200   /// For now, we provide a "stub" implementation.  It has the right
00201   /// methods and typedefs, but its constructors and methods all throw
00202   /// std::logic_error.  If you plan to use TSQR in Belos (e.g.,
00203   /// through TsqrOrthoManager) with SimpleMV, you must implement a
00204   /// functional TSQR adapter for SimpleMV.  Please refer to
00205   /// Epetra::TsqrAdapter (for Epetra_MultiVector) or
00206   /// Tpetra::TsqrAdaptor (for Tpetra::MultiVector) for examples of
00207   /// how to implement a TSQR adapter.
00208   typedef Belos::details::StubTsqrAdapter<SimpleMV> tsqr_adaptor_type;
00209 #endif // HAVE_BELOS_TSQR
00210 };
00211 
00212 
00213 /**
00214 
00215 */
00216 template <> 
00217 class OperatorTraits < double, SimpleMV, LinearOperator<double> >
00218 {
00219 public:
00220   typedef Anasazi::SimpleMV _MV;  
00221   typedef Anasazi::OperatorTraits<double, SimpleMV, LinearOperator<double> > AOPT;
00222   /**
00223   */    
00224   static void Apply ( 
00225     const LinearOperator< double >& Op, 
00226     const  _MV & x,  
00227     _MV & y,
00228     ETrans trans=NOTRANS
00229     )
00230     {
00231       if (trans==NOTRANS) AOPT::Apply(Op, x, y);
00232       else AOPT::Apply(Op.transpose(), x, y);
00233     }
00234     
00235 };
00236 
00237 
00238 
00239 }
00240 
00241 #endif

Site Contact