[petsc-users] Hypre AMS usage

Matthew Knepley knepley at gmail.com
Tue Sep 10 13:22:23 CDT 2024


On Tue, Sep 10, 2024 at 1:16 PM Karthikeyan Chockalingam - STFC UKRI <
karthikeyan.chockalingam at stfc.ac.uk> wrote:

> Hi Matt,
>
>
>
> I am not sure, if I understand how to read the source code, let’s take the
> below line
>
>
>
> PetscCall(PetscOptionsReal("-pc_hypre_ams_relax_weight", "Relaxation
> weight for AMS smoother", "None", jac->as_relax_weight,
> &jac->as_relax_weight, &flag3))
>
> (Q1) Am I doing it right by setting as in the line below? How does the
> line below match up to the above? Since the above line has more arguments
>
> petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_relax_weight", "1.0")
>

Yes. All options take a single string argument. The other arguments to the
function help present it to the user (say through a GUI).


> (Q2) How do I go about setting the four parameters below using
> PetscOptionsSetValue?
>

It looks like you have 5 parameters, and you would call

petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_amg_alpha_options",
"10,21,32,43,54")

  Thanks,

     Matt


> PetscCall(PetscOptionsIntArray("-pc_hypre_ams_amg_alpha_options", "AMG
> options for vector Poisson", "None", jac->as_amg_alpha_opts, &n, &flag2));
>
>   if (flag || flag2) {
>
>     PetscCallExternal(HYPRE_AMSSetAlphaAMGOptions, jac->hsolver,
> jac->as_amg_alpha_opts[0], /* AMG coarsen type */
>
>
> jac->as_amg_alpha_opts[1],                                            /*
> AMG agg_levels */
>
>                       jac->as_amg_alpha_opts[2],
>                       /* AMG relax_type */
>
>                       jac->as_amg_alpha_theta,
> jac->as_amg_alpha_opts[3],                   /* AMG interp_type */
>
>
> jac->as_amg_alpha_opts[4]);                                           /*
> AMG Pmax */
>
>   }
>
> Kind regards,
>
> Karthik.
>
>
>
>
>
>
>
> *From: *Matthew Knepley <knepley at gmail.com>
> *Date: *Monday, 9 September 2024 at 22:00
> *To: *Chockalingam, Karthikeyan (STFC,DL,HC) <
> karthikeyan.chockalingam at stfc.ac.uk>
> *Cc: *petsc-users at mcs.anl.gov <petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Hypre AMS usage
>
> On Mon, Sep 9, 2024 at 4:32 PM Karthikeyan Chockalingam - STFC UKRI <
> karthikeyan.chockalingam at stfc.ac.uk> wrote:
>
> Great. Thank you for letting me know.
>
> I got the reference to KSP as well from libmesh and I am not creating a
> new KSP.
>
> This time around, I didn’t have to set the PC type and it seems to work.
>
>
>
> Excellent.
>
>
>
>   Thanks,
>
>
>
>      Matt
>
>
>
>   petscErr = PetscOptionsSetValue(NULL,"-ksp_type", "gmres");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_type", "hypre");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_type", "ams");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_relax_type", "2");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_relax_times", "1");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_relax_weight",
> "1.0");
>
>   petscErr = PetscOptionsSetValue(NULL, "-pc_hypre_ams_omega", "1.0");
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_view", NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_monitor_true_residual",
> NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_converged_reason", NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-options_left", NULL);
>
>   petscErr = KSPSetFromOptions(ksp);
>
>
>
>   // Set discrete gradient
>
>   petscErr = PCHYPRESetDiscreteGradient(pc, par_G);
>
>
>
>   // Set vertex coordinates
>
>   petscErr = PCHYPRESetEdgeConstantVectors(pc, par_xvec, par_yvec,
> par_zvec);
>
>
>
>
>
> *From: *Matthew Knepley <knepley at gmail.com>
> *Date: *Monday, 9 September 2024 at 21:11
> *To: *Chockalingam, Karthikeyan (STFC,DL,HC) <
> karthikeyan.chockalingam at stfc.ac.uk>
> *Cc: *Stefano Zampini <stefano.zampini at gmail.com>, petsc-users at mcs.anl.gov
> <petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Hypre AMS usage
>
> On Mon, Sep 9, 2024 at 4:03 PM Karthikeyan Chockalingam - STFC UKRI <
> karthikeyan.chockalingam at stfc.ac.uk> wrote:
>
> I never explicitly called KSPGetPC(). I am embedding AMS in a libmesh
> (fem) example.  So, I only got the reference to pc from libmesh and set
>  petscErr = PCHYPRESetDiscreteGradient(pc, par_G);
>
>
>
> However, you are creating a _new_ KSP, so your PC will not match it. This
> is not doing what you want.
>
>
>
>   Thanks,
>
>
>
>      Matt
>
>
>
>
>
> *From: *Matthew Knepley <knepley at gmail.com>
> *Date: *Monday, 9 September 2024 at 20:57
> *To: *Chockalingam, Karthikeyan (STFC,DL,HC) <
> karthikeyan.chockalingam at stfc.ac.uk>
> *Cc: *Stefano Zampini <stefano.zampini at gmail.com>, petsc-users at mcs.anl.gov
> <petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Hypre AMS usage
>
> On Mon, Sep 9, 2024 at 3:38 PM Karthikeyan Chockalingam - STFC UKRI <
> karthikeyan.chockalingam at stfc.ac.uk> wrote:
>
> I didn’t know how to check the pc type but adding
>
> petscErr = PCSetType(pc, "hypre");
>
> before the two functions made it to work.
>
>
>
> But how does it work from the command line?
>
>
>
> For summary:
>
>
>
>   KSPCreate(mesh.comm().get(), &ksp);
>
>   petscErr = PetscOptionsSetValue(NULL,"-ksp_type", "gmres");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_type", "hypre");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_type", "ams");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_relax_type", "2");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_relax_times", "1");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_ams_relax_weigh", "1.0");
>
>   petscErr = PetscOptionsSetValue(NULL, "-pc_hypre_ams_omega", "1.0");
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_view", NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_monitor_true_residual",
> NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_converged_reason", NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-options_left", NULL);
>
>   petscErr = KSPSetFromOptions(ksp);
>
>
>
> Where is the KSPGetPC() call?
>
>
>
>   Thanks,
>
>
>
>      Matt
>
>
>
>   // Set pc type
>
>   petscErr = PCSetType(pc, "hypre");
>
>
>
>   // Set discrete gradient
>
>   petscErr = PCHYPRESetDiscreteGradient(pc, par_G);
>
>
>
>   // Set vertex coordinates
>
>   petscErr = PCHYPRESetEdgeConstantVectors(pc, par_xvec, par_yvec,
> par_zvec);
>
>
>
>
>
> *From: *Matthew Knepley <knepley at gmail.com>
> *Date: *Monday, 9 September 2024 at 19:46
> *To: *Chockalingam, Karthikeyan (STFC,DL,HC) <
> karthikeyan.chockalingam at stfc.ac.uk>
> *Cc: *Stefano Zampini <stefano.zampini at gmail.com>, petsc-users at mcs.anl.gov
> <petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Hypre AMS usage
>
> On Mon, Sep 9, 2024 at 2:18 PM Karthikeyan Chockalingam - STFC UKRI <
> karthikeyan.chockalingam at stfc.ac.uk> wrote:
>
> Calling the two functions after KSPSetFromOptions() did not work either.
>
>
>
> Can you check that the PC has the correct type when you call it?
>
>
>
>   Thanks,
>
>
>
>      Matt
>
>
>
> Everything works from the command line.
>
> I haven’t set KSPSetOperators, not sure if that is an issue.
>
>
>
> *From: *Matthew Knepley <knepley at gmail.com>
> *Date: *Monday, 9 September 2024 at 19:09
> *To: *Chockalingam, Karthikeyan (STFC,DL,HC) <
> karthikeyan.chockalingam at stfc.ac.uk>
> *Cc: *Stefano Zampini <stefano.zampini at gmail.com>, petsc-users at mcs.anl.gov
> <petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Hypre AMS usage
>
> On Mon, Sep 9, 2024 at 1:21 PM Karthikeyan Chockalingam - STFC UKRI <
> karthikeyan.chockalingam at stfc.ac.uk> wrote:
>
> Hi Stefano,
>
>
>
> Thank you. That was helpful. I tried the following:
>
>
>
>   petscErr = PetscOptionsSetValue(NULL,"-ksp_type", "gmres");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_type", "hypre");
>
>   petscErr = PetscOptionsSetValue(NULL,"-pc_hypre_type", "ams");
>
>
>
>   // Set discrete gradient
>
>   petscErr = PCHYPRESetDiscreteGradient(pc, par_G);
>
>
>
>   // Set vertex coordinates
>
>   petscErr = PCHYPRESetEdgeConstantVectors(pc, par_xvec, par_yvec,
> par_zvec);
>
>
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_view", NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_monitor_true_residual",
> NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-ksp_converged_reason", NULL);
>
>   petscErr = PetscOptionsSetValue(NULL, "-options_left", NULL);
>
>   petscErr = KSPSetFromOptions(ksp);
>
>
>
>
>
> It errored though I have set PCHYPRESetEdgeConstantVectors.
>
>
>
> My guess is that the "pc" was not yet Hypre, since that type was not set
> until KSPSetFromOptions() was called. So you need
>
> to call those two functions after that call.
>
>
>
>   Thanks,
>
>
>
>     Matt
>
>
>
> But my program works, without PetscOptionsSetValue and by passing
> everything on the command line.
>
>
>
> *[0]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------*
>
> [0]PETSC ERROR: HYPRE AMS preconditioner needs either the coordinate
> vectors via PCSetCoordinates() or the edge constant vectors via
> PCHYPRESetEdgeConstantVectors() or the interpolation matrix via
> PCHYPRESetInterpolations()
>
> [0]PETS C ERROR: WARNING! There are unused option(s) set! Could be the
> program crashed before usage or a spelling mistake, etc!
>
> [0]PETSC ERROR:   Option left: name:-ksp_converged_reason (no value)
> source: code
>
> [0]PETSC ERROR:   Option left: name:-options_left (no value) source: code
>
> [0]PETSC ERROR: See https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNFOiY6cO$  for trouble shooting.
>
> [0]PETSC ERROR: Petsc Release Version 3.20.3, unknown
>
> [0]PETSC ERROR: ./example-dbg on a arch-moose named HC20210312 by
> karthikeyan.chockalingam Mon Sep  9 18:10:13 2024
>
> [0]PETSC ERROR: Configure options --with-64-bit-indices
> --with-cxx-dialect=C++17 --with-debugging=no --with-fortran-bindings=0
> --with-mpi=1 --with-openmp=1 --with-shared-libraries=1 --with-sowing=0
> --download-fblaslapack=1 --download-hypre=1 --download-metis=1
> --download-mumps=1 --download-ptscotch=1 --download-parmetis=1
> --download-scalapack=1 --download-slepc=1 --download-strumpack=1
> --download-superlu_dist=1 --download-bison=1 --download-hdf5=1
> --with-hdf5-fortran-bindings=0
> --download-hdf5-configure-arguments="--with-zlib" --with-make-np=4
>
> [0]PETSC ERROR: #1 PCSetUp_HYPRE() at
> /Users/karthikeyan.chockalingam/moose/petsc/src/ksp/pc/impls/hypre/hypre.c:295
>
> [0]PETSC ERROR: #2 PCSetUp() at
> /Users/karthikeyan.chockalingam/moose/petsc/src/ksp/pc/interface/precon.c:1080
>
> [0]PETSC ERROR: #3 KSPSetUp() at
> /Users/karthikeyan.chockalingam/moose/petsc/src/ksp/ksp/interface/itfunc.c:415
>
> [0]PETSC ERROR: #4 KSPSolve_Private() at
> /Users/karthikeyan.chockalingam/moose/petsc/src/ksp/ksp/interface/itfunc.c:833
>
> [0]PETSC ERROR: #5 KSPSolve() at
> /Users/karthikeyan.chockalingam/moose/petsc/src/ksp/ksp/interface/itfunc.c:1080
>
> libMesh terminating:
>
> HYPRE AMS preconditioner needs either the coordinate vectors via
> PCSetCoordinates() or the edge constant vectors via
> PCHYPRESetEdgeConstantVectors() or the interpolation matrix via
> PCHYPRESetInterpolations()
>
>
>
>
>
> *From: *Stefano Zampini <stefano.zampini at gmail.com>
> *Date: *Monday, 9 September 2024 at 17:02
> *To: *Matthew Knepley <knepley at gmail.com>
> *Cc: *Chockalingam, Karthikeyan (STFC,DL,HC) <
> karthikeyan.chockalingam at stfc.ac.uk>, petsc-users at mcs.anl.gov <
> petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Hypre AMS usage
>
> I would say the best way is to look at the source code
>
>
>
>
> https://urldefense.us/v3/__https://gitlab.com/petsc/petsc/-/blob/main/src/ksp/pc/impls/hypre/hypre.c?ref_type=heads*L2061__;Iw!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNPqq4coE$ 
>
>
> https://urldefense.us/v3/__https://gitlab.com/petsc/petsc/-/blob/main/src/ksp/pc/impls/hypre/hypre.c?ref_type=heads*L1239__;Iw!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNJC3V4RB$ 
>
>
>
> Il giorno lun 9 set 2024 alle ore 18:50 Matthew Knepley <knepley at gmail.com>
> ha scritto:
>
> On Mon, Sep 9, 2024 at 10:17 AM Karthikeyan Chockalingam - STFC UKRI <
> karthikeyan.chockalingam at stfc.ac.uk> wrote:
>
> Hi Matt,
>
>
>
> You mentioned it doesn’t hurt to set the smoothing flags
>
>
>
>
> https://urldefense.us/v3/__https://github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L965__;Iw!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNNfmla18$ 
> <https://urldefense.us/v3/__https:/github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L965__;Iw!!G_uCfscf7eWS!biuCOXkV-qD-e8iTXiwaD9XJ0EWuqr2PidZJdPqbnFzCshQFo-mclGSvGBOPuGSTkjPvMKnIfJzWos7QShcSitpwyrhL55yOMB-v$>
>
> Do you know where I can look for the equivalent PETSc commands? Thank you.
>
>
>
> Yes, this is described here:
> https://urldefense.us/v3/__https://petsc.org/main/manualpages/PC/PCHYPRE/__;!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNJpowz1n$ 
> <https://urldefense.us/v3/__https:/petsc.org/main/manualpages/PC/PCHYPRE/__;!!G_uCfscf7eWS!dM_jKT_vX0NCVNcY9XOTTUaTSW-TF2CgOkB2UDeZ2ztYXapdmLlSBT4XhV7PknN3ewygGWDbKFEP9e0RfePC$>
>
> The best way is to run with -help or look at the source.
>
>
>
>   Thanks,
>
>
>
>      Matt
>
>
>
>
>
> Kind regards,
>
> Karthik.
>
>
>
> *From: *Matthew Knepley <knepley at gmail.com>
> *Date: *Friday, 6 September 2024 at 17:57
> *To: *Chockalingam, Karthikeyan (STFC,DL,HC) <
> karthikeyan.chockalingam at stfc.ac.uk>
> *Cc: *petsc-users at mcs.anl.gov <petsc-users at mcs.anl.gov>
> *Subject: *Re: [petsc-users] Hypre AMS usage
>
> On Fri, Sep 6, 2024 at 11:37 AM Karthikeyan Chockalingam - STFC UKRI via
> petsc-users <petsc-users at mcs.anl.gov> wrote:
>
> Hello,
>
>
>
> I am trying to use the Hypre AMS preconditioner for the first time.
>
>
>
> I am following the example problem from Hypre
>
>
> https://urldefense.us/v3/__https://github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L954__;Iw!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNNbhOKrb$ 
> <https://urldefense.us/v3/__https:/github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L954__;Iw!!G_uCfscf7eWS!biuCOXkV-qD-e8iTXiwaD9XJ0EWuqr2PidZJdPqbnFzCshQFo-mclGSvGBOPuGSTkjPvMKnIfJzWos7QShcSitpwyrhL5ySgA0XM$>
>
>
>
> I have so far successfully set the discrete gradient operator and vertex
> co-ordinates,
>
>
>
> // Set discrete gradient
>
>   petscErr = PCHYPRESetDiscreteGradient(pc, par_G);
>
>
>
>   // Set vertex coordinates
>
>   petscErr = PCHYPRESetEdgeConstantVectors(pc, par_xvec, par_yvec,
> par_zvec);
>
>
>
> Do I need to set the following smoothing options?
>
>
> https://urldefense.us/v3/__https://github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L965__;Iw!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNNfmla18$ 
> <https://urldefense.us/v3/__https:/github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L965__;Iw!!G_uCfscf7eWS!biuCOXkV-qD-e8iTXiwaD9XJ0EWuqr2PidZJdPqbnFzCshQFo-mclGSvGBOPuGSTkjPvMKnIfJzWos7QShcSitpwyrhL55yOMB-v$>
>
>
>
> It cannot hurt. I would set them to begin with.
>
>
>
> Also, do I need to convert from MATMPIAIJ to CSR?
>
>
> https://urldefense.us/v3/__https://github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L984__;Iw!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNEr3SNHH$ 
> <https://urldefense.us/v3/__https:/github.com/hypre-space/hypre/blob/3caa81955eb8d1b4e35d9b450e27cf6d07b50f6e/src/examples/ex15.c*L984__;Iw!!G_uCfscf7eWS!biuCOXkV-qD-e8iTXiwaD9XJ0EWuqr2PidZJdPqbnFzCshQFo-mclGSvGBOPuGSTkjPvMKnIfJzWos7QShcSitpwyrhL53Trv1si$>
>
>
>
> No.
>
>
>
> What are the other PETSc calls to invoke AMS? Is there an example problem
> I can look at?
>
>
>
> I do not know. I don't think we have an example.
>
>
>
>   Thanks,
>
>
>
>     Matt
>
>
>
> Thank you.
>
>
>
> Karthik.
>
>
>
> --
>
> *Karthik Chockalingam, Ph.D.*
>
> Senior Research Software Engineer
>
> High Performance Systems Engineering Group
>
> Hartree Centre | Science and Technology Facilities Council
>
> karthikeyan.chockalingam at stfc.ac.uk
>
>
>
>  [image: signature_3970890138]
>
>
>
>
>
>
> --
>
> 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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$ 
> <https://urldefense.us/v3/__http:/www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!dM_jKT_vX0NCVNcY9XOTTUaTSW-TF2CgOkB2UDeZ2ztYXapdmLlSBT4XhV7PknN3ewygGWDbKFEP9ZyhqOzc$>
>
>
>
>
> --
>
> 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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$ 
> <https://urldefense.us/v3/__http:/www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!dM_jKT_vX0NCVNcY9XOTTUaTSW-TF2CgOkB2UDeZ2ztYXapdmLlSBT4XhV7PknN3ewygGWDbKFEP9ZyhqOzc$>
>
>
>
>
> --
>
> Stefano
>
>
>
>
> --
>
> 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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$ 
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNP14Ygti$ >
>
>
>
>
> --
>
> 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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$ 
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNP14Ygti$ >
>
>
>
>
> --
>
> 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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$ 
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNP14Ygti$ >
>
>
>
>
> --
>
> 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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$ 
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNP14Ygti$ >
>
>
>
>
> --
>
> 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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$ 
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNP14Ygti$ >
>


-- 
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!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNHy1GP_5$  <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fzNdjKwlSE6S1482kG8uLzgyvhIjtNTWlib31rKGN57mOsU51yP6t9MkwWHKKx5J_BOV6XeIS4ibNP14Ygti$ >
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20240910/c82e6be1/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: image001.png
Type: image/png
Size: 18270 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20240910/c82e6be1/attachment-0001.png>


More information about the petsc-users mailing list