<div dir="ltr"><div>Hello,</div><div><br></div><div>I am teaching myself how to use Petsc for nonlinear equations and I've run into a problem that I can't quite figure out. I am trying to use the matrix coloring routines for the finite difference Jacobian approximation, and I've followed the steps in the manual to do this. </div><div>When I run the program with a MG preconditioner, I get back the error:</div><div><blockquote style="margin:0 0 0 40px;border:none;padding:0px"><div>[0]PETSC ERROR: --------------------- Error Message ----------------------------<br>[0]PETSC ERROR: Arguments are incompatible<br>[0]PETSC ERROR: Zero diagonal on row 0<br></div><div>.....</div></blockquote><br></div><div>What's interesting is that after I've added non-zero entries to the matrix with MatrixSetValues() and assembled the matrix with MatAssemblyBegin() + MatAssemblyEnd(), I can verify that every diagonal entry is non-zero with a call to MatGetValues. I've included a relevant code snippet below and I'm happy to send more. Any guidance is greatly appreciated.</div><div><br></div><div>Command line:</div><blockquote style="margin:0 0 0 40px;border:none;padding:0px"><div>mpirun -n 1 ./program -snes_view -snes_converged_reason -snes_monitor -ksp_monitor -ksp_converged_reason -pc_type mg</div></blockquote><div>Code snippet:<blockquote style="margin:0 0 0 40px;border:none;padding:0px"><div> ierr = FormJacobianColoring(snes,J);CHKERRQ(ierr); // this function set's matrix values and assembles the matrix<br><br>// removed the code, but this is where I've used MatGetValues() to ensure that the diagonal of J (as well as other entries) has been set to 1.0<br><br> ierr = MatColoringCreate(J,&coloring);CHKERRQ(ierr);<br>      ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);<br>      ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);<br>     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);<br>      ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);<br>       /*<br>            Create the data structure that SNESComputeJacobianDefaultColor() uses to compute the actual Jacobians via finite differences.<br> */<br>    ierr = MatFDColoringCreate(J,iscoloring,&fdcoloring);CHKERRQ(ierr);<br>       ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))FormFunction,NULL);CHKERRQ(ierr);<br>        ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);<br> ierr = MatFDColoringSetUp(J,iscoloring,fdcoloring);CHKERRQ(ierr);<br>     ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring);CHKERRQ(ierr);<br><br>  ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);<br><br>      ierr = VecSet(x,0.001);CHKERRQ(ierr);<br> ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);<br></div></blockquote><br></div><div><br></div><div>And here's the full error code:</div><div><blockquote style="margin:0 0 0 40px;border:none;padding:0px"><div>[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>[0]PETSC ERROR: Arguments are incompatible<br>[0]PETSC ERROR: Zero diagonal on row 0<br>[0]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html">http://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br>[0]PETSC ERROR: Petsc Release Version 3.10.3, Dec, 18, 2018<br>[0]PETSC ERROR: ./petsc_flowroute on a  named i0n22 by alexprescott Thu Apr  2 20:14:01 2020<br>[0]PETSC ERROR: Configure options --prefix=/cm/shared/uaapps/petsc/3.10.3 --download-fblaslapack --download-metis --download-parmetis --download-hypre PETSC_ARCH=linux-gnu --with-debugging=no COPTFLAGS=-O3 CXXOPTFLAGS=-O3<br>[0]PETSC ERROR: #1 MatInvertDiagonal_SeqAIJ() line 1662 in /cm/local/uabuild/petsc/petsc-3.10.3/src/mat/impls/aij/seq/aij.c<br>[0]PETSC ERROR: #2 MatSOR_SeqAIJ() line 1693 in /cm/local/uabuild/petsc/petsc-3.10.3/src/mat/impls/aij/seq/aij.c<br>[0]PETSC ERROR: #3 MatSOR() line 3932 in /cm/local/uabuild/petsc/petsc-3.10.3/src/mat/interface/matrix.c<br>[0]PETSC ERROR: #4 PCApply_SOR() line 31 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/pc/impls/sor/sor.c<br>[0]PETSC ERROR: #5 PCApply() line 462 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/pc/interface/precon.c<br>[0]PETSC ERROR: #6 KSP_PCApply() line 281 in /cm/local/uabuild/petsc/petsc-3.10.3/include/petsc/private/kspimpl.h<br>[0]PETSC ERROR: #7 KSPInitialResidual() line 67 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/interface/itres.c<br>[0]PETSC ERROR: #8 KSPSolve_GMRES() line 233 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/impls/gmres/gmres.c<br>[0]PETSC ERROR: #9 KSPSolve() line 780 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/interface/itfunc.c<br>[0]PETSC ERROR: #10 KSPSolve_Chebyshev() line 367 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/impls/cheby/cheby.c<br>[0]PETSC ERROR: #11 KSPSolve() line 780 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/interface/itfunc.c<br>[0]PETSC ERROR: #12 PCMGMCycle_Private() line 20 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/pc/impls/mg/mg.c<br>[0]PETSC ERROR: #13 PCApply_MG() line 377 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/pc/impls/mg/mg.c<br>[0]PETSC ERROR: #14 PCApply() line 462 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/pc/interface/precon.c<br>[0]PETSC ERROR: #15 KSP_PCApply() line 281 in /cm/local/uabuild/petsc/petsc-3.10.3/include/petsc/private/kspimpl.h<br>[0]PETSC ERROR: #16 KSPInitialResidual() line 67 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/interface/itres.c<br>[0]PETSC ERROR: #17 KSPSolve_GMRES() line 233 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/impls/gmres/gmres.c<br>[0]PETSC ERROR: #18 KSPSolve() line 780 in /cm/local/uabuild/petsc/petsc-3.10.3/src/ksp/ksp/interface/itfunc.c<br>[0]PETSC ERROR: #19 SNESSolve_NEWTONLS() line 224 in /cm/local/uabuild/petsc/petsc-3.10.3/src/snes/impls/ls/ls.c<br>[0]PETSC ERROR: #20 SNESSolve() line 4397 in /cm/local/uabuild/petsc/petsc-3.10.3/src/snes/interface/snes.c<br>[0]PETSC ERROR: #21 main() line 705 in /home/u16/alexprescott/petsc_example/flowroute/petsc_flowroute.c<br>[0]PETSC ERROR: PETSc Option Table entries:<br>[0]PETSC ERROR: -ksp_converged_reason<br>[0]PETSC ERROR: -ksp_monitor<br>[0]PETSC ERROR: -pc_type mg<br>[0]PETSC ERROR: -snes_converged_reason<br>[0]PETSC ERROR: -snes_monitor<br>[0]PETSC ERROR: -snes_view<br>[0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint@mcs.anl.gov----------<br>application called MPI_Abort(MPI_COMM_WORLD, 75) - process 0<br></div></blockquote></div><div><br></div><div>Printed matrix for a 25 x 25 example</div><div><blockquote style="margin:0 0 0 40px;border:none;padding:0px"><div>1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 0 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>1 0 0 0 0 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 1 0 0 0 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 0 1 0 0 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 0 0 1 0 0 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 0 0 0 1 0 0 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>1 0 0 0 0 1 0 0 0 0 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 1 0 0 0 0 1 0 0 0 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0<br>0 0 1 0 0 0 0 1 0 0 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0<br>0 0 0 1 0 0 0 0 1 0 0 1 1 1 0 0 0 0 0 0 0 0 0 0 0<br>0 0 0 0 1 0 0 0 0 1 0 0 1 1 1 0 0 0 0 0 0 0 0 0 0<br>0 0 0 0 0 1 0 0 0 0 1 0 0 0 0 1 0 0 0 0 0 0 0 0 0<br>0 0 0 0 0 0 1 0 0 0 0 1 0 0 0 1 1 0 0 0 0 0 0 0 0<br>0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 1 1 1 0 0 0 0 0 0 0<br>0 0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 1 1 1 0 0 0 0 0 0<br>0 0 0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 1 1 1 0 0 0 0 0<br>0 0 0 0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 0 0 1 0 0 0 0<br>0 0 0 0 0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 0 1 1 0 0 0<br>0 0 0 0 0 0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 1 1 1 0 0<br>0 0 0 0 0 0 0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 1 1 1 0<br>0 0 0 0 0 0 0 0 0 0 0 0 0 0 1 0 0 0 0 1 0 0 1 1 1<br></div></blockquote></div><div><br></div><div><br></div><div>Best,</div><div>Alexander</div><div><br></div><div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature" data-smartmail="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div style="font-size:12.8px"><span style="font-family:arial,sans-serif">Alexander</span><span style="font-family:arial,sans-serif"> </span><span style="font-family:arial,sans-serif">Prescott</span><br></div><div style="font-size:12.8px"><span style="font-family:arial,sans-serif"><a href="mailto:alexprescott@email.arizona.edu" target="_blank">alexprescott@email.arizona.edu</a></span></div><div style="font-size:12.8px"><div style="font-family:arial,sans-serif"><span>PhD</span> <span>Candidate</span>, The University of Arizona</div><div style="font-family:arial,sans-serif">Department of Geosciences</div><div style="font-family:arial,sans-serif">1040 E. 4th Street</div><div style="font-family:arial,sans-serif">Tucson, AZ, 85721</div></div></div></div></div></div></div></div></div></div>