Ifpack_ex_BlockRelaxation.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_MultiVector.h"
00038 #include "Epetra_LinearProblem.h"
00039 #include "Galeri_Maps.h"
00040 #include "Galeri_CrsMatrices.h"
00041 #include "Teuchos_ParameterList.hpp"
00042 #include "AztecOO.h"
00043 #include "Ifpack_AdditiveSchwarz.h"
00044 #include "Ifpack_PointRelaxation.h"
00045 #include "Ifpack_BlockRelaxation.h"
00046 #include "Ifpack_SparseContainer.h"
00047 #include "Ifpack_Amesos.h"
00048 
00049 int main(int argc, char *argv[])
00050 {
00051   // initialize MPI and Epetra communicator
00052 #ifdef HAVE_MPI
00053   MPI_Init(&argc,&argv);
00054   Epetra_MpiComm Comm( MPI_COMM_WORLD );
00055 #else
00056   Epetra_SerialComm Comm;
00057 #endif
00058 
00059   Teuchos::ParameterList GaleriList;
00060 
00061   // The problem is defined on a 2D grid, global size is nx * nx.
00062   int nx = 30; 
00063   GaleriList.set("nx", nx);
00064   GaleriList.set("ny", nx * Comm.NumProc());
00065   GaleriList.set("mx", 1);
00066   GaleriList.set("my", Comm.NumProc());
00067   Epetra_Map* Map = Galeri::CreateMap("Cartesian2D", Comm, GaleriList);
00068   Epetra_RowMatrix* A = Galeri::CreateCrsMatrix("Laplace2D", Map, GaleriList);
00069 
00070   // =============================================================== //
00071   // B E G I N N I N G   O F   I F P A C K   C O N S T R U C T I O N //
00072   // =============================================================== //
00073 
00074   Teuchos::ParameterList List;
00075 
00076   // builds an Ifpack_AdditiveSchwarz. This is templated with
00077   // the local solvers, in this case Ifpack_BlockRelaxation.
00078   // Ifpack_BlockRelaxation requires as a templated a container
00079   // class. A container defines
00080   // how to store the diagonal blocks. Two choices are avaiable:
00081   // Ifpack_DenseContainer (to store them as dense block,
00082   // than use LAPACK' factorization to apply the inverse of
00083   // each block), of Ifpack_SparseContainer (to store
00084   // the diagonal block as Epetra_CrsMatrix's). 
00085   // 
00086   // Here, we use Ifpack_SparseContainer, which in turn is
00087   // templated with the class to use to apply the inverse
00088   // of each block. For example, we can use Ifpack_Amesos.
00089  
00090   // We still have to decide the overlap among the processes,
00091   // and the overlap among the blocks. The two values
00092   // can be different. The overlap among the blocks is
00093   // considered only if block Jacobi is used.
00094   int OverlapProcs = 2;
00095   int OverlapBlocks = 0;
00096 
00097   // define the block below to use dense containers
00098 #if 0
00099   Ifpack_AdditiveSchwarz<Ifpack_BlockRelaxation<Ifpack_DenseContainer> > Prec(A, OverlapProcs);
00100 #else
00101   Ifpack_AdditiveSchwarz<Ifpack_BlockRelaxation<Ifpack_SparseContainer<Ifpack_Amesos> > > Prec(A, OverlapProcs);
00102 #endif
00103 
00104   List.set("relaxation: type", "symmetric Gauss-Seidel");
00105   List.set("partitioner: overlap", OverlapBlocks);
00106 #ifdef HAVE_IFPACK_METIS
00107   // use METIS to create the blocks. This requires --enable-ifpack-metis.
00108   // If METIS is not installed, the user may select "linear". 
00109   List.set("partitioner: type", "metis");
00110 #else
00111   // or a simple greedy algorithm is METIS is not enabled
00112   List.set("partitioner: type", "greedy");
00113 #endif
00114   // defines here the number of local blocks. If 1,
00115   // and only one process is used in the computation, then
00116   // the preconditioner must converge in one iteration. 
00117   List.set("partitioner: local parts", 4);
00118 
00119   // sets the parameters
00120   IFPACK_CHK_ERR(Prec.SetParameters(List));
00121 
00122   // initialize the preconditioner. 
00123   IFPACK_CHK_ERR(Prec.Initialize());
00124 
00125   // Builds the preconditioners.
00126   IFPACK_CHK_ERR(Prec.Compute());
00127 
00128   // =================================================== //
00129   // E N D   O F   I F P A C K   C O N S T R U C T I O N //
00130   // =================================================== //
00131 
00132   // At this point, we need some additional objects
00133   // to define and solve the linear system.
00134 
00135   // defines LHS and RHS
00136   Epetra_Vector LHS(A->OperatorDomainMap());
00137   Epetra_Vector RHS(A->OperatorDomainMap());
00138 
00139   LHS.PutScalar(0.0);
00140   RHS.Random();
00141 
00142   // need an Epetra_LinearProblem to define AztecOO solver
00143   Epetra_LinearProblem Problem(A,&LHS,&RHS);
00144 
00145   // now we can allocate the AztecOO solver
00146   AztecOO Solver(Problem);
00147 
00148   // specify solver
00149   Solver.SetAztecOption(AZ_solver,AZ_cg);
00150   Solver.SetAztecOption(AZ_output,32);
00151 
00152   // HERE WE SET THE IFPACK PRECONDITIONER
00153   Solver.SetPrecOperator(&Prec);
00154 
00155   // .. and here we solve
00156   // NOTE: with one process, the solver must converge in
00157   // one iteration.
00158   Solver.Iterate(1550,1e-5);
00159 
00160   delete A;
00161   delete Map;
00162 
00163 #ifdef HAVE_MPI
00164   MPI_Finalize() ; 
00165 #endif
00166 
00167   return(EXIT_SUCCESS);
00168 }

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