[Trilinos-Users] FW: Two solvers (direct and iterative) to solve the matrices created through Newton method

Heroux, Michael A maherou at sandia.gov
Thu Sep 11 11:23:34 MDT 2008


Forwarding on behalf of Isa.  For some reason the message bounced.

Mike

------ Forwarded Message
From: Isabel Gil <misabel at litec.csic.es>
Date: Thu, 11 Sep 2008 10:36:03 -0600
To: <trilinos-users-bounces at software.sandia.gov>
Conversation: Two solvers (direct and iterative) to solve the matrices created through Newton method
Subject: Two solvers (direct and iterative) to solve the matrices created through Newton method

Hi everybody,

I install the Trilinos package using the following arguments with the
configure command:

***********************************************************************************************
../configure --prefix=/home//LIBRERIA_TRILINOS/trilinos-8.0.8/LINUX_MPI
--enable-mpi --with-mpi-compilers --disable-default-packages
--enable-teuchos --enable-epetra --enable-aztecoo --enable-triutils
--enable-amesos --enable-ifpack --enable-ml --enable-nox
--enable-didasko --enable-sacado --enable-sacado-tests --enable-tests
--enable-shared
***********************************************************************************************

Later, I include the following .h files:

***********************************************************************************************
#include <Epetra_SerialComm.h>
#include <Epetra_Map.h>
#include <Epetra_CrsGraph.h>
#include <Epetra_CrsMatrix.h>
#include <Epetra_Vector.h>
#include <Teuchos_ParameterList.hpp>
#include <AztecOO.h>
#include <AztecOO_Operator.h>
#include <Amesos.h>
#include <Sacado.hpp>
***********************************************************************************************

I define the following variables:

***********************************************************************************************
     Epetra_SerialComm               communicator;
     std::auto_ptr<Epetra_Map>       Map;
     std::auto_ptr<Epetra_CrsMatrix> Matrix;
***********************************************************************************************

In order to solver the system, I have two possibilities: to use a direct
solver or an iterative one:

***********************************************************************************************
   Epetra_Vector x(View, *Map, newton_update.begin());
   Epetra_Vector b(View, *Map, right_hand_side.begin());


   switch (resolviendo)
     {
       case 0: //using a direct solver
       {
        Epetra_LinearProblem prob;
        prob.SetOperator (Matrix.get());

        Amesos_BaseSolver *solver = Amesos().Create ("Amesos_Klu", prob);
        Assert (solver != NULL, ExcInternalError());

        solver->SymbolicFactorization();

        solver->NumericFactorization();

        prob.SetRHS(&b);
        prob.SetLHS(&x);

        solver->Solve();

        delete solver;

        return std::make_pair<unsigned int, double> (0, 0);
       }

       case 1: //using an iterative solver: gmres
       {
        AztecOO solver;
        solver.SetAztecOption(AZ_output,

                               AZ_none
                               );
        solver.SetAztecOption(AZ_solver, AZ_gmres);
        solver.SetRHS(&b);
        solver.SetLHS(&x);

        solver.SetAztecOption(AZ_precond,         AZ_dom_decomp);
        solver.SetAztecOption(AZ_subdomain_solve, AZ_ilut);
        solver.SetAztecOption(AZ_overlap,         0);
        solver.SetAztecOption(AZ_reorder,         0);

        solver.SetAztecParam(AZ_drop,      1e-6);
        solver.SetAztecParam(AZ_ilut_fill, 1.5);
        solver.SetAztecParam(AZ_athresh,   1e-6);
        solver.SetAztecParam(AZ_rthresh,   1.0);

        solver.SetUserMatrix(Matrix.get());

        solver.Iterate(300, 1e-10);

        return std::make_pair<unsigned int, double> (solver.NumIters(),
                                                     solver.TrueResidual());
       }
     }
***********************************************************************************************

Result I have using the direct solver:

***********************************************************************************************
   8.103e-02        0000        0.00e+00
   1.294e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
   1.280e-01        0000        0.00e+00
.
.
.
***********************************************************************************************
It can be seen that the residual doesn't reduce.

Result I have using the iterative solver:

***********************************************************************************************
AztecOO::SetProcConfig ERROR, failed to dynamic_cast Comm to Epetra_MpiComm.
AztecOO::SetProcConfig ERROR, failed to dynamic_cast Comm to Epetra_MpiComm.
***********************************************************************************************

And then, the program breaks.

I don't know where the problem can be.
Any advice will be welcome.

Thanks in advance
Best
Isa




------ End of Forwarded Message
-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://software.sandia.gov/mailman/private/trilinos-users/attachments/20080911/c30f19de/attachment-0001.html 


More information about the Trilinos-Users mailing list