[Trilinos-Users] AztecOO residual vector

Heroux, Michael A maherou at sandia.gov
Tue Feb 2 19:04:25 EST 2016


Danielle,

If I understand your question correctly, there is no way to have StatusTestResNorm give you more than the final residual.  I see two possibilities:


  1.  You could loop around the Iterate() method repeatedly, doing only enough iterations to complete one Arnoldi block and then call the GetResNormValue() after each of those.
  2.  If you were only interested in getting the history of the residual values after the solve was complete, you could create your own StatusTest derived from the AztecOO_StatusTest base class.  In your class you could capture the history of the residual norm values in an array.  Then the values would be available to you afterwards.  You could use the AztecOO_StatusTestResNorm class as a resource for creating your own status test class.  In fact, I think you would only need a little extra code to record the residual values.

Mike

From: Danielle Catherine Maddix <dcmaddix at stanford.edu<mailto:dcmaddix at stanford.edu>>
Date: Tuesday, February 2, 2016 at 10:47 AM
To: Michael A Heroux <maherou at sandia.gov<mailto:maherou at sandia.gov>>
Cc: Trilinos Users <trilinos-users at trilinos.org<mailto:trilinos-users at trilinos.org>>
Subject: [EXTERNAL] Re: [Trilinos-Users] AztecOO residual vector



Hi Michael,

Yes, the estimated residual value from iterations after restarting would be all I would need.  When I call GetResNorm(), it only returns the norm after the final iteration and convergence.  Is there a way I can make the calls to checkStatus within the Iterate function more often, such as each time after restarting to get this value.  From the documentation, I am a little unclear on how to do this. Any example would be appreciated. Right now, I have:


AztecOO_StatusTestResNorm restartResNorm(*Trilinos::K, *Trilinos::X,

(Epetra_Vector&) Trilinos::F[0], tol);

restartResNorm.DefineResForm(AztecOO_StatusTestResNorm::Implicit,

AztecOO_StatusTestResNorm::TwoNorm);

Solver.SetStatusTest(&restartResNorm);
int status = Solver.Iterate(numMaxIterations, tol);
restartResNorm.GetResNormValue();

Thank you,
Danielle

________________________________
From: Heroux, Michael A <maherou at sandia.gov<mailto:maherou at sandia.gov>>
Sent: Tuesday, February 2, 2016 7:01 AM
To: Duk-Soon Oh; Danielle Catherine Maddix
Cc: trilinos-users at trilinos.org<mailto:trilinos-users at trilinos.org>
Subject: Re: [Trilinos-Users] AztecOO residual vector

Danielle,

One thing to note is that GMRES solvers do not generate an update to the solution or residual vector at each iteration.  Algorithmically there is scalar value that in exact arithmetic represents the norm of the residual at each iteration, which is what we track for monitoring convergence, but it is not efficient to compute the residual until after the solution vector is updated at the end of the Arnoldi process.

Asking for anything but the estimated residual value from GMRES at each iteration is prohibitively expensive.

Mike

From: Trilinos-Users <trilinos-users-bounces at trilinos.org<mailto:trilinos-users-bounces at trilinos.org>> on behalf of Duk-Soon Oh <duksoon at cims.nyu.edu<mailto:duksoon at cims.nyu.edu>>
Date: Tuesday, February 2, 2016 at 3:00 AM
To: Danielle Catherine Maddix <dcmaddix at stanford.edu<mailto:dcmaddix at stanford.edu>>
Cc: Trilinos Users <trilinos-users at trilinos.org<mailto:trilinos-users at trilinos.org>>
Subject: [EXTERNAL] Re: [Trilinos-Users] AztecOO residual vector

Hi,

I just used my custom AztecOO_StatusTest object based on the AztecOO_StatusTestResNorm object. At that time, I had to access the intermediate vectors to debug my custom preconditioner object. I think you can just add few lines in the CheckStatus method.

Best,

Duk-Soon


On Mon, Feb 1, 2016 at 4:31 PM, Danielle Catherine Maddix <dcmaddix at stanford.edu<mailto:dcmaddix at stanford.edu>> wrote:

I got the initialization working AztecOO_StatusTestResNorm restartResNorm(*Trilinos::K, *Trilinos::X,

(Epetra_Vector&) Trilinos::F[0], tol);

restartResNorm.DefineResForm(AztecOO_StatusTestResNorm::Implicit,

AztecOO_StatusTestResNorm::TwoNorm);

Solver.SetStatusTest(&restartResNorm);


I am unclear on how to get the intermediate norms though because it only returns the final one on the last call to check status.

double  GetResNormValue<https://trilinos.org/docs/r8.0/packages/aztecoo/doc/html/classAztecOO__StatusTestResNorm.html#d1c5d2fed099c119384bd7c7dc21649f> () const
        Returns the residual norm value, [$ \|r\| $] , computed in most recent call to CheckStatus.

How do we call checkStatus within the status = Solver.Iterate(numMaxIterations, tol); the inner iteration loops.
Any example of pseudocode would be helpful.

Thanks,
Danielle


________________________________
From: Duk-Soon Oh <duksoon at cims.nyu.edu<mailto:duksoon at cims.nyu.edu>>
Sent: Sunday, January 31, 2016 2:03 PM
To: Danielle Catherine Maddix
Subject: Re: [Trilinos-Users] AztecOO residual vector

Hi,

I had a similar issue previously. Please take a look at the AztecOO_StatusTest object. You may need to define your own AztecOO_StatusTest object to access the intermediate residual vectors. I hope it would be helpful.

Best,

Duk-Soon


On Sun, Jan 31, 2016 at 11:31 AM, Danielle Catherine Maddix <dcmaddix at stanford.edu<mailto:dcmaddix at stanford.edu>> wrote:

Hello,


I was wondering if it would be possible to access the residual vector at each iteration of the AztecOO linear solver.  I know that it returns the final residual.  I want to pass back to our in house finite element code the ratio of the final residual and the residual at the beginning of the last restart iteration in gmres.


Thank you,

Danielle Maddix

_______________________________________________
Trilinos-Users mailing list
Trilinos-Users at trilinos.org<mailto:Trilinos-Users at trilinos.org>
https://trilinos.org/mailman/listinfo/trilinos-users



-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://trilinos.org/pipermail/trilinos-users/attachments/20160203/dd4e569c/attachment.html>


More information about the Trilinos-Users mailing list