[Trilinos-Users] AztecOO recursion

Mike Heroux maherou at sandia.gov
Mon Nov 7 05:21:37 MST 2005


Bora,

You should not need to use RecursiveIterate().  I cannot see anything else
out of the ordinary.  Can you run Purify or valgrind to check if there is
some errors in your program?  I cannot see any other way that inX would
change.

Mike 

> -----Original Message-----
> From: trilinos-users-bounces at software.sandia.gov 
> [mailto:trilinos-users-bounces at software.sandia.gov] On Behalf 
> Of Bora Ucar
> Sent: Sunday, November 06, 2005 4:54 PM
> To: trilinos-users at software.sandia.gov
> Subject: [Trilinos-Users] AztecOO recursion
> 
> Hi,
> 
> I have to solve a linear system. the system is defined using 
> an Epetra_Operator AOp as
> 
>   Epetra_LinearProblem probDef(AOp, solution, rhs);
>   AztecOO probSolver(probDef);
>   probSolver.SetPrecOperator(P);
>   probSolver.SetAztecOption(AZ_kspace,outerNIter);
>   probSolver.SetAztecOption(AZ_solver, AZ_gmres);
>   probSolver.Iterate(outerNIter, 1.0E-8);
> 
> where P is an Epetra_Operator which will be used as a preconditioner.
> The ApplyInverse of P requires two linear system solutions 
> with matrices P1 and P2.
> (The codes are given below). the problem is that the input 
> inX of ApplyInverse changes during the solve_P2 function; the calls to
> inX.Norm2() at lines L1 and L2 below result in different 
> values. I have tried all combinations of recursiveIterate and 
> iterate in the
> solve_P1 and solve_P2 functions, none worked.
> 
> Any suggestions?
> 
> Bora
> --------------------------------------------------------------
> ---------------------------------------------------------------------
> int myPre::ApplyInverse(const Epetra_MultiVector &inX, 
> Epetra_MultiVector &outY ) const {
> 
>   int s1, s2;
>   double nn;
>   Epetra_MultiVector y1(outY);
> 
>   y1.PutScalar(0.0);
>   s1 = solve_P1(y1, inX);
> 
>   inX.Norm2(&nn);
> /*L1*/  cout <<"Before second solve norm2(inX)= "<<nn<<endl;
> 
>   outY.PutScalar(0.0);
>   s2 = solve_P2(outY, y1);
> 
>   inX.Norm2(&nn);
> /*L2*/  cout <<"After second solve  norm2(inX)="<<nn<<endl;
>   return s1+s2;
> }
> --------------------------------------------------------------
> ---------------------------------------------------------------------
> where solve_P1 and solve_P2 solves linear systems using AzteccOO
> (solve_P2 is the same as above where P1 is replaced by P2):
> 
> 
> int myPre::solve_P1(Epetra_MultiVector &initialGuess, const 
> Epetra_MultiVector &rhs ) const {
> 
> int outerNIter = 50;
>  int rv;
>  Epetra_MultiVector tmpRhs (rhs);
> 
>   Epetra_LinearProblem probDef(P1, &initialGuess, &tmpRhs);
>   AztecOO probSolver(probDef);
> 
>   probSolver.SetAztecOption(AZ_precond, AZ_none);
>   probSolver.SetAztecOption(AZ_output, AZ_all);
>   probSolver.SetAztecOption(AZ_diagnostics, AZ_all);
>   probSolver.SetAztecOption(AZ_kspace, outerNIter);
>   probSolver.SetAztecOption(AZ_solver, AZ_gmres);
>   probSolver.recursiveIterate(outerNIter, 1.0E-5);
> 
>  return rv //rv is set according to the status of probSolver; }
> 
> 
> _______________________________________________
> Trilinos-Users mailing list
> Trilinos-Users at software.sandia.gov
> http://software.sandia.gov/mailman/listinfo/trilinos-users
> 





More information about the Trilinos-Users mailing list