test/Container/cxx_main.cpp

Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 // 
00004 //                IFPACK
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 "Ifpack_ConfigDefs.h"
00030 
00031 #ifdef HAVE_MPI
00032 #include "Epetra_MpiComm.h"
00033 #endif
00034 #include "Epetra_SerialComm.h"
00035 #include "Epetra_CrsMatrix.h"
00036 #include "Epetra_Vector.h"
00037 #include "Epetra_LinearProblem.h"
00038 #include "Galeri_Maps.h"
00039 #include "Galeri_CrsMatrices.h"
00040 #include "Galeri_Utils.h"
00041 #include "Teuchos_ParameterList.hpp"
00042 #include "Ifpack_DenseContainer.h"
00043 #include "Ifpack_SparseContainer.h"
00044 #include "Ifpack_Amesos.h"
00045 #include "Ifpack_LocalFilter.h"
00046 
00047 static bool verbose = false;
00048 
00049 // ======================================================================
00050 bool TestContainer(string Type, Epetra_RowMatrix* A)
00051 {
00052   int NumVectors = 3;
00053   int NumMyRows = A->NumMyRows();
00054 
00055   Epetra_MultiVector LHS_exact(A->RowMatrixRowMap(), NumVectors);
00056   Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
00057   Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
00058   LHS_exact.Random(); LHS.PutScalar(0.0); 
00059   A->Multiply(false, LHS_exact, RHS);
00060 
00061   Epetra_LinearProblem Problem(A, &LHS, &RHS);
00062 
00063   if (verbose) {
00064     cout << "Container type = " << Type << endl;
00065     cout << "NumMyRows = " << NumMyRows << ", NumVectors = " << NumVectors << endl;
00066   }
00067   LHS.PutScalar(0.0);
00068   
00069   Ifpack_Container* Container;
00070 
00071   if (Type == "dense")
00072     Container = new Ifpack_DenseContainer(A->NumMyRows(), NumVectors);
00073   else
00074     Container = new Ifpack_SparseContainer<Ifpack_Amesos>(A->NumMyRows(), NumVectors);
00075 
00076   assert (Container != 0);
00077 
00078   IFPACK_CHK_ERR(Container->Initialize());
00079   // set as ID all the local rows of A
00080   for (int i = 0 ; i < A->NumMyRows() ; ++i)
00081     Container->ID(i) = i;
00082 
00083   // extract submatrix (in this case, the entire matrix)
00084   // and complete setup
00085   IFPACK_CHK_ERR(Container->Compute(*A));
00086 
00087   // set the RHS and LHS
00088   for (int i = 0 ; i < A->NumMyRows() ; ++i)
00089     for (int j = 0 ; j < NumVectors ; ++j) {
00090       Container->RHS(i,j) = RHS[j][i];
00091       Container->LHS(i,j) = LHS[j][i];
00092     }
00093   
00094   // set parameters (empty for dense containers)
00095   Teuchos::ParameterList List;
00096   List.set("amesos: solver type", Type);
00097   IFPACK_CHK_ERR(Container->SetParameters(List));
00098 
00099   // solve the linear system
00100   IFPACK_CHK_ERR(Container->ApplyInverse());
00101 
00102   // get the computed solution, store it in LHS
00103   for (int i = 0 ; i < A->NumMyRows() ; ++i)
00104     for (int j = 0 ; j < NumVectors ; ++j) {
00105        LHS[j][i] = Container->LHS(i,j);
00106     }
00107 
00108   double residual = Galeri::ComputeNorm(&LHS, &LHS_exact);
00109 
00110   if (A->Comm().MyPID() == 0 && verbose) {
00111     cout << "||x_exact - x||_2 = " << residual << endl;
00112     cout << *Container;
00113   }
00114 
00115   bool passed = false;
00116   if (residual < 1e-5)
00117     passed = true;
00118 
00119   delete Container;
00120 
00121   return(passed);
00122 }
00123 
00124 // ======================================================================
00125 int main(int argc, char *argv[])
00126 {
00127 
00128 #ifdef HAVE_MPI
00129   MPI_Init(&argc,&argv);
00130   Epetra_MpiComm Comm(MPI_COMM_WORLD);
00131 #else
00132   Epetra_SerialComm Comm;
00133 #endif
00134   Epetra_SerialComm SerialComm;
00135 
00136   verbose = (Comm.MyPID() == 0);
00137 
00138   const int nx = 30;
00139   Teuchos::ParameterList GaleriList;
00140   GaleriList.set("n", nx * nx);
00141   GaleriList.set("nx", nx);
00142   GaleriList.set("ny", nx);
00143 
00144   Epetra_Map* Map = Galeri::CreateMap("Linear", Comm, GaleriList);
00145   Epetra_RowMatrix* Matrix = Galeri::CreateCrsMatrix("Laplace2D", Map, GaleriList);
00146   
00147   Ifpack_LocalFilter* LocalMatrix = new Ifpack_LocalFilter(Matrix);
00148   int TestPassed = true;
00149 
00150   if (!TestContainer("dense",LocalMatrix))  TestPassed = false;
00151   if (!TestContainer("sparse",LocalMatrix)) TestPassed = false;
00152 
00153   if (TestPassed)
00154   {
00155     if (verbose)
00156       cout << "Test `TestContainer.exe' passed!" << endl;
00157   }
00158   else {
00159     cout << "Test `TestContainer.exe' FAILED!" << endl;
00160     exit(EXIT_FAILURE);
00161   }
00162   
00163   delete LocalMatrix;
00164   delete Matrix;
00165   delete Map;
00166 
00167 #ifdef HAVE_MPI
00168   MPI_Finalize(); 
00169 #endif
00170 
00171   return(EXIT_SUCCESS);
00172 }

Generated on Thu Sep 18 12:37:21 2008 for Ifpack Package Browser (Single Doxygen Collection) by doxygen 1.3.9.1