<html xmlns:v="urn:schemas-microsoft-com:vml" xmlns:o="urn:schemas-microsoft-com:office:office" xmlns:w="urn:schemas-microsoft-com:office:word" xmlns:m="http://schemas.microsoft.com/office/2004/12/omml" xmlns="http://www.w3.org/TR/REC-html40"><head><meta http-equiv=Content-Type content="text/html; charset=utf-8"><meta name=Generator content="Microsoft Word 15 (filtered medium)"><!--[if !mso]><style>v\:* {behavior:url(#default#VML);}
o\:* {behavior:url(#default#VML);}
w\:* {behavior:url(#default#VML);}
.shape {behavior:url(#default#VML);}
</style><![endif]--><style><!--
/* Font Definitions */
@font-face
{font-family:Helvetica;
panose-1:2 11 6 4 2 2 2 2 2 4;}
@font-face
{font-family:"Cambria Math";
panose-1:2 4 5 3 5 4 6 3 2 4;}
@font-face
{font-family:DengXian;
panose-1:2 1 6 0 3 1 1 1 1 1;}
@font-face
{font-family:Calibri;
panose-1:2 15 5 2 2 2 4 3 2 4;}
@font-face
{font-family:Consolas;
panose-1:2 11 6 9 2 2 4 3 2 4;}
@font-face
{font-family:"\@DengXian";
panose-1:2 1 6 0 3 1 1 1 1 1;}
/* Style Definitions */
p.MsoNormal, li.MsoNormal, div.MsoNormal
{margin:0in;
margin-bottom:.0001pt;
font-size:11.0pt;
font-family:"Calibri",sans-serif;}
a:link, span.MsoHyperlink
{mso-style-priority:99;
color:blue;
text-decoration:underline;}
a:visited, span.MsoHyperlinkFollowed
{mso-style-priority:99;
color:purple;
text-decoration:underline;}
pre
{mso-style-priority:99;
mso-style-link:"HTML Preformatted Char";
margin:0in;
margin-bottom:.0001pt;
font-size:10.0pt;
font-family:"Courier New";}
p.msonormal0, li.msonormal0, div.msonormal0
{mso-style-name:msonormal;
mso-margin-top-alt:auto;
margin-right:0in;
mso-margin-bottom-alt:auto;
margin-left:0in;
font-size:11.0pt;
font-family:"Calibri",sans-serif;}
span.apple-converted-space
{mso-style-name:apple-converted-space;}
span.HTMLPreformattedChar
{mso-style-name:"HTML Preformatted Char";
mso-style-priority:99;
mso-style-link:"HTML Preformatted";
font-family:Consolas;}
span.hljs-keyword
{mso-style-name:hljs-keyword;}
span.line
{mso-style-name:line;}
span.hljs-number
{mso-style-name:hljs-number;}
span.hljs-builtin
{mso-style-name:hljs-builtin;}
span.hljs-string
{mso-style-name:hljs-string;}
span.EmailStyle27
{mso-style-type:personal-reply;
font-family:"Calibri",sans-serif;
color:windowtext;}
span.hljs-comment
{mso-style-name:hljs-comment;}
span.hljs-builtin0
{mso-style-name:hljs-built_in;}
.MsoChpDefault
{mso-style-type:export-only;
font-size:10.0pt;}
@page WordSection1
{size:8.5in 11.0in;
margin:1.0in 1.0in 1.0in 1.0in;}
div.WordSection1
{page:WordSection1;}
--></style><!--[if gte mso 9]><xml>
<o:shapedefaults v:ext="edit" spidmax="1026" />
</xml><![endif]--><!--[if gte mso 9]><xml>
<o:shapelayout v:ext="edit">
<o:idmap v:ext="edit" data="1" />
</o:shapelayout></xml><![endif]--></head><body lang=en-DE link=blue vlink=purple><div class=WordSection1><p class=MsoNormal><span lang=EN-US>Thanks for the clarification. It is more clear to me now about the global to local processes after checking the examples, e.g. ksp/ksp/tutorials/ex14f.F90. <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><span lang=EN-US>And for using Vec locally, I followed your advice of VecGet.. and VecRestore… In fact I used DMDAVecGetArrayReadF90() and some other relevant subroutines. <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><span lang=EN-US>For your comment on DMSNESSetJacobianLocal(). It seems that I need to use both SNESSetJacobian() and then DMSNESSetJacobianLocal() to get things working. When I do only SNESSetJacobian(), it does not work, meaning the following does not work<o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><span lang=EN-US>…… <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US style='font-size:10.0pt;font-family:"Courier New"'> </span><span style='font-size:10.0pt;font-family:"Courier New"'>call DMDASNESsetFunctionLocal(DM_mech,INSERT_VALUES,formResidual,PETSC_NULL_SNES,err_PETSc) </span><span lang=EN-US style='font-size:10.0pt;font-family:"Courier New"'><o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,&<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> 9*product(cells(1:2))*cells3,9*product(cells(1:2))*cells3,&<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> 0,Jac_PETSc,err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call MatShellSetOperation(Jac_PETSc,MATOP_MULT,GK_op,err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call SNESSetJacobian(SNES_mech,Jac_PETSc,Jac_PETSc,PETSC_NULL_FUNCTION,0,err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> !call DMSNESsetJacobianLocal(DM_mech,formJacobian,PETSC_NULL_SNES,err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> !CHKERRQ(err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call SNESsetConvergenceTest(SNES_mech,converged,PETSC_NULL_SNES,PETSC_NULL_FUNCTION,err_PETSc) <o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call SNESSetDM(SNES_mech,DM_mech,err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)<o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>……<o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><span lang=EN-US>It gives me the message <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: No support for this operation for this object type <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: Code not yet written for matrix type shell <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: See https://petsc.org/release/faq/ for trouble shooting. <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: Petsc Release Version 3.16.4, Feb 02, 2022<o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: Configure options PETSC_ARCH=linux-gnu --with-fortran-bindings --with-mpi-f90module-visibility=0 --download-fftw --download-hdf5 --download-hdf5-fortran-bindings --download-fblaslapack --download-ml --download-zlib <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #1 MatFDColoringCreate() at /home/yi/app/petsc-3.16.4/src/mat/matfd/fdmatrix.c:471 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #2 SNESComputeJacobian_DMDA() at /home/yi/app/petsc-3.16.4/src/snes/utils/dmdasnes.c:173[0]PETSC ERROR: #3 SNESComputeJacobian() at /home/yi/app/petsc-3.16.4/src/snes/interface/snes.c:2864 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #4 SNESSolve_NEWTONLS() at /home/yi/app/petsc-3.16.4/src/snes/impls/ls/ls.c:222 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #5 SNESSolve() at /home/yi/app/petsc-3.16.4/src/snes/interface/snes.c:4809 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #6 User provided function() at User file:0 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #7 VecSetErrorIfLocked() at /home/yi/app/petsc-3.16.4/include/petscvec.h:623 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #8 VecGetArray() at /home/yi/app/petsc-3.16.4/src/vec/vec/interface/rvector.c:1769 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #9 User provided function() at User file:0 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #10 MatFDColoringCreate() at /home/yi/app/petsc-3.16.4/src/mat/matfd/fdmatrix.c:471 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #11 SNESComputeJacobian_DMDA() at /home/yi/app/petsc-3.16.4/src/snes/utils/dmdasnes.c:173 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #12 SNESComputeJacobian() at /home/yi/app/petsc-3.16.4/src/snes/interface/snes.c:2864 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #13 SNESSolve_NEWTONLS() at /home/yi/app/petsc-3.16.4/src/snes/impls/ls/ls.c:222 <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>[0]PETSC ERROR: #14 SNESSolve() at /home/yi/app/petsc-3.16.4/src/snes/interface/snes.c:4809<o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><span lang=EN-US>It seems that I have to use a DMSNESSetJacobianLocal() to “activate” the use of my shell matrix, although the formJacobian() in the DMSNESSetJacobianLocal() is doing nothing. <o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><span lang=EN-US>Best wishes,<o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US>Yi<o:p></o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><span lang=EN-US><o:p> </o:p></span></p><p class=MsoNormal><o:p> </o:p></p><div><div style='border:none;border-top:solid #E1E1E1 1.0pt;padding:3.0pt 0in 0in 0in'><p class=MsoNormal><b><span lang=EN-US>From:</span></b><span lang=EN-US> Barry Smith <bsmith@petsc.dev> <br><b>Sent:</b> Tuesday, January 9, 2024 4:49 PM<br><b>To:</b> Yi Hu <y.hu@mpie.de><br><b>Cc:</b> petsc-users <petsc-users@mcs.anl.gov><br><b>Subject:</b> Re: [petsc-users] SNES seems not use my matrix-free operation<o:p></o:p></span></p></div></div><p class=MsoNormal><o:p> </o:p></p><blockquote style='margin-top:5.0pt;margin-bottom:5.0pt'><div><p class=MsoNormal><span lang=EN-US>However, my GK_op (which is the reloaded MATOP_MULT) gives me some problem. It is entered but then crashed with Segmentation Violation error. So I guess my indices may be wrong. I wonder do I need to use the local Vec (of dF), and should my output Vec also be in the correct shape (i.e. after calculation I need to transform back into a Vec)? As you can see here, my dF is a tensor defined on every grid point. </span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></blockquote><div><p class=MsoNormal><o:p> </o:p></p></div><p class=MsoNormal> The input for the matrix-vector product is a global vector, as is the result. (Not like the arguments to DMSNESSetJacobianLocal).<o:p></o:p></p><div><p class=MsoNormal><o:p> </o:p></p></div><div><p class=MsoNormal> This means that your MATOP_MULT function needs to do the DMGlobalToLocal() vector operation first then the "unwrapping" from the vector to the array format at the beginning of the routine. Similarly it needs to "unwrap" the result vector as an array. See src/snes/tutorials/ex14f.F90 and in particular the code block<o:p></o:p></p></div><div><p class=MsoNormal><o:p> </o:p></p></div><div><div><p class=MsoNormal> PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr))<o:p></o:p></p></div><div><p class=MsoNormal> PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr))<o:p></o:p></p></div><div><p class=MsoNormal><o:p> </o:p></p></div><div><p class=MsoNormal>! Get pointers to vector data<o:p></o:p></p></div><div><p class=MsoNormal><o:p> </o:p></p></div><div><p class=MsoNormal> PetscCall(VecGetArrayReadF90(localX,xx,ierr))<o:p></o:p></p></div><div><p class=MsoNormal> PetscCall(VecGetArrayF90(F,ff,ierr))<o:p></o:p></p></div></div><div><p class=MsoNormal><o:p> </o:p></p></div><div><p class=MsoNormal> Barry<o:p></o:p></p></div><div><p class=MsoNormal><o:p> </o:p></p></div><div><p class=MsoNormal>You really shouldn't be using DMSNESSetJacobianLocal() for your code. Basically all the DMSNESSetJacobianLocal() gives you is that it automatically handles the global to local mapping and unwrapping of the vector to an array, but it doesn't work for shell matrices.<o:p></o:p></p></div><div><p class=MsoNormal><o:p> </o:p></p></div><div><div><p class=MsoNormal><br><br><o:p></o:p></p><blockquote style='margin-top:5.0pt;margin-bottom:5.0pt'><div><p class=MsoNormal>On Jan 9, 2024, at 6:30 AM, Yi Hu <<a href="mailto:y.hu@mpie.de">y.hu@mpie.de</a>> wrote:<o:p></o:p></p></div><p class=MsoNormal><o:p> </o:p></p><div><div><p class=MsoNormal><span lang=EN-US>Dear Barry,</span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US>Thanks for your help.<span class=apple-converted-space> </span></span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US>It works when doing first SNESSetJacobian() with my created shell matrix Jac in the main (or module) and then DMSNESSetJacobianLocal() to associate with my DM and an dummy formJacobian callback (which is doing nothing). My SNES can now recognize my shell matrix and do my customized operation.</span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US>However, my GK_op (which is the reloaded MATOP_MULT) gives me some problem. It is entered but then crashed with Segmentation Violation error. So I guess my indices may be wrong. I wonder do I need to use the local Vec (of dF), and should my output Vec also be in the correct shape (i.e. after calculation I need to transform back into a Vec)? As you can see here, my dF is a tensor defined on every grid point.<span class=apple-converted-space> </span></span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US>Best wishes,</span><o:p></o:p></p></div><div><p class=MsoNormal><span lang=EN-US>Yi</span><o:p></o:p></p></div><div><p class=MsoNormal> <o:p></o:p></p></div><div><div style='border:none;border-top:solid windowtext 1.0pt;padding:3.0pt 0in 0in 0in;border-color:currentcolor currentcolor;border-image: none'><div><p class=MsoNormal><b><span lang=EN-US>From:</span></b><span class=apple-converted-space><span lang=EN-US> </span></span><span lang=EN-US>Barry Smith <<a href="mailto:bsmith@petsc.dev"><span style='color:purple'>bsmith@petsc.dev</span></a>><span class=apple-converted-space> </span><br><b>Sent:</b><span class=apple-converted-space> </span>Monday, January 8, 2024 6:41 PM<br><b>To:</b><span class=apple-converted-space> </span>Yi Hu <<a href="mailto:y.hu@mpie.de"><span style='color:purple'>y.hu@mpie.de</span></a>><br><b>Cc:</b><span class=apple-converted-space> </span><a href="mailto:petsc-users@mcs.anl.gov"><span style='color:purple'>petsc-users@mcs.anl.gov</span></a><br><b>Subject:</b><span class=apple-converted-space> </span>Re: [petsc-users] SNES seems not use my matrix-free operation</span><o:p></o:p></p></div></div></div><div><p class=MsoNormal> <o:p></o:p></p></div><div><div><p class=MsoNormal> <o:p></o:p></p></div></div><div><p class=MsoNormal> "formJacobian" should not be __creating__ the matrices. Here "form" means computing the numerical values in the matrix (or when using a shell matrix it means keeping a copy of X so that your custom matrix-free multiply knows the base location where the matrix free Jacobian-vector products are computed.)<o:p></o:p></p></div><div><div><p class=MsoNormal> <o:p></o:p></p></div></div><div><div><p class=MsoNormal> You create the shell matrices up in your main program and pass them in with SNESSetJacobian(). <o:p></o:p></p></div></div><div><div><p class=MsoNormal> <o:p></o:p></p></div></div><div><div><p class=MsoNormal> Try first calling SNESSetJacobian() to provide the matrices (provide a dummy function argument) and then call DMSNESSetJacobianLocal() to provide your "formjacobian" function (that does not create the matrices).<o:p></o:p></p></div></div><div><div><div><p class=MsoNormal> <o:p></o:p></p></div></div><div><div><p class=MsoNormal> Barry<o:p></o:p></p></div></div><div><div><p class=MsoNormal> <o:p></o:p></p></div></div><div><div><p class=MsoNormal> <o:p></o:p></p></div></div><div><div><p class=MsoNormal> Yes, "form" is a bad word that should not have been used in our code.<o:p></o:p></p></div></div><div><div><p class=MsoNormal> <o:p></o:p></p></div></div><div><div><p class=MsoNormal> <o:p></o:p></p></div><div><div><p class=MsoNormal><br><br><br><o:p></o:p></p></div><blockquote style='margin-top:5.0pt;margin-bottom:5.0pt'><div><div><p class=MsoNormal>On Jan 8, 2024, at 12:24 PM, Yi Hu <<a href="mailto:y.hu@mpie.de"><span style='color:purple'>y.hu@mpie.de</span></a>> wrote:<o:p></o:p></p></div></div><div><p class=MsoNormal> <o:p></o:p></p></div><div><div><div><p class=MsoNormal><span lang=EN-US>Dear PETSc Experts,</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>I am implementing a matrix-free jacobian for my SNES solver in Fortran. (command line option -snes_type newtonls -ksp_type gmres)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>In the main program, I define my residual and jacobian and matrix-free jacobian like the following,</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>…</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>call DMDASNESSetFunctionLocal(DM_mech, INSERT_VALUES, formResidual, PETSC_NULL_SNES, err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>call DMSNESSetJacobianLocal(DM_mech, formJacobian, PETSC_NULL_SNES, err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>…</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'>subroutine formJacobian(residual_subdomain,F,Jac_pre,Jac,dummy,err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> <span class=apple-converted-space> </span></span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'>#include <petsc/finclude/petscmat.h></span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> use petscmat</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> implicit None</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: &</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> residual_subdomain !< DMDA info (needs to be named "in" for macros like XRANGE to work)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> real(pREAL), dimension(3,3,cells(1),cells(2),cells3), intent(in) :: &</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> F !< deformation gradient field</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> Mat :: Jac, Jac_pre</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> PetscObject :: dummy</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> PetscErrorCode :: err_PETSc</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> PetscInt :: N_dof ! global number of DoF, maybe only a placeholder</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> N_dof = 9*product(cells(1:2))*cells3<span class=apple-converted-space> </span></span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> print*, 'in my jac'</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> <span class=apple-converted-space> </span></span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N_dof,N_dof,0,Jac,err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call MatShellSetOperation(Jac,MATOP_MULT,GK_op,err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> <span class=apple-converted-space> </span></span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> print*, 'in my jac'</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> ! for jac preconditioner</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N_dof,N_dof,0,Jac_pre,err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> call MatShellSetOperation(Jac_pre,MATOP_MULT,GK_op,err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> CHKERRQ(err_PETSc)</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> print*, 'in my jac'</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span style='font-size:10.0pt;font-family:"Courier New"'>end subroutine formJacobian</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><pre><span class=hljs-keyword>subroutine</span><span class=line> GK_op(Jac,dF,output,err_PETSc)</span><o:p></o:p></pre><pre><span class=line> </span><span class=hljs-keyword>real</span><span class=line>(pREAL), </span><span class=hljs-keyword>dimension</span><span class=line>(</span><span class=hljs-number>3</span><span class=line>,</span><span class=hljs-number>3</span><span class=line>,cells(</span><span class=hljs-number>1</span><span class=line>),cells(</span><span class=hljs-number>2</span><span class=line>),cells3), </span><span class=hljs-keyword>intent</span><span class=line>(</span><span class=hljs-keyword>in</span><span class=line>) :: &</span><o:p></o:p></pre><pre><span class=line> dF </span><o:p></o:p></pre><pre><span class=line> </span><span class=hljs-keyword>real</span><span class=line>(pREAL), </span><span class=hljs-keyword>dimension</span><span class=line>(</span><span class=hljs-number>3</span><span class=line>,</span><span class=hljs-number>3</span><span class=line>,cells(</span><span class=hljs-number>1</span><span class=line>),cells(</span><span class=hljs-number>2</span><span class=line>),cells3), </span><span class=hljs-keyword>intent</span><span class=line>(</span><span class=hljs-keyword>out</span><span class=line>) :: &</span><o:p></o:p></pre><pre><span class=line> output </span><o:p></o:p></pre><pre><span class=line> </span><span class=hljs-keyword>real</span><span class=line>(pREAL), </span><span class=hljs-keyword>dimension</span><span class=line>(</span><span class=hljs-number>3</span><span class=line>,</span><span class=hljs-number>3</span><span class=line>) :: &</span><o:p></o:p></pre><pre><span class=line> deltaF_aim = </span><span class=hljs-number>0.0_pREAL</span><o:p></o:p></pre><pre> <o:p></o:p></pre><pre><span class=line> Mat :: Jac</span><o:p></o:p></pre><pre><span class=line> PetscErrorCode :: err_PETSc</span><o:p></o:p></pre><pre> <o:p></o:p></pre><pre><span class=line> </span><span class=hljs-keyword>integer</span><span class=line> :: i, j, k, e</span><o:p></o:p></pre><pre> <o:p></o:p></pre><pre><span class=line> </span><span class=line><span lang=EN-US>… a lot of calculations …</span></span><o:p></o:p></pre><pre><span class=line><span lang=EN-US> </span></span><o:p></o:p></pre><pre><span class=line> </span><span class=hljs-builtin>print</span><span class=line>*, </span><span class=hljs-string>'in GK op'</span><o:p></o:p></pre><pre><span class=line> </span><o:p></o:p></pre><pre><span class=hljs-keyword>end</span><span class=line> </span><span class=hljs-keyword>subroutine</span><span class=line> GK_op</span><o:p></o:p></pre><pre> <o:p></o:p></pre><div><div><p class=MsoNormal><span lang=EN-US>The first question is that: it seems I still need to explicitly define the interface of MatCreateShell() and MatShellSetOperation() to properly use them, even though I include them via “use petscmat”. It is a little bit strange to me, since some examples do not perform this step.<span class=apple-converted-space> </span></span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>Then the main issue is that I can build my own Jacobian from my call back function formJacobian, and confirm my Jacobian is a shell matrix (by MatView). However, my customized operator GK_op is not called when solving the nonlinear system (not print my “in GK op”). When I try to monitor my SNES, it gives me some conventional output not mentioning my matrix-free operations. So I guess my customized MATOP_MULT may be not associated with Jacobian. Or my configuration is somehow wrong. Could you help me solve this issue?<span class=apple-converted-space> </span></span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>Thanks,</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US>Yi</span><o:p></o:p></p></div></div><div><div><p class=MsoNormal><span lang=EN-US> </span><o:p></o:p></p></div></div><div><p class=MsoNormal><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'><br><br><br></span><o:p></o:p></p></div><div class=MsoNormal align=center style='text-align:center'><hr size=2 width="100%" align=center></div><div><p class=MsoNormal><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'>-------------------------------------------------<br>Stay up to date and follow us on LinkedIn, Twitter and YouTube.<br><br>Max-Planck-Institut für Eisenforschung GmbH<br>Max-Planck-Straße 1<br>D-40237 Düsseldorf<br> <br>Handelsregister B 2533 <br>Amtsgericht Düsseldorf<br> <br>Geschäftsführung<br>Prof. Dr. Gerhard Dehm<br>Prof. Dr. Jörg Neugebauer<br>Prof. Dr. Dierk Raabe<br>Dr. Kai de Weldige<br> <br>Ust.-Id.-Nr.: DE 11 93 58 514 <br>Steuernummer: 105 5891 1000<br><br><br>Please consider that invitations and e-mails of our institute are <br>only valid if they end with …@mpie.de. <br>If you are not sure of the validity please contact </span><a href="mailto:rco@mpie.de"><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif;color:#954F72'>rco@mpie.de</span></a><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'><br><br>Bitte beachten Sie, dass Einladungen zu Veranstaltungen und E-Mails<br>aus unserem Haus nur mit der Endung …@mpie.de gültig sind. <br>In Zweifelsfällen wenden Sie sich bitte an </span><a href="mailto:rco@mpie.de"><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif;color:#954F72'>rco@mpie.de</span></a><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'><br>-------------------------------------------------</span><o:p></o:p></p></div></div></blockquote></div><div><p class=MsoNormal> <o:p></o:p></p></div></div></div><p class=MsoNormal><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'><br style='caret-color: rgb(0, 0, 0);font-variant-caps: normal;text-align:start;-webkit-text-stroke-width: 0px;word-spacing:0px'><br></span><o:p></o:p></p><div class=MsoNormal align=center style='text-align:center'><hr size=2 width="100%" align=center></div><p class=MsoNormal><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'>-------------------------------------------------<br>Stay up to date and follow us on LinkedIn, Twitter and YouTube.<br><br>Max-Planck-Institut für Eisenforschung GmbH<br>Max-Planck-Straße 1<br>D-40237 Düsseldorf<br> <br>Handelsregister B 2533 <br>Amtsgericht Düsseldorf<br> <br>Geschäftsführung<br>Prof. Dr. Gerhard Dehm<br>Prof. Dr. Jörg Neugebauer<br>Prof. Dr. Dierk Raabe<br>Dr. Kai de Weldige<br> <br>Ust.-Id.-Nr.: DE 11 93 58 514 <br>Steuernummer: 105 5891 1000<br><br><br>Please consider that invitations and e-mails of our institute are <br>only valid if they end with …@mpie.de. <br>If you are not sure of the validity please contact </span><a href="mailto:rco@mpie.de"><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif;color:purple'>rco@mpie.de</span></a><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'><br><br>Bitte beachten Sie, dass Einladungen zu Veranstaltungen und E-Mails<br>aus unserem Haus nur mit der Endung …@mpie.de gültig sind. <br>In Zweifelsfällen wenden Sie sich bitte an </span><a href="mailto:rco@mpie.de"><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif;color:purple'>rco@mpie.de</span></a><span style='font-size:13.5pt;font-family:"Helvetica",sans-serif'><br>-------------------------------------------------</span><o:p></o:p></p></div></blockquote></div><p class=MsoNormal><o:p> </o:p></p></div></div><BR />
<BR />
<HR />
-------------------------------------------------<BR />
Stay up to date and follow us on LinkedIn, Twitter and YouTube.<BR />
<BR />
Max-Planck-Institut für Eisenforschung GmbH<BR />
Max-Planck-Straße 1<BR />
D-40237 Düsseldorf<BR />
<BR />
Handelsregister B 2533 <BR />
Amtsgericht Düsseldorf<BR />
<BR />
Geschäftsführung<BR />
Prof. Dr. Gerhard Dehm<BR />
Prof. Dr. Jörg Neugebauer<BR />
Prof. Dr. Dierk Raabe<BR />
Dr. Kai de Weldige<BR />
<BR />
Ust.-Id.-Nr.: DE 11 93 58 514 <BR />
Steuernummer: 105 5891 1000<BR />
<BR />
<BR />
Please consider that invitations and e-mails of our institute are <BR />
only valid if they end with …@mpie.de. <BR />
If you are not sure of the validity please contact rco@mpie.de<BR />
<BR />
Bitte beachten Sie, dass Einladungen zu Veranstaltungen und E-Mails<BR />
aus unserem Haus nur mit der Endung …@mpie.de gültig sind. <BR />
In Zweifelsfällen wenden Sie sich bitte an rco@mpie.de<BR />
-------------------------------------------------<BR />
</body></html>