[petsc-users] [petsc-maint] norm L2 problemQuestion about changing the norm used in nonlinear solvers (L2 Euclidean vs. L2 Lebesgue)

Matthew Knepley knepley at gmail.com
Wed Aug 13 14:52:09 CDT 2025


Ah, then you can just use

  https://urldefense.us/v3/__https://petsc.org/main/manualpages/SNES/SNESSetConvergenceTest/__;!!G_uCfscf7eWS!b5Rm6Uca7nRuIRvz5loe-2If8OgOq1pxZCIhkk6IOO8PsiK7P7ooXpEK053J03WuANNw-ZftUDLu1NeNMT8G$ 

and compute your weighted norm yourself.

  Thanks,

      Matt

On Wed, Aug 13, 2025 at 3:45 PM Ali ALI AHMAD <ali.ali_ahmad at utt.fr> wrote:

> * As I mentioned earlier, the weighted norm I am referring to is:
>
> ||f||^2_2 = \int_{\Omega} |f(x)|^2 \, dx \approx \sum_{K \in \Omega}
> \sum_q |f(G_q)|^2 \, area(K)
>
> where G_q is the centroid of triangle K, computed using a first-order
> Gaussian quadrature rule.
>
>
> Yes, I am using the inexact Newton method, where the method stops when
>
>
> ||F(x_k) + J_{x_k} \, d_k|| \leq \nu_k \, ||F(x_k)||.
>
>
> Here, I would also like to change the norm to see whether this
> modification affects the quality of the descent direction. I know that
> these norms are theoretically equivalent, but I would like to test them to
> observe any difference between the two. My initial expectation is that the
> method may converge faster when using the L^2 norm.
>
>
> Best regards,
>
> Ali ALI AHMAD
> ------------------------------
> *De: *"Matthew Knepley" <knepley at gmail.com>
> *À: *"Ali ALI AHMAD" <ali.ali_ahmad at utt.fr>
> *Cc: *"Barry Smith" <bsmith at petsc.dev>, "petsc-users" <
> petsc-users at mcs.anl.gov>, "petsc-maint" <petsc-maint at mcs.anl.gov>
> *Envoyé: *Lundi 28 Juillet 2025 18:55:04
> *Objet: *Re: [petsc-maint] norm L2 problemQuestion about changing the
> norm used in nonlinear solvers (L2 Euclidean vs. L2 Lebesgue)
>
> On Mon, Jul 28, 2025 at 11:00 AM Ali ALI AHMAD <ali.ali_ahmad at utt.fr>
> wrote:
>
>> I’m sorry for getting back to you so late. Thank you for your patience
>> and understanding.
>>
>>    -
>>
>>    For example, when using L2 algorithms where a different norm is
>>    applied in the line search, see *Line_search_L2.png* as an example
>>    from this reference: https://urldefense.us/v3/__https://arxiv.org/abs/1607.04254__;!!G_uCfscf7eWS!b5Rm6Uca7nRuIRvz5loe-2If8OgOq1pxZCIhkk6IOO8PsiK7P7ooXpEK053J03WuANNw-ZftUDLu1CcZFu6I$ 
>>    <https://urldefense.us/v3/__https://arxiv.org/abs/1607.04254__;!!G_uCfscf7eWS!ep85s6aXge00Uc1_kixR9yGEU0FBO9W8HTZ6LsVYzy9n9Jel5_r26mPMP138k1Ilhi2CAN6cR4DGe4kKZG997W1E0q3AIw$>.
>>    Here, we can change the norm from the L2 Euclidean to the L2 Lebesgue norm.
>>
>> What does "L2 Lebesgue" mean here? I know what I mean by L_2. It is
>
>   ||f||^2_2 = \int_Omega |f(x)|^2 dx \approx \sum_q |f(x_q)|^2 w_q
>
> Oh, from below it seems you want a weight function wt(x) in the norm. So
> we would have
>
>   \sum_q |f(x_q)|^2 wt(x_q) w_q
>
>>
>>    - For GMRES, we need to replace NORM_2 (L2 Euclidean) in your code
>>    with weighted_NormL2 (L2 Lebesgue) everywhere, including all the
>>    details such as in the Arnoldi algorithm...
>>
>> I do not quite understand here, because GMRES does not use the L:_2 norm.
> It uses the l_2 norm, which is
>
>   ||v||^2_2 = \sum |v_i|^2
>
> The things in the vector are coefficients of basis functions. I guess, if
> you had an interpolatory element, you could interpret this as quadrature
> rule, meaning you would have
>
>   \sum wt(x_i) |v_i|^2
>
> where x_i were the coordinates of the dual basis evaluation functionals.
>
>>
>>    -
>>
>>    For the convergence test as well, in the linear system for finding
>>    the direction *d* (A · d = b),
>>
>> You want to use the inner product that generates your weighted L_2 norm?
>
>>
>>    -
>>
>>    and also when we search for a good step using the line search
>>    formula: xk+1=xk+step×d.
>>
>> You want to minimize your weighted L_2 norm. That should work in the same
> way.
>
>   Thanks,
>
>      Matt
>
>
>> I hope this explanation is clear for you.
>>
>> Best regards,
>> Ali ALI AHMAD
>>
>> ------------------------------
>> *De: *"Barry Smith" <bsmith at petsc.dev>
>> *À: *"Ali ALI AHMAD" <ali.ali_ahmad at utt.fr>
>> *Cc: *"petsc-users" <petsc-users at mcs.anl.gov>, "petsc-maint" <
>> petsc-maint at mcs.anl.gov>
>> *Envoyé: *Vendredi 20 Juin 2025 15:50:17
>> *Objet: *Re: [petsc-maint] norm L2 problemQuestion about changing the
>> norm used in nonlinear solvers (L2 Euclidean vs. L2 Lebesgue)
>>
>>
>>
>> On Jun 20, 2025, at 5:10 AM, Ali ALI AHMAD <ali.ali_ahmad at utt.fr> wrote:
>>
>>
>> * Yes, I am indeed using an inexact Newton method in my code. The descent
>> direction is computed by solving a linear system involving the Jacobian, so
>> the update follows the classical formula "J(un)^{-1}d(un)=-F(un)"  I'm also
>> trying to use a line search strategy based on a weighted L2 norm (in the
>> Lebesgue sense), which a priori should lead to better accuracy and faster
>> convergence in anisotropic settings.
>>
>>
>>    Ok, could you point to sample code (any language) or written
>> algorithms where a different norm is used in the line search?
>>
>>
>> * During the subsequent iterations, I apply the Eisenstat–Walker method
>> to adapt the tolerance, which should also involve modifying the norm used
>> in the algorithm.
>>
>> * The current implementation still uses the standard Euclidean L2 norm in
>> PETSc's linear solver and in GMRES. I believe this should ideally be
>> replaced by a weighted L2 norm consistent with the discretization. However,
>> I haven't yet succeeded in modifying the norm used internally by the linear
>> solver in PETSc, so, I'm not yet sure how much impact this change would
>> have on the overall convergence, but I suspect it could improve robustness,
>> especially for highly anisotropic problems. I would greatly appreciate any
>> guidance on how to implement this properly in PETSc.
>>
>>
>>    Norms are used in multiple ways in GMRES.
>>
>>    1) defining convergence
>>
>>    2) as part of preconditioning
>>
>>   Again can you point to sample code (any language) or written algorithms
>> that describe exactly what you would like to accomplish.
>>
>>    Barry
>>
>>
>> Do not hesitate to contact me again if anything remains unclear or if you
>> need further information.
>>
>> Best regards,
>> Ali ALI AHMAD
>>
>> ------------------------------
>> *De: *"Barry Smith" <bsmith at petsc.dev>
>> *À: *"Ali ALI AHMAD" <ali.ali_ahmad at utt.fr>
>> *Cc: *"petsc-users" <petsc-users at mcs.anl.gov>, "petsc-maint" <
>> petsc-maint at mcs.anl.gov>
>> *Envoyé: *Samedi 14 Juin 2025 01:06:52
>> *Objet: *Re: [petsc-maint] norm L2 problemQuestion about changing the
>> norm used in nonlinear solvers (L2 Euclidean vs. L2 Lebesgue)
>>
>> I appreciate the clarification. I would call 3) preconditioning.
>>    To increase my understanding, you are already using Newton's method?
>> That is, you compute the Jacobian of the function and use  - J^{-1}(u^n)
>> F(u^n) as your update direction?
>>
>>     When you switch the inner product (or precondition) how will the
>> search direction be different?
>>
>>    Thanks
>>
>>    Barry
>>
>> The case you need support for is becoming important to PETSc so we need
>> to understand it well and support it well which is why I am asking these
>> (perhaps to you) trivial questions.
>>
>>
>>
>> On Jun 13, 2025, at 4:55 AM, Ali ALI AHMAD <ali.ali_ahmad at utt.fr> wrote:
>>
>> Thank you for your message.
>>
>> To answer your question: I would like to use the L2 norm in the sense of
>> Lebesgue for *all three purposes*, especially the *third one*.
>>
>> *1- For displaying residuals* during the nonlinear iterations, I would
>> like to observe the convergence behavior using a norm that better reflects
>> the physical properties of the problem.
>>
>> *2- For convergence testing*, I would like the stopping criterion to be
>> based on a weighted L2 norm that accounts for the geometry of the mesh
>> (since I am working with unstructured, anisotropic triangular meshes).
>>
>> *3 - Most importantly*, I would like to modify the *inner product* used
>> in the algorithm so that it aligns with the weighted L2 norm (since I am
>> working with unstructured, anisotropic triangular meshes).
>>
>> Best regards,
>> Ali ALI AHMAD
>> ------------------------------
>> *De: *"Barry Smith" <bsmith at petsc.dev>
>> *À: *"Ali ALI AHMAD" <ali.ali_ahmad at utt.fr>
>> *Cc: *"petsc-users" <petsc-users at mcs.anl.gov>, "petsc-maint" <
>> petsc-maint at mcs.anl.gov>
>> *Envoyé: *Vendredi 13 Juin 2025 03:14:06
>> *Objet: *Re: [petsc-maint] norm L2 problemQuestion about changing the
>> norm used in nonlinear solvers (L2 Euclidean vs. L2 Lebesgue)
>>
>> You haven't answered my question. Where (conceptually) and for what
>> purpose do you want to use the L2 norm.
>>     1) displaying norms to observe the convergence behavior
>>
>>      2) in the convergence testing to determine when to stop
>>
>>      3) changing the "inner product" in the algorithm which amounts to
>> preconditioning.
>>
>>   Barry
>>
>>
>> On Jun 12, 2025, at 9:42 AM, Ali ALI AHMAD <ali.ali_ahmad at utt.fr> wrote:
>>
>> Thank you for your answer.
>>
>> I am currently working with the nonlinear solvers *newtonls* (with bt, l2,
>> etc.) and *newtontr* (using newton, cauchy, and dogleg strategies)
>> combined with the linear solver *gmres* and the *ILU* preconditioner,
>> since my Jacobian matrix is nonsymmetric.
>>
>> I also use the Eisenstat-Walker method for newtonls, as my initial guess
>> is often very far from the exact solution.
>>
>> What I would like to do now is to *replace the standard Euclidean L2
>> norm* with the *L2 norm in the Lebesgue sense in the above numerical
>> algorithm*, because my problem is defined on an *unstructured,
>> anisotropic triangular mesh* where a weighted norm would be more
>> physically appropriate.
>>
>> Would you be able to advise me on how to implement this change properly?
>>
>> I would deeply appreciate any guidance or suggestions you could provide.
>>
>> Thank you in advance for your help.
>>
>> Best regards,
>> Ali ALI AHMAD
>>
>> ------------------------------
>> *De: *"Ali ALI AHMAD" <ali.ali_ahmad at utt.fr>
>> *À: *"Barry Smith" <bsmith at petsc.dev>
>> *Cc: *"petsc-users" <petsc-users at mcs.anl.gov>, "petsc-maint" <
>> petsc-maint at mcs.anl.gov>
>> *Envoyé: *Jeudi 12 Juin 2025 15:28:02
>> *Objet: *Re: [petsc-maint] norm L2 problemQuestion about changing the
>> norm used in nonlinear solvers (L2 Euclidean vs. L2 Lebesgue)
>>
>> Thank you for your answer.
>>
>> I am currently working with the nonlinear solvers *newtonls* (with bt, l2,
>> etc.) and *newtontr* (using newton, cauchy, and dogleg strategies)
>> combined with the linear solver *gmres* and the *ILU* preconditioner,
>> since my Jacobian matrix is nonsymmetric.
>>
>> I also use the Eisenstat-Walker method for newtonls, as my initial guess
>> is often very far from the exact solution.
>>
>> What I would like to do now is to *replace the standard Euclidean L2
>> norm* with the *L2 norm in the Lebesgue sense*, because my problem is
>> defined on an *unstructured, anisotropic triangular mesh* where a
>> weighted norm would be more physically appropriate.
>>
>> Would you be able to advise me on how to implement this change properly?
>>
>> I would deeply appreciate any guidance or suggestions you could provide.
>>
>> Thank you in advance for your help.
>>
>> Best regards,
>> Ali ALI AHMAD
>> ------------------------------
>> *De: *"Barry Smith" <bsmith at petsc.dev>
>> *À: *"Ali ALI AHMAD" <ali.ali_ahmad at utt.fr>
>> *Cc: *"petsc-users" <petsc-users at mcs.anl.gov>, "petsc-maint" <
>> petsc-maint at mcs.anl.gov>
>> *Envoyé: *Jeudi 12 Juin 2025 14:57:40
>> *Objet: *Re: [petsc-maint] norm L2 problemQuestion about changing the
>> norm used in nonlinear solvers (L2 Euclidean vs. L2 Lebesgue)
>>
>>   Do you wish to use a different norm
>>
>>    1) ONLY for displaying (printing out) the residual norms to track
>> progress
>>
>>    2) in the convergence testing
>>
>>    3) to change the numerical algorithm (for example using the L2 inner
>> product instead of the usual linear algebra R^N l2 inner product).
>>
>>    For 1) use SNESMonitorSet() and in your monitor function use
>> SNESGetSolution() to grab the solution and then VecGetArray(). Now you can
>> compute any weighted norm you want on the solution.
>>
>>    For 2) similar but you need to use SNESSetConvergenceTest
>>
>>    For 3) yes, but you need to ask us specifically.
>>
>>   Barry
>>
>>
>> On Jun 11, 2025, at 4:45 AM, Ali ALI AHMAD <ali.ali_ahmad at utt.fr> wrote:
>>
>> Dear PETSc team,
>>
>> I hope this message finds you well.
>>
>> I am currently using PETSc in a C++, where I rely on the nonlinear
>> solvers `SNES` with either `newtonls` or `newtontr` methods. I would like
>> to ask if it is possible to change the default norm used (typically the L2
>> Euclidean norm) to a custom norm, specifically the L2 norm in the sense of
>> Lebesgue (e.g., involving cell-wise weighted integrals over the domain).
>>
>> My main goal is to define a custom residual norm that better reflects the
>> physical quantities of interest in my simulation.
>>
>> Would this be feasible within the PETSc framework? If so, could you point
>> me to the recommended approach (e.g., redefining the norm manually, using
>> specific PETSc hooks or options)?
>>
>> Thank you very much in advance for your help and for the great work on
>> PETSc!
>>
>> Best regards,
>>
>> ------------------------------
>> Ali ALI AHMAD
>> PhD Student
>> University of Technology of Troyes - UTT - France
>> GAMMA3 Project - Office H008 - Phone No: +33 7 67 44 68 18
>> 12 rue Marie Curie - CS 42060 10004 TROYES Cedex
>>
>>
>>
>>
>>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
> https://urldefense.us/v3/__https://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!b5Rm6Uca7nRuIRvz5loe-2If8OgOq1pxZCIhkk6IOO8PsiK7P7ooXpEK053J03WuANNw-ZftUDLu1NXIqWUQ$ 
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!b5Rm6Uca7nRuIRvz5loe-2If8OgOq1pxZCIhkk6IOO8PsiK7P7ooXpEK053J03WuANNw-ZftUDLu1OjIfo34$ >
>
>

-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener

https://urldefense.us/v3/__https://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!b5Rm6Uca7nRuIRvz5loe-2If8OgOq1pxZCIhkk6IOO8PsiK7P7ooXpEK053J03WuANNw-ZftUDLu1NXIqWUQ$  <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!b5Rm6Uca7nRuIRvz5loe-2If8OgOq1pxZCIhkk6IOO8PsiK7P7ooXpEK053J03WuANNw-ZftUDLu1OjIfo34$ >
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20250813/fd382caa/attachment-0001.html>


More information about the petsc-users mailing list