00001
00002
00003
00004
00005 #ifndef BELOS_PLAYA_ADAPTER_HPP
00006 #define BELOS_PLAYA_ADAPTER_HPP
00007
00008
00009 #include "PlayaDefs.hpp"
00010 #include "PlayaAnasaziAdapter.hpp"
00011 #include "BelosMultiVecTraits.hpp"
00012 #include "BelosOperatorTraits.hpp"
00013
00014 namespace Belos
00015 {
00016
00017
00018 using Playa::Vector;
00019 using Teuchos::RCP;
00020 using Teuchos::Array;
00021 using Anasazi::SimpleMV;
00022
00023
00024 template<>
00025 class MultiVecTraits<double, SimpleMV>
00026 {
00027 public:
00028 typedef SimpleMV _MV;
00029 typedef Teuchos::ScalarTraits<double> SCT;
00030 typedef Anasazi::MultiVecTraits<double, _MV> AMVT;
00031
00032 static double one() {static double rtn = SCT::one(); return rtn;}
00033 static double zero() {static double rtn = SCT::zero(); return rtn;}
00034
00035
00036 static RCP<_MV> Clone( const _MV & mv, const int numvecs )
00037 {return AMVT::Clone(mv, numvecs);}
00038
00039
00040 static RCP< _MV > CloneCopy( const _MV & mv )
00041 {return AMVT::CloneCopy(mv);}
00042
00043
00044 static RCP< _MV > CloneCopy( const _MV & mv, const std::vector<int>& index )
00045 {return AMVT::CloneCopy(mv, index);}
00046
00047
00048 static RCP< _MV > CloneViewNonConst( _MV & mv,
00049 const std::vector<int>& index )
00050 {return AMVT::CloneViewNonConst(mv, index);}
00051
00052
00053 static RCP<_MV>
00054 CloneViewNonConst (_MV & mv,
00055 const Teuchos::Range1D& index)
00056 {
00057 return AMVT::CloneViewNonConst (mv, index);
00058 }
00059
00060
00061 static RCP<const _MV > CloneView( const _MV & mv, const std::vector<int>& index )
00062 {return AMVT::CloneView(mv, index);}
00063
00064
00065 static RCP<const _MV > CloneView( const _MV & mv, const Teuchos::Range1D& index )
00066 {return AMVT::CloneView(mv, index);}
00067
00068
00069 static int GetVecLength( const _MV & mv )
00070 {return AMVT::GetVecLength(mv);}
00071
00072
00073 static int GetNumberVecs( const _MV & mv )
00074 {return AMVT::GetNumberVecs(mv);}
00075
00076
00077
00078 static void MvTimesMatAddMv( const double alpha, const _MV & A,
00079 const Teuchos::SerialDenseMatrix<int,double>& B,
00080 const double beta, _MV & mv )
00081 {AMVT::MvTimesMatAddMv(alpha, A, B, beta, mv);}
00082
00083
00084
00085 static void MvAddMv( const double alpha, const _MV & A,
00086 const double beta, const _MV & B, _MV & mv )
00087 {AMVT::MvAddMv(alpha, A, beta, B, mv);}
00088
00089
00090
00091
00092 static void MvTransMv( const double alpha, const _MV & A, const _MV & mv,
00093 Teuchos::SerialDenseMatrix<int,double>& B )
00094 {AMVT::MvTransMv(alpha, A, mv, B);}
00095
00096
00097
00098
00099
00100 static void MvDot( const _MV & mv, const _MV & A, std::vector<double> &b )
00101 {AMVT::MvDot(mv, A, b);}
00102
00103
00104
00105
00106 static void MvScale ( _MV & mv, const double alpha )
00107 {AMVT::MvScale(mv, alpha);}
00108
00109
00110
00111
00112
00113 static void MvScale ( _MV & mv, const std::vector<double>& alpha )
00114 {AMVT::MvScale(mv, alpha);}
00115
00116
00117 static void MvNorm( const _MV & mv,
00118 std::vector<Teuchos::ScalarTraits<double>::magnitudeType> &normvec,
00119 NormType type = TwoNorm)
00120 {
00121 normvec.resize(mv.size());
00122 for (int i=0; i<mv.size(); i++)
00123 {
00124 if (type==OneNorm)
00125 normvec[i] = mv[i].norm1();
00126 else if (type==TwoNorm)
00127 normvec[i] = mv[i].norm2();
00128 else
00129 normvec[i] = mv[i].normInf();
00130 }
00131 }
00132
00133
00134
00135 static void SetBlock( const _MV & A, const std::vector<int>& index,
00136 _MV & mv )
00137 {AMVT::SetBlock(A, index, mv);}
00138
00139
00140 static void Assign( const _MV& A, _MV& mv ) {
00141 AMVT::Assign (A, mv);
00142 }
00143
00144
00145
00146 static void MvRandom( _MV & mv )
00147 {AMVT::MvRandom(mv);}
00148
00149
00150
00151 static void MvInit( _MV & mv,
00152 double alpha = Teuchos::ScalarTraits<double>::zero() )
00153 { AMVT::MvInit(mv, alpha);}
00154
00155
00156 static void MvPrint( const _MV & mv, std::ostream& os )
00157 { AMVT::MvPrint(mv, os);}
00158
00159 #ifdef HAVE_BELOS_TSQR
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171 typedef Belos::details::StubTsqrAdapter<SimpleMV> tsqr_adaptor_type;
00172 #endif // HAVE_BELOS_TSQR
00173 };
00174
00175
00176
00177
00178
00179 template <>
00180 class OperatorTraits < double, SimpleMV, LinearOperator<double> >
00181 {
00182 public:
00183 typedef Anasazi::SimpleMV _MV;
00184 typedef Anasazi::OperatorTraits<double, SimpleMV, LinearOperator<double> > AOPT;
00185
00186
00187 static void Apply (
00188 const LinearOperator< double >& Op,
00189 const _MV & x,
00190 _MV & y,
00191 ETrans trans=NOTRANS
00192 )
00193 {
00194 if (trans==NOTRANS) AOPT::Apply(Op, x, y);
00195 else AOPT::Apply(Op.transpose(), x, y);
00196 }
00197
00198 };
00199
00200
00201
00202 }
00203
00204 #endif