Ifpack_ex_Amesos.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_Map.h"
00039 #include "Epetra_BlockMap.h"
00040 #include "Epetra_LinearProblem.h"
00041 #include "Galeri_Maps.h"
00042 #include "Galeri_CrsMatrices.h"
00043 #include "Teuchos_ParameterList.hpp"
00044 #include "AztecOO.h"
00045 #include "Ifpack.h"
00046 
00047 int main(int argc, char *argv[])
00048 {
00049 
00050   // initialize MPI and Epetra communicator
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   Teuchos::ParameterList GaleriList;
00059 
00060   // The problem is defined on a 2D grid, global size is nx * nx.
00061   int nx = 30; 
00062   GaleriList.set("nx", nx);
00063   GaleriList.set("ny", nx * Comm.NumProc());
00064   GaleriList.set("mx", 1);
00065   GaleriList.set("my", Comm.NumProc());
00066   Epetra_Map* Map = Galeri::CreateMap("Cartesian2D", Comm, GaleriList);
00067   Epetra_RowMatrix* A = Galeri::CreateCrsMatrix("Laplace2D", Map, GaleriList);
00068 
00069   // =============================================================== //
00070   // 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 //
00071   // =============================================================== //
00072 
00073   Teuchos::ParameterList List;
00074 
00075   // allocates an IFPACK factory. No data is associated 
00076   // to this object (only method Create()).
00077   Ifpack Factory;
00078 
00079   // create the preconditioner. For valid PrecType values,
00080   // please check the documentation
00081   string PrecType = "Amesos";
00082   int OverlapLevel = 2; // must be >= 0. If Comm.NumProc() == 1,
00083                         // it is ignored.
00084 
00085   Ifpack_Preconditioner* Prec = Factory.Create(PrecType, A, OverlapLevel);
00086   assert(Prec != 0);
00087 
00088   // specify the Amesos solver to be used. 
00089   // If the selected solver is not available,
00090   // IFPACK will try to use Amesos' KLU (which is usually always
00091   // compiled). Amesos' serial solvers are:
00092   // "Amesos_Klu", "Amesos_Umfpack", "Amesos_Superlu"
00093   List.set("amesos: solver type", "Amesos_Klu");
00094 
00095   // sets the parameters
00096   IFPACK_CHK_ERR(Prec->SetParameters(List));
00097 
00098   // initialize the preconditioner. At this point the matrix must
00099   // have been FillComplete()'d, but actual values are ignored.
00100   // At this call, Amesos will perform the symbolic factorization.
00101   IFPACK_CHK_ERR(Prec->Initialize());
00102 
00103   // Builds the preconditioners, by looking for the values of 
00104   // the matrix. At this call, Amesos will perform the
00105   // numeric factorization.
00106   IFPACK_CHK_ERR(Prec->Compute());
00107 
00108   // =================================================== //
00109   // E N D   O F   I F P A C K   C O N S T R U C T I O N //
00110   // =================================================== //
00111 
00112   // At this point, we need some additional objects
00113   // to define and solve the linear system.
00114 
00115   // defines LHS and RHS
00116   Epetra_Vector LHS(A->OperatorDomainMap());
00117   Epetra_Vector RHS(A->OperatorDomainMap());
00118 
00119   // solution is constant
00120   LHS.PutScalar(1.0);
00121   // now build corresponding RHS
00122   A->Apply(LHS,RHS);
00123 
00124   // now randomize the solution
00125   RHS.Random();
00126 
00127   // need an Epetra_LinearProblem to define AztecOO solver
00128   Epetra_LinearProblem Problem(A,&LHS,&RHS);
00129 
00130   // now we can allocate the AztecOO solver
00131   AztecOO Solver(Problem);
00132 
00133   // specify solver
00134   Solver.SetAztecOption(AZ_solver,AZ_gmres);
00135   Solver.SetAztecOption(AZ_output,32);
00136 
00137   // HERE WE SET THE IFPACK PRECONDITIONER
00138   Solver.SetPrecOperator(Prec);
00139 
00140   // .. and here we solve
00141   // NOTE: with one process, the solver must converge in
00142   // one iteration.
00143   Solver.Iterate(1550,1e-8);
00144 
00145   // free memory
00146   delete Prec;
00147   delete A;
00148   delete Map;
00149   
00150 #ifdef HAVE_MPI
00151   MPI_Finalize() ; 
00152 #endif
00153 
00154     return(EXIT_SUCCESS);
00155 }

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