00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
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
00106 static int GetVecLength( const _MV & mv )
00107 {return AMVT::GetVecLength(mv);}
00108
00109
00110 static int GetNumberVecs( const _MV & mv )
00111 {return AMVT::GetNumberVecs(mv);}
00112
00113
00114
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
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
00127
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
00136
00137 static void MvDot( const _MV & mv, const _MV & A, std::vector<double> &b )
00138 {AMVT::MvDot(mv, A, b);}
00139
00140
00141
00142
00143 static void MvScale ( _MV & mv, const double alpha )
00144 {AMVT::MvScale(mv, alpha);}
00145
00146
00147
00148
00149
00150 static void MvScale ( _MV & mv, const std::vector<double>& alpha )
00151 {AMVT::MvScale(mv, alpha);}
00152
00153
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
00171
00172 static void SetBlock( const _MV & A, const std::vector<int>& index,
00173 _MV & mv )
00174 {AMVT::SetBlock(A, index, mv);}
00175
00176
00177 static void Assign( const _MV& A, _MV& mv ) {
00178 AMVT::Assign (A, mv);
00179 }
00180
00181
00182
00183 static void MvRandom( _MV & mv )
00184 {AMVT::MvRandom(mv);}
00185
00186
00187
00188 static void MvInit( _MV & mv,
00189 double alpha = Teuchos::ScalarTraits<double>::zero() )
00190 { AMVT::MvInit(mv, alpha);}
00191
00192
00193 static void MvPrint( const _MV & mv, std::ostream& os )
00194 { AMVT::MvPrint(mv, os);}
00195
00196 #ifdef HAVE_BELOS_TSQR
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
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