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