Anasazi Version of the Day
AnasaziSpecializedEpetraAdapter.cpp
Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 //
00004 //                 Anasazi: Block Eigensolvers Package
00005 //                 Copyright (2004) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00025 //
00026 // ***********************************************************************
00027 // @HEADER
00028 
00029 #include "AnasaziSpecializedEpetraAdapter.hpp"
00030 #include "Teuchos_ScalarTraits.hpp"
00031 
00036 namespace Anasazi {
00037 
00039   //
00040   //--------Anasazi::EpetraOpMultiVec Implementation-------------------------------
00041   //
00043 
00044   // Construction/Destruction
00045  
00046    EpetraOpMultiVec::EpetraOpMultiVec(const Teuchos::RCP<Epetra_Operator> &Op, const Epetra_BlockMap& Map_in, const int numvecs)
00047      : Epetra_OP( Op )
00048    {
00049      Epetra_MV = Teuchos::rcp( new Epetra_MultiVector(Map_in, numvecs) );
00050      Epetra_MV_Temp = Teuchos::rcp( new Epetra_MultiVector(Map_in, numvecs) );
00051    }    
00052    
00053    EpetraOpMultiVec::EpetraOpMultiVec(const Teuchos::RCP<Epetra_Operator> &Op, 
00054                                       const Epetra_BlockMap& Map_in, double * array, const int numvecs, const int stride)
00055      : Epetra_OP( Op )
00056    {
00057      Epetra_MV = Teuchos::rcp( new Epetra_MultiVector(Copy, Map_in, array, stride, numvecs) ); 
00058      Epetra_MV_Temp = Teuchos::rcp( new Epetra_MultiVector(Map_in, numvecs) ); 
00059    }
00060     
00061    EpetraOpMultiVec::EpetraOpMultiVec(const Teuchos::RCP<Epetra_Operator> &Op, 
00062                                       Epetra_DataAccess CV, const Epetra_MultiVector& P_vec, const std::vector<int>& index)
00063      : Epetra_OP( Op )
00064    {
00065      Epetra_MV = Teuchos::rcp( new Epetra_MultiVector(CV, P_vec, &(const_cast<std::vector<int> &>(index))[0], index.size()) );
00066      Epetra_MV_Temp = Teuchos::rcp( new Epetra_MultiVector( P_vec.Map(), index.size()) );
00067    }
00068 
00069    EpetraOpMultiVec::EpetraOpMultiVec(const EpetraOpMultiVec& P_vec)
00070      : Epetra_OP( P_vec.Epetra_OP )
00071    {
00072      Epetra_MV = Teuchos::rcp( new Epetra_MultiVector( *(P_vec.Epetra_MV) ) );
00073      Epetra_MV_Temp = Teuchos::rcp( new Epetra_MultiVector( *(P_vec.Epetra_MV_Temp) ) );
00074    }
00075 
00076   //
00077   //  member functions inherited from Anasazi::MultiVec
00078   //
00079   //
00080   //  Simulating a virtual copy constructor. If we could rely on the co-variance
00081   //  of virtual functions, we could return a pointer to EpetraOpMultiVec
00082   //  (the derived type) instead of a pointer to the pure virtual base class.
00083   //
00084   
00085   MultiVec<double>* EpetraOpMultiVec::Clone ( const int numvecs ) const
00086   {
00087     EpetraOpMultiVec *ptr_apv = new EpetraOpMultiVec( Epetra_OP, Epetra_MV->Map(), numvecs );
00088     return ptr_apv; // safe upcast.
00089   }
00090   //
00091   //  the following is a virtual copy constructor returning
00092   //  a pointer to the pure virtual class. vector values are
00093   //  copied.
00094   //
00095   
00096   MultiVec<double>* EpetraOpMultiVec::CloneCopy() const
00097   {
00098     EpetraOpMultiVec *ptr_apv = new EpetraOpMultiVec(*this);
00099     return ptr_apv; // safe upcast
00100   }
00101   
00102   
00103   MultiVec<double>* EpetraOpMultiVec::CloneCopy ( const std::vector<int>& index ) const
00104   {
00105     EpetraOpMultiVec * ptr_apv = new EpetraOpMultiVec( Epetra_OP, Copy, *Epetra_MV, index);
00106     return ptr_apv; // safe upcast.
00107   }
00108   
00109   
00110   MultiVec<double>* EpetraOpMultiVec::CloneViewNonConst ( const std::vector<int>& index ) 
00111   {
00112     EpetraOpMultiVec * ptr_apv = new EpetraOpMultiVec( Epetra_OP, View, *Epetra_MV, index);
00113     return ptr_apv; // safe upcast.
00114   }
00115   
00116   const MultiVec<double>* EpetraOpMultiVec::CloneView ( const std::vector<int>& index ) const
00117   {
00118     EpetraOpMultiVec * ptr_apv = new EpetraOpMultiVec( Epetra_OP, View, *Epetra_MV, index);
00119     return ptr_apv; // safe upcast.
00120   }
00121   
00122   
00123   void EpetraOpMultiVec::SetBlock( const MultiVec<double>& A, const std::vector<int>& index ) 
00124   {
00125     // this should be revisited to e
00126     EpetraOpMultiVec temp_vec(Epetra_OP, View, *Epetra_MV, index);
00127 
00128     int numvecs = index.size();
00129     if ( A.GetNumberVecs() != numvecs ) {
00130       std::vector<int> index2( numvecs );
00131       for(int i=0; i<numvecs; i++)
00132         index2[i] = i;
00133       EpetraOpMultiVec *tmp_vec = dynamic_cast<EpetraOpMultiVec *>(&const_cast<MultiVec<double> &>(A)); 
00134       TEUCHOS_TEST_FOR_EXCEPTION( tmp_vec==NULL, std::invalid_argument, "Anasazi::EpetraOpMultiVec::SetBlocks() cast of MultiVec<double> to EpetraOpMultiVec failed.");
00135       EpetraOpMultiVec A_vec(Epetra_OP, View, *(tmp_vec->GetEpetraMultiVector()), index2);
00136       temp_vec.MvAddMv( 1.0, A_vec, 0.0, A_vec );
00137     }
00138     else {
00139       temp_vec.MvAddMv( 1.0, A, 0.0, A );
00140     }
00141   }
00142 
00143   //-------------------------------------------------------------
00144   //
00145   // *this <- alpha * A * B + beta * (*this)
00146   //
00147   //-------------------------------------------------------------
00148   
00149   void EpetraOpMultiVec::MvTimesMatAddMv ( double alpha, const MultiVec<double>& A, 
00150       const Teuchos::SerialDenseMatrix<int,double>& B, double beta ) 
00151   {
00152     Epetra_LocalMap LocalMap(B.numRows(), 0, Epetra_MV->Map().Comm());
00153     Epetra_MultiVector B_Pvec(View, LocalMap, B.values(), B.stride(), B.numCols());
00154     
00155     EpetraOpMultiVec *A_vec = dynamic_cast<EpetraOpMultiVec *>(&const_cast<MultiVec<double> &>(A)); 
00156     TEUCHOS_TEST_FOR_EXCEPTION( A_vec==NULL,  std::invalid_argument, "Anasazi::EpetraOpMultiVec::SetBlocks() cast of MultiVec<double> to EpetraOpMultiVec failed.");
00157     
00158     TEUCHOS_TEST_FOR_EXCEPTION( 
00159         Epetra_MV->Multiply( 'N', 'N', alpha, *(A_vec->GetEpetraMultiVector()), B_Pvec, beta ) != 0,
00160         EpetraSpecializedMultiVecFailure, "Anasazi::EpetraOpMultiVec::MvTimesMatAddMv() call to Epetra_MultiVec::Multiply() returned a nonzero value.");
00161   }
00162 
00163   //-------------------------------------------------------------
00164   //
00165   // *this <- alpha * A + beta * B
00166   //
00167   //-------------------------------------------------------------
00168   
00169   void EpetraOpMultiVec::MvAddMv ( double alpha , const MultiVec<double>& A, 
00170                                  double beta, const MultiVec<double>& B) 
00171   {
00172     EpetraOpMultiVec *A_vec = dynamic_cast<EpetraOpMultiVec *>(&const_cast<MultiVec<double> &>(A)); 
00173     TEUCHOS_TEST_FOR_EXCEPTION( A_vec==NULL,  std::invalid_argument, "Anasazi::EpetraOpMultiVec::MvAddMv() cast of MultiVec<double> to EpetraOpMultiVec failed.");
00174     EpetraOpMultiVec *B_vec = dynamic_cast<EpetraOpMultiVec *>(&const_cast<MultiVec<double> &>(B)); 
00175     TEUCHOS_TEST_FOR_EXCEPTION( B_vec==NULL,  std::invalid_argument, "Anasazi::EpetraOpMultiVec::MvAddMv() cast of MultiVec<double> to EpetraOpMultiVec failed.");
00176     
00177     TEUCHOS_TEST_FOR_EXCEPTION( 
00178         Epetra_MV->Update( alpha, *(A_vec->GetEpetraMultiVector()), beta, *(B_vec->GetEpetraMultiVector()), 0.0 ) != 0,
00179         EpetraSpecializedMultiVecFailure, "Anasazi::EpetraOpMultiVec::MvAddMv() call to Epetra_MultiVec::Update() returned a nonzero value.");
00180   }
00181 
00182   //-------------------------------------------------------------
00183   //
00184   // dense B <- alpha * A^T * OP * (*this)
00185   //
00186   //-------------------------------------------------------------
00187   
00188   void EpetraOpMultiVec::MvTransMv ( double alpha, const MultiVec<double>& A,
00189                                    Teuchos::SerialDenseMatrix<int,double>& B
00190 #ifdef HAVE_ANASAZI_EXPERIMENTAL
00191                                    , ConjType conj
00192 #endif
00193                                   ) const
00194   {    
00195     EpetraOpMultiVec *A_vec = dynamic_cast<EpetraOpMultiVec *>(&const_cast<MultiVec<double> &>(A));
00196     
00197     if (A_vec) {
00198       Epetra_LocalMap LocalMap(B.numRows(), 0, Epetra_MV->Map().Comm());
00199       Epetra_MultiVector B_Pvec(View, LocalMap, B.values(), B.stride(), B.numCols());
00200      
00201       int info = Epetra_OP->Apply( *Epetra_MV, *Epetra_MV_Temp );
00202       TEUCHOS_TEST_FOR_EXCEPTION( info != 0, EpetraSpecializedMultiVecFailure, 
00203         "Anasazi::EpetraOpMultiVec::MvTransMv(): Error returned from Epetra_Operator::Apply()" );
00204 
00205       TEUCHOS_TEST_FOR_EXCEPTION( 
00206         B_Pvec.Multiply( 'T', 'N', alpha, *(A_vec->GetEpetraMultiVector()), *Epetra_MV_Temp, 0.0 ) != 0,
00207         EpetraSpecializedMultiVecFailure, "Anasazi::EpetraOpMultiVec::MvTransMv() call to Epetra_MultiVector::Multiply() returned a nonzero value.");
00208     }
00209   }
00210   
00211   //-------------------------------------------------------------
00212   //
00213   // b[i] = A[i]^T * OP * this[i]
00214   // 
00215   //-------------------------------------------------------------
00216   
00217   void EpetraOpMultiVec::MvDot ( const MultiVec<double>& A, std::vector<double> & b
00218 #ifdef HAVE_ANASAZI_EXPERIMENTAL
00219                                , ConjType conj
00220 #endif
00221                              ) const
00222   {
00223     EpetraOpMultiVec *A_vec = dynamic_cast<EpetraOpMultiVec *>(&const_cast<MultiVec<double> &>(A)); 
00224     TEUCHOS_TEST_FOR_EXCEPTION( A_vec==NULL,  std::invalid_argument, "Anasazi::EpetraOpMultiVec::MvDot() cast of MultiVec<double> to EpetraOpMultiVec failed.");
00225 
00226     int info = Epetra_OP->Apply( *Epetra_MV, *Epetra_MV_Temp );
00227     TEUCHOS_TEST_FOR_EXCEPTION( info != 0, EpetraSpecializedMultiVecFailure, 
00228       "Anasazi::EpetraOpMultiVec::MvDot(): Error returned from Epetra_Operator::Apply()" );
00229 
00230     if (( (int)b.size() >= A_vec->GetNumberVecs() ) ) {
00231       TEUCHOS_TEST_FOR_EXCEPTION( 
00232           Epetra_MV_Temp->Dot( *(A_vec->GetEpetraMultiVector()), &b[0] ) != 0,
00233           EpetraSpecializedMultiVecFailure, "Anasazi::EpetraOpMultiVec::MvDot() call to Epetra_MultiVector::Dot() returned an error.");
00234     }
00235   }
00236 
00237   //-------------------------------------------------------------
00238   //
00239   // normvec[i] = || this[i] ||_OP
00240   // 
00241   //-------------------------------------------------------------
00242 
00243   void EpetraOpMultiVec::MvNorm ( std::vector<double> & normvec ) const
00244   {
00245     int info = Epetra_OP->Apply( *Epetra_MV, *Epetra_MV_Temp );
00246     TEUCHOS_TEST_FOR_EXCEPTION( info != 0, EpetraSpecializedMultiVecFailure,
00247       "Anasazi::EpetraOpMultiVec::MvNorm(): Error returned from Epetra_Operator::Apply()" );
00248 
00249     if (( (int)normvec.size() >= Epetra_MV->NumVectors() ) ) {
00250       TEUCHOS_TEST_FOR_EXCEPTION(
00251           Epetra_MV_Temp->Dot( *Epetra_MV, &normvec[0] ) != 0,
00252           EpetraSpecializedMultiVecFailure, "Anasazi::EpetraOpMultiVec::MvNorm() call to Epetra_MultiVector::Dot() returned an error.");
00253     }
00254 
00255     for (int i=0; i<Epetra_MV->NumVectors(); ++i)
00256       normvec[i] = Teuchos::ScalarTraits<double>::squareroot( normvec[i] );
00257   }
00258 
00259   //-------------------------------------------------------------
00260   //
00261   // this[i] = alpha[i] * this[i]
00262   // 
00263   //-------------------------------------------------------------
00264   void EpetraOpMultiVec::MvScale ( const std::vector<double>& alpha )
00265   {
00266     // Check to make sure the vector is as long as the multivector has columns.
00267     int numvecs = this->GetNumberVecs();
00268     TEUCHOS_TEST_FOR_EXCEPTION( (int)alpha.size() != numvecs, std::invalid_argument, 
00269         "Anasazi::EpetraOpMultiVec::MvScale() alpha argument size was inconsistent with number of vectors in mv.");
00270     
00271     std::vector<int> tmp_index( 1, 0 );
00272     for (int i=0; i<numvecs; i++) {
00273       Epetra_MultiVector temp_vec(View, *Epetra_MV, &tmp_index[0], 1);
00274       TEUCHOS_TEST_FOR_EXCEPTION( 
00275           temp_vec.Scale( alpha[i] ) != 0,
00276           EpetraSpecializedMultiVecFailure, "Anasazi::EpetraOpMultiVec::MvScale() call to Epetra_MultiVector::Scale() returned a nonzero value.");
00277       tmp_index[0]++;
00278     }
00279   }
00280 
00281 } // namespace Anasazi
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends