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 #ifndef BELOS_EPETRA_ADAPTER_HPP
00030 #define BELOS_EPETRA_ADAPTER_HPP
00031
00036 #include "Epetra_MultiVector.h"
00037 #include "Epetra_Operator.h"
00038 #include "Epetra_Map.h"
00039 #include "Epetra_LocalMap.h"
00040 #include "Teuchos_SerialDenseMatrix.hpp"
00041
00042 #include "BelosConfigDefs.hpp"
00043 #include "BelosMultiVec.hpp"
00044 #include "BelosOperator.hpp"
00045
00046 namespace Belos {
00047
00048
00049 class EpetraMultiVec : public MultiVec<double>, public Epetra_MultiVector {
00050 public:
00051
00052 EpetraMultiVec(const Epetra_BlockMap& Map, double * array, const int numvecs, const int stride=0);
00053 EpetraMultiVec(const Epetra_BlockMap& Map, const int numvecs);
00054 EpetraMultiVec(Epetra_DataAccess CV, const Epetra_MultiVector& P_vec, const std::vector<int>& index);
00055 EpetraMultiVec& operator=(const EpetraMultiVec& pv) { Epetra_MultiVector::operator=(pv); return *this; }
00056 EpetraMultiVec(const Epetra_MultiVector & P_vec);
00057 ~EpetraMultiVec();
00058
00059
00060
00061
00062
00063
00064
00065
00066 MultiVec<double> * Clone ( const int numvecs ) const;
00067
00068
00069
00070
00071
00072
00073 MultiVec<double> * CloneCopy () const;
00074
00075
00076
00077
00078
00079
00080 MultiVec<double> * CloneCopy ( const std::vector<int>& index ) const;
00081
00082
00083
00084
00085
00086 MultiVec<double> * CloneView ( const std::vector<int>& index );
00087
00088
00089
00090
00091 void SetBlock ( const MultiVec<double>& A, const std::vector<int>& index );
00092
00093 int GetNumberVecs () const { return NumVectors(); }
00094 int GetVecLength () const { return GlobalLength(); }
00095
00096
00097
00098 void MvTimesMatAddMv ( const double alpha, const MultiVec<double>& A,
00099 const Teuchos::SerialDenseMatrix<int,double>& B, const double beta );
00100
00101
00102
00103 void MvAddMv ( const double alpha, const MultiVec<double>& A, const double beta,
00104 const MultiVec<double>& B);
00105
00106
00107
00108 void MvTransMv ( const double alpha, const MultiVec<double>& A, Teuchos::SerialDenseMatrix<int,double>& B ) const;
00109
00110
00111
00112 void MvDot ( const MultiVec<double>& A, std::vector<double>* b ) const;
00113
00114
00115
00116 void MvNorm ( std::vector<double>* normvec, NormType norm_type = TwoNorm ) const;
00117
00118
00119
00120 void MvRandom() { assert( Random() == 0 ); };
00121
00122
00123
00124 void MvInit ( const double alpha ) { assert( PutScalar( alpha ) == 0 ); };
00125
00126
00127
00128 void MvPrint( ostream& os ) const { os << *this << endl; };
00129 private:
00130 };
00131
00132
00134
00135
00136 class EpetraOp : public virtual Operator<double> {
00137 public:
00138 EpetraOp( const Teuchos::RefCountPtr<Epetra_Operator> &Op );
00139 ~EpetraOp() {};
00140 ReturnType Apply ( const MultiVec<double>& x, MultiVec<double>& y, ETrans trans=NOTRANS ) const;
00141 private:
00142 Teuchos::RefCountPtr<Epetra_Operator> Epetra_Op;
00143 };
00144
00146
00147
00148 class EpetraPrecOp : public virtual Operator<double> {
00149 public:
00150 EpetraPrecOp( const Teuchos::RefCountPtr<Epetra_Operator> &Op );
00151 ~EpetraPrecOp() {};
00152 ReturnType Apply ( const MultiVec<double>& x, MultiVec<double>& y, ETrans trans=NOTRANS ) const;
00153 private:
00154 Teuchos::RefCountPtr<Epetra_Operator> Epetra_Op;
00155 };
00156
00157
00159
00160
00161
00163
00164 template<>
00165 class MultiVecTraits<double, Epetra_MultiVector>
00166 {
00167 public:
00169 static Teuchos::RefCountPtr<Epetra_MultiVector> Clone( const Epetra_MultiVector& mv, const int numvecs )
00170 { return Teuchos::rcp( new Epetra_MultiVector(mv.Map(), numvecs) ); }
00172 static Teuchos::RefCountPtr<Epetra_MultiVector> CloneCopy( const Epetra_MultiVector& mv )
00173 { return Teuchos::rcp( new Epetra_MultiVector( mv ) ); }
00175 static Teuchos::RefCountPtr<Epetra_MultiVector> CloneCopy( const Epetra_MultiVector& mv, const std::vector<int>& index )
00176 {
00177 std::vector<int>& tmp_index = const_cast<std::vector<int> &>( index );
00178 return Teuchos::rcp( new Epetra_MultiVector(Copy, mv, &tmp_index[0], index.size()) );
00179 }
00181 static Teuchos::RefCountPtr<Epetra_MultiVector> CloneView( Epetra_MultiVector& mv, const std::vector<int>& index )
00182 {
00183 std::vector<int>& tmp_index = const_cast<std::vector<int> &>( index );
00184 return Teuchos::rcp( new Epetra_MultiVector(View, mv, &tmp_index[0], index.size()) );
00185 }
00187 static Teuchos::RefCountPtr<const Epetra_MultiVector> CloneView( const Epetra_MultiVector& mv, const std::vector<int>& index )
00188 {
00189 std::vector<int>& tmp_index = const_cast<std::vector<int> &>( index );
00190 return Teuchos::rcp( new Epetra_MultiVector(View, mv, &tmp_index[0], index.size()) );
00191 }
00193 static int GetVecLength( const Epetra_MultiVector& mv )
00194 { return mv.GlobalLength(); }
00196 static int GetNumberVecs( const Epetra_MultiVector& mv )
00197 { return mv.NumVectors(); }
00199 static void MvTimesMatAddMv( const double alpha, const Epetra_MultiVector& A,
00200 const Teuchos::SerialDenseMatrix<int,double>& B,
00201 const double beta, Epetra_MultiVector& mv )
00202 {
00203 Epetra_LocalMap LocalMap(B.numRows(), 0, mv.Map().Comm());
00204 Epetra_MultiVector B_Pvec(Copy, LocalMap, B.values(), B.stride(), B.numCols());
00205
00206 assert( mv.Multiply( 'N', 'N', alpha, A, B_Pvec, beta ) == 0 );
00207 }
00209 static void MvAddMv( const double alpha, const Epetra_MultiVector& A, const double beta, const Epetra_MultiVector& B, Epetra_MultiVector& mv )
00210 {
00211 assert( mv.Update( alpha, A, beta, B, 0.0 ) == 0 );
00212 }
00214 static void MvTransMv( const double alpha, const Epetra_MultiVector& A, const Epetra_MultiVector& mv, Teuchos::SerialDenseMatrix<int,double>& B )
00215 {
00216 Epetra_LocalMap LocalMap(B.numRows(), 0, mv.Map().Comm());
00217 Epetra_MultiVector B_Pvec(View, LocalMap, B.values(), B.stride(), B.numCols());
00218
00219 assert( B_Pvec.Multiply( 'T', 'N', alpha, A, mv, 0.0 ) == 0 );
00220 }
00222 static void MvDot( const Epetra_MultiVector& mv, const Epetra_MultiVector& A, std::vector<double>* b )
00223 {
00224 assert( mv.Dot( A, &(*b)[0] ) == 0 );
00225 }
00227 static void MvNorm( const Epetra_MultiVector& mv, std::vector<double>* normvec, NormType type = TwoNorm )
00228 {
00229 if (normvec && ((int)normvec->size() >= mv.NumVectors())) {
00230 switch( type ) {
00231 case ( OneNorm ) :
00232 assert( mv.Norm1(&(*normvec)[0]) == 0 );
00233 break;
00234 case ( TwoNorm ) :
00235 assert( mv.Norm2(&(*normvec)[0]) == 0 );
00236 break;
00237 case ( InfNorm ) :
00238 assert( mv.NormInf(&(*normvec)[0]) == 0 );
00239 break;
00240 default:
00241 break;
00242 }
00243 }
00244 }
00246 static void SetBlock( const Epetra_MultiVector& A, const std::vector<int>& index, Epetra_MultiVector& mv )
00247 {
00248
00249 int numvecs = index.size();
00250 std::vector<int>& tmp_index = const_cast<std::vector<int> &>( index );
00251 Epetra_MultiVector temp_vec(View, mv, &tmp_index[0], numvecs);
00252
00253 if ( A.NumVectors() != numvecs ) {
00254 std::vector<int> index2( numvecs );
00255 for(int i=0; i<numvecs; i++)
00256 index2[i] = i;
00257 Epetra_MultiVector A_vec(View, A, &index2[0], numvecs);
00258 assert( temp_vec.Update( 1.0, A_vec, 0.0, A_vec, 0.0 ) == 0 );
00259 }
00260 else
00261 assert( temp_vec.Update( 1.0, A, 0.0, A, 0.0 ) == 0 );
00262 }
00264 static void MvRandom( Epetra_MultiVector& mv )
00265 { assert( mv.Random() == 0 ); }
00267 static void MvInit( Epetra_MultiVector& mv, double alpha = Teuchos::ScalarTraits<double>::zero() )
00268 { assert( mv.PutScalar(alpha) == 0); }
00270 static void MvPrint( const Epetra_MultiVector& mv, ostream& os )
00271 { os << mv << endl; }
00272
00273 };
00274
00276
00277
00278
00280
00281 template <>
00282 class OperatorTraits < double, Epetra_MultiVector, Epetra_Operator >
00283 {
00284 public:
00285
00287 static ReturnType Apply ( const Epetra_Operator& Op,
00288 const Epetra_MultiVector& x,
00289 Epetra_MultiVector& y,
00290 ETrans trans=NOTRANS )
00291 {
00292 int info = 0;
00293 if ( trans )
00294 const_cast<Epetra_Operator &>(Op).SetUseTranspose( true );
00295 info = Op.Apply( x, y );
00296 if ( trans )
00297 const_cast<Epetra_Operator &>(Op).SetUseTranspose( false );
00298
00299 return ( info == 0 ? Ok : Error );
00300 }
00301
00302 };
00303
00304 }
00305
00306 #endif
00307