[petsc-users] SNESSetFunctionDomainError with NEWTONTR method

Barry Smith bsmith at petsc.dev
Mon Jul 28 22:20:08 CDT 2025


   Stefano will need to provide the full explanation with trust regions I will just try to provide some background based on looking at the code.

   Something that is not made clear in the manual page is that for SNESSetFunctionDomainError() to work (in line search codes) properly, you must also set some Infinity or NaN into the computed vector.  I have attempted to make this clear in https://urldefense.us/v3/__https://gitlab.com/petsc/petsc/-/merge_requests/8603__;!!G_uCfscf7eWS!c2YPoljsHfCf3J1-tV6cQV7kuT9ezJmi1Kuv-pGMFz6-6ZuXvyEXCAbHOTEZ-AX9a_LlrYH2HTafcS-AU2WCXLw$ 


   In looking at SNESSolve_NEWTONTR() 

   if the initial "guess" is not in the domain; that is you called SNESSetFunctionDomainError() and set some vector entry to NaN or infinity) then it hits the lines

   if (!snes->vec_func_init_set) {
    PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
  } else snes->vec_func_init_set = PETSC_FALSE;

  PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
  SNESCheckFunctionNorm(snes, fnorm);

  and will return from SNESSolve() with SNES_DIVERGED_FUNCTION_DOMAIN (I think this is likely what we want to happen).

   Later in the algorithm, the nonlinear function gets called with

       /* Compute new objective function */
    PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1));
    if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1;

   So if part of the function vector has an Inf or NaN it sets rho to neP->eta1.
   But rho must be > neP->eta1 so the step is rejected and presumably
   the trust region will be shrunk or the algorithm will give up. Note that any value you gave for SNESSetFunctionDomainError() is ignored.

   Barry


    

   






> On Jul 28, 2025, at 12:52 PM, Ali ALI AHMAD <ali.ali_ahmad at utt.fr> wrote:
> 
> Hi,
> 
> I am currently using the NEWTONTR method for compressible flow simulations. In some cases, I obtain non-physical solutions (for example, negative pressure). Can I use SNESSetFunctionDomainError in this case? I tried, but I don’t fully understand how it works with NEWTONTR.
> 
> Thank you in advance
> 
> best regrads,
> Ali ALI AHMAD

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20250728/6f119af0/attachment.html>


More information about the petsc-users mailing list