Ifpack_ex_Factory.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.h"
00044 #include "Ifpack_AdditiveSchwarz.h"
00045 
00046 int main(int argc, char *argv[])
00047 {
00048 
00049 #ifdef HAVE_MPI
00050   MPI_Init(&argc,&argv);
00051   Epetra_MpiComm Comm( MPI_COMM_WORLD );
00052 #else
00053   Epetra_SerialComm Comm;
00054 #endif
00055 
00056   Teuchos::ParameterList GaleriList;
00057 
00058   // The problem is defined on a 2D grid, global size is nx * nx.
00059   int nx = 30; 
00060   GaleriList.set("n", nx * nx);
00061   GaleriList.set("nx", nx);
00062   GaleriList.set("ny", nx);
00063   Epetra_Map* Map = Galeri::CreateMap("Linear", Comm, GaleriList);
00064   Epetra_RowMatrix* A = Galeri::CreateCrsMatrix("Laplace2D", Map, GaleriList);
00065 
00066   // =============================================================== //
00067   // 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 //
00068   // =============================================================== //
00069 
00070   Teuchos::ParameterList List;
00071 
00072   // allocates an IFPACK factory. No data is associated 
00073   // to this object (only method Create()).
00074   Ifpack Factory;
00075 
00076   // create the preconditioner. For valid PrecType values,
00077   // please check the documentation
00078   string PrecType = "ILU"; // incomplete LU
00079   int OverlapLevel = 1; // must be >= 0. If Comm.NumProc() == 1,
00080                         // it is ignored.
00081 
00082   Ifpack_Preconditioner* Prec = Factory.Create(PrecType, A, OverlapLevel);
00083   assert(Prec != 0);
00084 
00085   // specify parameters for ILU
00086   List.set("fact: drop tolerance", 1e-9);
00087   List.set("fact: level-of-fill", 1);
00088   // the combine mode is on the following:
00089   // "Add", "Zero", "Insert", "InsertAdd", "Average", "AbsMax"
00090   // Their meaning is as defined in file Epetra_CombineMode.h   
00091   List.set("schwarz: combine mode", "Add");
00092   // sets the parameters
00093   IFPACK_CHK_ERR(Prec->SetParameters(List));
00094 
00095   // initialize the preconditioner. At this point the matrix must
00096   // have been FillComplete()'d, but actual values are ignored.
00097   IFPACK_CHK_ERR(Prec->Initialize());
00098 
00099   // Builds the preconditioners, by looking for the values of 
00100   // the matrix.
00101   IFPACK_CHK_ERR(Prec->Compute());
00102 
00103   // =================================================== //
00104   // E N D   O F   I F P A C K   C O N S T R U C T I O N //
00105   // =================================================== //
00106 
00107   // At this point, we need some additional objects
00108   // to define and solve the linear system.
00109 
00110   // defines LHS and RHS
00111   Epetra_Vector LHS(A->OperatorDomainMap());
00112   Epetra_Vector RHS(A->OperatorDomainMap());
00113 
00114   // solution is constant
00115   LHS.PutScalar(1.0);
00116   // now build corresponding RHS
00117   A->Apply(LHS,RHS);
00118 
00119   // now randomize the solution
00120   RHS.Random();
00121 
00122   // need an Epetra_LinearProblem to define AztecOO solver
00123   Epetra_LinearProblem Problem(A,&LHS,&RHS);
00124 
00125   // now we can allocate the AztecOO solver
00126   AztecOO Solver(Problem);
00127 
00128   // specify solver
00129   Solver.SetAztecOption(AZ_solver,AZ_gmres);
00130   Solver.SetAztecOption(AZ_output,32);
00131 
00132   // HERE WE SET THE IFPACK PRECONDITIONER
00133   Solver.SetPrecOperator(Prec);
00134 
00135   // .. and here we solve
00136   Solver.Iterate(1550,1e-8);
00137 
00138   cout << *Prec;
00139 
00140   // delete memory
00141   delete Prec;
00142   delete Map;
00143   delete A;
00144   
00145 #ifdef HAVE_MPI
00146   MPI_Finalize() ; 
00147 #endif
00148 
00149   return(EXIT_SUCCESS);
00150 }

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