On Fri, Jul 12, 2024 at 7:25 AM Blauth, Sebastian <
sebastian.bla...@itwm.fraunhofer.de> wrote:

> Dear Barry,
>
>
>
> thanks for the clarification. Oh, that’s unfortunate, that not all TAO
> algorithms use the supplied matrix for the norm (and then probably also not
> for computing inner products in, e.g., the limited memory formulas).
>
>
>
> I fear that I don’t have sufficient time at the moment to make a MR. I
> could, however, provide some “minimal” example where the behavior is shown.
> However, that example would be using petsc4py as I am only familiar with
> that and I would use the fenics FEM package to define the matrices. Would
> this be okay?
>

Yes


> And if that’s the case, should I post the example here or at the petsc
> gitlab?
>

Either place is fine. Gitlab makes it easier to track.

  Thanks,

     Matt


> Best regards,
>
> Sebastian
>
>
>
> --
>
> Dr. Sebastian Blauth
>
> Fraunhofer-Institut für
>
> Techno- und Wirtschaftsmathematik ITWM
>
> Abteilung Transportvorgänge
>
> Fraunhofer-Platz 1, 67663 Kaiserslautern
>
> Telefon: +49 631 31600-4968
>
> sebastian.bla...@itwm.fraunhofer.de
>
> https://urldefense.us/v3/__https://www.itwm.fraunhofer.de__;!!G_uCfscf7eWS!aSNMLW5E3ZNlZwOAfQJMnwBM_4sDBJ8jXhjmfgZFQJJJFSw-7QaqXFgtPRpKudragLDd8ONBV5pJxKWFVczr$
>  
>
>
>
> *From:* Barry Smith <bsm...@petsc.dev>
> *Sent:* Tuesday, July 9, 2024 3:32 PM
> *To:* Blauth, Sebastian <sebastian.bla...@itwm.fraunhofer.de>; Munson,
> Todd <tmun...@mcs.anl.gov>; toby Isaac <toby.is...@gmail.com>
> *Cc:* petsc-users@mcs.anl.gov
> *Subject:* Re: [petsc-users] Questions on TAO and gradient norm / inner
> products
>
>
>
>
>
>    From
>
>
>
> $ git grep TaoGradientNorm
>
> bound/impls/blmvm/blmvm.c:  PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> bound/impls/blmvm/blmvm.c:    PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> bound/impls/bnk/bnk.c:  PetscCall(*TaoGradientNorm*(tao, tao->gradient,
> NORM_2, &bnk->gnorm));
>
> bound/impls/bnk/bnk.c:          PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &bnk->gnorm));
>
> bound/impls/bnk/bnls.c:      PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &bnk->gnorm));
>
> bound/impls/bnk/bntl.c:      PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &bnk->gnorm));
>
> bound/impls/bnk/bntl.c:        PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &bnk->gnorm));
>
> bound/impls/bnk/bntr.c:        PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &bnk->gnorm));
>
> interface/taosolver.c:.seealso: [](ch_tao), `Tao`,
> `TaoGetGradientNorm()`, `*TaoGradientNorm*()`
>
> interface/taosolver.c:.seealso: [](ch_tao), `Tao`,
> `TaoSetGradientNorm()`, `*TaoGradientNorm*()`
>
> interface/taosolver.c:  *TaoGradientNorm* - Compute the norm using the
> `NormType`, the user has selected
>
> interface/taosolver.c:PetscErrorCode *TaoGradientNorm*(Tao tao, Vec
> gradient, NormType type, PetscReal *gnorm)
>
> unconstrained/impls/lmvm/lmvm.c:  PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> unconstrained/impls/lmvm/lmvm.c:      PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> unconstrained/impls/nls/nls.c:  PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> unconstrained/impls/nls/nls.c:          PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> unconstrained/impls/nls/nls.c:    PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> unconstrained/impls/ntr/ntr.c:  PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> unconstrained/impls/ntr/ntr.c:        PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
> unconstrained/impls/ntr/ntr.c:      PetscCall(*TaoGradientNorm*(tao,
> tao->gradient, NORM_2, &gnorm));
>
>
>
> it appears only some of the algorithm implementations use the norm you
> provide. While
>
>
>
>  git grep VecNorm
>
>
>
> indicates many places where it is not used. Likely some of the other
> algorithm implementations could be easily "fixed" to support by
>
> changing the norm computed. But I am not an expert on the algorithms and
> don't know if all algorithms can mathematically support a user provided
> norm.
>
>
>
> You are welcome to take a stab at making the change in an MR, or do you
> have a simple test problem with a mass matrix we can use to fix the
>
> "missing" implementations?
>
>
>
>   Barry
>
>
>
>
>
>
>
>
>
>
>
> On Jul 9, 2024, at 3:47 AM, Blauth, Sebastian <
> sebastian.bla...@itwm.fraunhofer.de> wrote:
>
>
>
> Hello,
>
>
>
> I have some questions regarding TAO and the use the gradient norm.
>
>
>
> First, I want to use a custom inner product for the optimization in TAO
> (for computing the gradient norm and, e.g., in the double loop of a
> quasi-Newton method). I have seen that there is the method
> TAOSetGradientNorm
> https://urldefense.us/v3/__https://petsc.org/release/manualpages/Tao/TaoSetGradientNorm/__;!!G_uCfscf7eWS!aSNMLW5E3ZNlZwOAfQJMnwBM_4sDBJ8jXhjmfgZFQJJJFSw-7QaqXFgtPRpKudragLDd8ONBV5pJxDTBKkKq$
>   which seems
> to do this. According to the petsc4py docs
> https://urldefense.us/v3/__https://petsc.org/release/petsc4py/reference/petsc4py.PETSc.TAO.html*petsc4py.PETSc.TAO.setGradientNorm__;Iw!!G_uCfscf7eWS!aSNMLW5E3ZNlZwOAfQJMnwBM_4sDBJ8jXhjmfgZFQJJJFSw-7QaqXFgtPRpKudragLDd8ONBV5pJxNddlrzE$
>  ,
> this should do what I want. However, the method does not always seem to
> perform correctly: When I use it with “-tao_type lmvm”, it really seems to
> work and gives the correct scaling of the residual in the default TAO
> monitor. However, when I use, e.g., “-tao_type bqnls”, “-tao_type cg”, or
> “-tao_type bncg”, the (initial) residual is the same as it is when I do not
> use the TAOSetGradientNorm. However, there seem to be some slight internal
> changes (at least for the bqnls), as the number of iterations to reach the
> tolerance changes from 15 without TAOSetGradientNorm to 17 with
> TAOSetGradientNorm.
>
>
>
> For the context: Here, I am trying to solve a PDE constrained optimal
> control problem, which I tackle in a reduced fashion (using a reduced cost
> functional which results in an unconstrained optimization using the adjoint
> approach). For this, I would like to use the L2 inner product induced by
> the FEM discretization, so the L2 mass matrix.
>
>
>
> Moreover, I noticed that the performance of “-tao_type lmvm” and
> “-tao_type bqnls” as well as “-tao_type cg” and “-tao_type bncg” are
> drastically different for the same unconstrained problem. I would have
> expected that the algorithms are (more or less) identical for that case. Is
> this to be expected?
>
>
>
> Finally, I would like to use TAO for solving PDE constrained shape
> optimization problems. To do so, I would need to be able to specify the
> inner product used in the solver (see the above part) and this inner
> product would need to change in each iteration. Is it possible to do this
> with TAO? And could anyone give me some hints how to do so in python with
> petsc4py?
>
>
>
> Thanks a lot in advance,
>
> Sebastian
>
>
>
>
>
>
>
> --
>
> Dr. Sebastian Blauth
>
> Fraunhofer-Institut für
>
> Techno- und Wirtschaftsmathematik ITWM
>
> Abteilung Transportvorgänge
>
> Fraunhofer-Platz 1, 67663 Kaiserslautern
>
> Telefon: +49 631 31600-4968
>
> sebastian.bla...@itwm.fraunhofer.de
>
> https://urldefense.us/v3/__https://www.itwm.fraunhofer.de__;!!G_uCfscf7eWS!aSNMLW5E3ZNlZwOAfQJMnwBM_4sDBJ8jXhjmfgZFQJJJFSw-7QaqXFgtPRpKudragLDd8ONBV5pJxKWFVczr$
>  
>
>
>


-- 
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!aSNMLW5E3ZNlZwOAfQJMnwBM_4sDBJ8jXhjmfgZFQJJJFSw-7QaqXFgtPRpKudragLDd8ONBV5pJxCQ5b78V$
  
<https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!aSNMLW5E3ZNlZwOAfQJMnwBM_4sDBJ8jXhjmfgZFQJJJFSw-7QaqXFgtPRpKudragLDd8ONBV5pJxFbGYEze$
 >

Reply via email to