Ifpack_ex_Reordering.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_Reordering.h"
00038 #include "Ifpack_RCMReordering.h"
00039 #include "Ifpack_ReorderFilter.h"
00040 #include "Ifpack_Utils.h"
00041 
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   // only one process
00054   if (Comm.NumProc() != 1) {
00055 #ifdef HAVE_MPI
00056     MPI_Finalize();
00057 #endif
00058     if (Comm.MyPID() == 0)
00059       cout << "Please run this test with one process only" << endl;
00060     // return success not to break the tests
00061     exit(EXIT_SUCCESS);
00062   }
00063 
00064   // ======================================================== //
00065   // now create the famous "upper arrow" matrix, which        //
00066   // should be reordered as a "lower arrow". Sparsity pattern //
00067   // will be printed on screen.                               //
00068   // ======================================================== //
00069   
00070   int NumPoints = 16;
00071   
00072   Epetra_Map Map(-1,NumPoints,0,Comm);
00073   
00074   vector<int> Indices(NumPoints);
00075   vector<double> Values(NumPoints);
00076 
00077   Epetra_CrsMatrix A(Copy,Map,0);
00078   for (int i = 0 ; i < NumPoints ; ++i) {
00079     
00080     int NumEntries;
00081     if (i == 0) {
00082       NumEntries = NumPoints;
00083       for (int j = 0 ; j < NumPoints ; ++j) {
00084   Indices[j] = j;
00085   Values[j] = 1.0;
00086       }
00087     }
00088     else {
00089       NumEntries = 2;
00090       Indices[0] = 0;
00091       Indices[1] = i;
00092       Values[0] = 1.0;
00093       Values[1] = 1.0;
00094     }
00095 
00096     A.InsertGlobalValues(i, NumEntries, &Values[0], &Indices[0]);
00097   }
00098 
00099   A.FillComplete();
00100 
00101 #ifdef HAVE_IFPACK_TEUCHOS
00102   // print the sparsity to file, postscript format
00104 
00105   // create the reordering...
00106   Ifpack_RCMReordering Reorder;
00107   // and compute is on A
00108   IFPACK_CHK_ERR(Reorder.Compute(A));
00109 
00110   // cout information
00111   cout << Reorder;
00112 
00113   // create a reordered matrix
00114   Ifpack_ReorderFilter ReordA(&A,&Reorder);
00115 
00116   // print the sparsity to file, postscript format
00118 #endif
00119 
00120 #ifdef HAVE_MPI
00121   MPI_Finalize(); 
00122 #endif
00123   return(EXIT_SUCCESS);
00124 
00125 }

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