Ifpack_ex_Filtering.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 #ifdef HAVE_MPI
00031 #include "Epetra_MpiComm.h"
00032 #else
00033 #include "Epetra_SerialComm.h"
00034 #endif
00035 #include "Epetra_Map.h"
00036 #include "Epetra_CrsMatrix.h"
00037 #include "Ifpack_DropFilter.h"
00038 #include "Ifpack_SparsityFilter.h"
00039 #include "Ifpack_SingletonFilter.h"
00040 #include "Ifpack_Utils.h"
00041 
00042 int main(int argc, char *argv[])
00043 {
00044 
00045 #ifdef HAVE_MPI
00046   MPI_Init(&argc,&argv);
00047   Epetra_MpiComm Comm( MPI_COMM_WORLD );
00048 #else
00049   Epetra_SerialComm Comm;
00050 #endif
00051 
00052   if (Comm.NumProc() != 1) {
00053     cerr << "This example must be run with one process only." << endl;
00054     // exit with success not to break the test harness
00055 #ifdef HAVE_MPI
00056     MPI_Finalize();
00057 #endif
00058     exit(EXIT_SUCCESS);
00059   }
00060   
00061   int NumPoints = 5;
00062   Epetra_Map Map(NumPoints,0,Comm);
00063 
00064   Epetra_CrsMatrix Matrix(Copy,Map,0);
00065 
00066   vector<int> Indices(NumPoints);
00067   vector<double> Values(NumPoints);
00068   double Diag = 0.0;
00069 
00070   for (int i = 0 ; i < NumPoints ; ++i) {
00071     // add a diagonal
00072     Matrix.InsertGlobalValues(i,1,&Diag,&i);
00073 
00074     // add off-diagonals
00075     int NumEntries = 0;
00076     for (int j = i + 1 ; j < NumPoints ; ++j) {
00077       Indices[NumEntries] = j;
00078       Values[NumEntries] = 1.0 * (j - i);
00079       ++NumEntries;
00080     }
00081     Matrix.InsertGlobalValues(i,NumEntries,&Values[0],&Indices[0]);
00082   }
00083   Matrix.FillComplete();
00084 
00085   // ================================= //
00086   // print sparsity of original matrix //
00087   // ================================= //
00088  
00089   cout << "Sparsity, non-dropped matrix" << endl;
00090   Ifpack_PrintSparsity_Simple(Matrix);
00091 
00092   // ====================================== //
00093   // create a new matrix, dropping by value //
00094   // ====================================== //
00095   //
00096   // drop all elements below 4.0. Only the upper-right element
00097   // is maintained, plus all the diagonals that are not
00098   // considering in dropping.
00099   Ifpack_DropFilter DropA(&Matrix,4.0);
00100   assert (DropA.MaxNumEntries() == 2);
00101 
00102   cout << "Sparsity, dropping by value" << endl;
00103   Ifpack_PrintSparsity_Simple(DropA);
00104 
00105   // ========================================= //
00106   // create a new matrix, dropping by sparsity //
00107   // ========================================= //
00108   //
00109   // Mantain 2 off-diagonal elements.
00110   Ifpack_SparsityFilter SparsityA(&Matrix,2);
00111 
00112   cout << "Sparsity, dropping by sparsity" << endl;
00113   Ifpack_PrintSparsity_Simple(SparsityA);
00114   assert (SparsityA.MaxNumEntries() == 3);
00115 
00116   // ======================================== //
00117   // create new matrices, dropping singletons //
00118   // ======================================== //
00119   //
00120   // If we apply this filter NumPoints - 1 times, 
00121   // we end up with a one-row matrix
00122   Ifpack_SingletonFilter Filter1(&Matrix);
00123   Ifpack_SingletonFilter Filter2(&Filter1);
00124   Ifpack_SingletonFilter Filter3(&Filter2);
00125   Ifpack_SingletonFilter Filter4(&Filter3);
00126 
00127   cout << "Sparsity, dropping singletons 4 times" << endl;
00128   Ifpack_PrintSparsity_Simple(Filter4);
00129   assert (Filter4.NumMyRows() == 1);
00130 
00131 #ifdef HAVE_MPI
00132   MPI_Finalize();
00133 #endif
00134   return(EXIT_SUCCESS);
00135 }
00136 

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