test/OverlappingRowMatrix/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 #else
00034 #include "Epetra_SerialComm.h"
00035 #endif
00036 #include "Epetra_CrsMatrix.h"
00037 #include "Epetra_Vector.h"
00038 #include "Epetra_LinearProblem.h"
00039 #include "Epetra_Map.h"
00040 #include "Epetra_Import.h"
00041 #include "Epetra_Time.h"
00042 #include "Galeri_Maps.h"
00043 #include "Galeri_CrsMatrices.h"
00044 #include "Teuchos_ParameterList.hpp"
00045 #include "Ifpack_OverlappingRowMatrix.h"
00046 #include "Ifpack_LocalFilter.h"
00047 #include "Ifpack_Utils.h"
00048 
00049 int main(int argc, char *argv[])
00050 {
00051 #ifdef HAVE_MPI
00052   MPI_Init(&argc,&argv);
00053   Epetra_MpiComm Comm( MPI_COMM_WORLD );
00054 #else
00055   Epetra_SerialComm Comm;
00056 #endif
00057 
00058   if (Comm.NumProc() == 1)
00059   {
00060 #ifdef HAVE_MPI
00061     MPI_Finalize();
00062 #endif
00063     exit(EXIT_SUCCESS);
00064   }
00065 
00066   Teuchos::ParameterList GaleriList;
00067   int nx = 100; 
00068   GaleriList.set("n", nx * nx);
00069   GaleriList.set("nx", nx);
00070   GaleriList.set("ny", nx);
00071   Epetra_Map* Map = Galeri::CreateMap("Linear", Comm, GaleriList);
00072   Epetra_CrsMatrix* A = Galeri::CreateCrsMatrix("Laplace2D", Map, GaleriList);
00073 
00074   int OverlapLevel = 5;
00075   Epetra_Time Time(Comm);
00076 
00077   // ======================================== //
00078   // Build the overlapping matrix using class //
00079   // Ifpack_OverlappingRowMatrix.             //
00080   // ======================================== //
00081  
00082   Time.ResetStartTime();
00083   Ifpack_OverlappingRowMatrix B(A,OverlapLevel);
00084   if (Comm.MyPID() == 0)
00085     cout << "Time to create B = " << Time.ElapsedTime() << endl;
00086 
00087   int NumGlobalRowsB = B.NumGlobalRows();
00088   int NumGlobalNonzerosB = B.NumGlobalNonzeros();
00089 
00090   Epetra_Vector X(A->RowMatrixRowMap());
00091   Epetra_Vector Y(A->RowMatrixRowMap());
00092   for (int i = 0 ; i < A->NumMyRows() ; ++i) 
00093     X[i] = 1.0* A->RowMatrixRowMap().GID(i);
00094   Y.PutScalar(0.0);
00095 
00096   Epetra_Vector ExtX_B(B.RowMatrixRowMap());
00097   Epetra_Vector ExtY_B(B.RowMatrixRowMap());
00098   ExtY_B.PutScalar(0.0);
00099 
00100   IFPACK_CHK_ERR(B.ImportMultiVector(X,ExtX_B));
00101   IFPACK_CHK_ERR(B.Multiply(false,ExtX_B,ExtY_B));
00102   IFPACK_CHK_ERR(B.ExportMultiVector(ExtY_B,Y,Add));
00103 
00104   double Norm_B;
00105   Y.Norm2(&Norm_B);
00106   if (Comm.MyPID() == 0)
00107     cout << "Norm of Y using B = " << Norm_B << endl;
00108   
00109   // ================================================== //
00110   //Build the overlapping matrix as an Epetra_CrsMatrix //
00111   // ================================================== //
00112 
00113   Time.ResetStartTime();
00114   Epetra_CrsMatrix& C = 
00115     *(Ifpack_CreateOverlappingCrsMatrix(A,OverlapLevel));
00116   if (Comm.MyPID() == 0)
00117     cout << "Time to create C = " << Time.ElapsedTime() << endl;
00118 
00119   // simple checks on global quantities
00120   int NumGlobalRowsC = C.NumGlobalRows();
00121   int NumGlobalNonzerosC = C.NumGlobalNonzeros();
00122   assert (NumGlobalRowsB == NumGlobalRowsC);
00123   assert (NumGlobalNonzerosB == NumGlobalNonzerosC);
00124 
00125   Epetra_Vector ExtX_C(C.RowMatrixRowMap());
00126   Epetra_Vector ExtY_C(C.RowMatrixRowMap());
00127   ExtY_C.PutScalar(0.0);
00128   Y.PutScalar(0.0);
00129 
00130   IFPACK_CHK_ERR(C.Multiply(false,X,Y));
00131 
00132   double Norm_C;
00133   Y.Norm2(&Norm_C);
00134   if (Comm.MyPID() == 0)
00135     cout << "Norm of Y using C = " << Norm_C << endl;
00136 
00137   if (IFPACK_ABS(Norm_B - Norm_C) > 1e-5)
00138     IFPACK_CHK_ERR(-1);
00139 
00140   // ======================= //
00141   // now localize the matrix //
00142   // ======================= //
00143 
00144   Ifpack_LocalFilter D(&B);
00145 
00146   // free memory
00147   delete A;
00148   delete Map;
00149 
00150 #ifdef HAVE_MPI
00151   MPI_Finalize() ; 
00152 #endif
00153 
00154   if (Comm.MyPID() == 0)
00155     cout << "Test `TestOverlappingRowMatrix.exe' passed!" << endl;
00156 
00157   return(EXIT_SUCCESS);
00158 }

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