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);
From: Matthew Knepley <knep...@gmail.com> Date: Monday, 9 September 2024 at 20:57 To: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalin...@stfc.ac.uk> Cc: Stefano Zampini <stefano.zamp...@gmail.com>, petsc-users@mcs.anl.gov <petsc-users@mcs.anl.gov> Subject: Re: [petsc-users] Hypre AMS usage On Mon, Sep 9, 2024 at 3:38 PM Karthikeyan Chockalingam - STFC UKRI <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@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 <knep...@gmail.com<mailto:knep...@gmail.com>> Date: Monday, 9 September 2024 at 19:46 To: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@stfc.ac.uk>> Cc: Stefano Zampini <stefano.zamp...@gmail.com<mailto:stefano.zamp...@gmail.com>>, petsc-users@mcs.anl.gov<mailto:petsc-users@mcs.anl.gov> <petsc-users@mcs.anl.gov<mailto:petsc-users@mcs.anl.gov>> Subject: Re: [petsc-users] Hypre AMS usage On Mon, Sep 9, 2024 at 2:18 PM Karthikeyan Chockalingam - STFC UKRI <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@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 <knep...@gmail.com<mailto:knep...@gmail.com>> Date: Monday, 9 September 2024 at 19:09 To: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@stfc.ac.uk>> Cc: Stefano Zampini <stefano.zamp...@gmail.com<mailto:stefano.zamp...@gmail.com>>, petsc-users@mcs.anl.gov<mailto:petsc-users@mcs.anl.gov> <petsc-users@mcs.anl.gov<mailto:petsc-users@mcs.anl.gov>> Subject: Re: [petsc-users] Hypre AMS usage On Mon, Sep 9, 2024 at 1:21 PM Karthikeyan Chockalingam - STFC UKRI <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBlaDxg3B$ 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.zamp...@gmail.com<mailto:stefano.zamp...@gmail.com>> Date: Monday, 9 September 2024 at 17:02 To: Matthew Knepley <knep...@gmail.com<mailto:knep...@gmail.com>> Cc: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@stfc.ac.uk>>, petsc-users@mcs.anl.gov<mailto:petsc-users@mcs.anl.gov> <petsc-users@mcs.anl.gov<mailto:petsc-users@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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBiuQiYeB$ 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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBpXp4zIJ$ Il giorno lun 9 set 2024 alle ore 18:50 Matthew Knepley <knep...@gmail.com<mailto:knep...@gmail.com>> ha scritto: On Mon, Sep 9, 2024 at 10:17 AM Karthikeyan Chockalingam - STFC UKRI <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBs_nEdNZ$ <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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBt3EToc2$ <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 <knep...@gmail.com<mailto:knep...@gmail.com>> Date: Friday, 6 September 2024 at 17:57 To: Chockalingam, Karthikeyan (STFC,DL,HC) <karthikeyan.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@stfc.ac.uk>> Cc: petsc-users@mcs.anl.gov<mailto:petsc-users@mcs.anl.gov> <petsc-users@mcs.anl.gov<mailto:petsc-users@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@mcs.anl.gov<mailto:petsc-users@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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBp6LJT5B$ <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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBs_nEdNZ$ <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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBlPLfOch$ <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.chockalin...@stfc.ac.uk<mailto:karthikeyan.chockalin...@stfc.ac.uk> [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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBoVYxr6r$ <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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBoVYxr6r$ <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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBoVYxr6r$ <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBtNwizK8$ > -- 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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBoVYxr6r$ <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBtNwizK8$ > -- 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!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBoVYxr6r$ <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!dkm28NStZmKXBlf2nSMgCMkJC8CQ8T2tnmhAwsND0Bz74mxG_LyGyHWDGWF6E85XLcebgv75DeJV7kbMV7Pq8IjFgcgXBtNwizK8$ >