src/ksp/pc/examples/tests/ex8f.F is a simple example, but I made the modification below and it runs fine.<div><div><br></div><div><div>diff --git a/src/ksp/pc/examples/tests/ex8f.F b/src/ksp/pc/examples/tests/ex8f.F</div><div>
--- a/src/ksp/pc/examples/tests/ex8f.F</div><div>+++ b/src/ksp/pc/examples/tests/ex8f.F</div><div>@@ -35,7 +35,7 @@</div><div>       Vec              x,b,u</div><div>       PC               pc</div><div>       PetscInt  n,dim,istart,iend</div>
<div>-      PetscInt  i,j,jj,ii,one,zero</div><div>+      PetscInt  i,j,jj,ii,one,two,zero</div><div>       PetscErrorCode ierr</div><div>       PetscScalar v</div><div>       external         MyResidual</div><div>@@ -50,6 +50,7 @@</div>
<div>       pfive = .5d0</div><div>       n      = 6</div><div>       dim    = n*n</div><div>+      two = 2</div><div>       one    = 1</div><div>       zero   = 0</div><div> </div><div>@@ -138,7 +139,7 @@</div><div>       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)</div>
<div>       call KSPGetPC(ksp,pc,ierr)</div><div>       call PCSetType(pc,PCMG,ierr)</div><div>-      call PCMGSetLevels(pc,one,PETSC_NULL_OBJECT,ierr)</div><div>+      call PCMGSetLevels(pc,two,PETSC_NULL_OBJECT,ierr)</div>
<div>       call PCMGSetResidual(pc,zero,PCMGDefaultResidual,A,ierr)</div><div> </div><div>       call PCMGSetResidual(pc,zero,MyResidual,A,ierr)</div><div>@@ -147,6 +148,7 @@</div><div> !  also serves as the preconditioning matrix.</div>
<div> </div><div>       call KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN,ierr)</div><div>+      call KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD,ierr)</div><div> </div><div> </div><div>       call KSPDestroy(ksp,ierr)</div>
<div><br><div>$ mpirun.hydra -n 2 ./ex8f</div><div>KSP Object: 2 MPI processes</div><div>  type: gmres</div><div>    GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div>
<div>    GMRES: happy breakdown tolerance 1e-30</div><div>  maximum iterations=10000, initial guess is zero</div><div>  tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>  left preconditioning</div><div>
  using DEFAULT norm type for convergence test</div><div>PC Object: 2 MPI processes</div><div>  type: mg</div><div>    MG: type is MULTIPLICATIVE, levels=2 cycles=v</div><div>      Cycles per PCApply=1</div><div>      Not using Galerkin computed coarse grid matrices</div>
<div>  Coarse grid solver -- level -------------------------------</div><div>    KSP Object:    (mg_coarse_)     2 MPI processes</div><div>      type: preonly</div><div>      maximum iterations=1, initial guess is zero</div>
<div>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>      left preconditioning</div><div>      using NONE norm type for convergence test</div><div>    PC Object:    (mg_coarse_)     2 MPI processes</div>
<div>      type: redundant</div><div>        Redundant preconditioner: Not yet setup</div><div>  Down solver (pre-smoother) on level 1 -------------------------------</div><div>    KSP Object:    (mg_levels_1_)     2 MPI processes</div>
<div>      type: chebyshev</div><div>        Chebyshev: eigenvalue estimates:  min = 0, max = 0</div><div>      maximum iterations=2, initial guess is zero</div><div>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div>
<div>      left preconditioning</div><div>      using NONE norm type for convergence test</div><div>    PC Object:    (mg_levels_1_)     2 MPI processes</div><div>      type: sor</div><div>        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1</div>
<div>  Up solver (post-smoother) same as down solver (pre-smoother)</div><div>  linear system matrix = precond matrix:</div><div>  Matrix Object:   2 MPI processes</div><div>    type: mpiaij</div><div>    rows=36, cols=36</div>
<div>    total: nonzeros=156, allocated nonzeros=252</div><div>    total number of mallocs used during MatSetValues calls =0</div><div>      not using I-node (on process 0) routines</div><div><br></div><br></div><br><div class="gmail_quote">
On Fri, Aug 31, 2012 at 7:27 PM, Gaetan Kenway <span dir="ltr"><<a href="mailto:kenway@utias.utoronto.ca" target="_blank">kenway@utias.utoronto.ca</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
Hi<div><br></div><div>I'm having an issue with the PCMG preconditioner in PETSc. I'm trying to use PCMG to precondition a linear system of equations resulting from the discretization of the EUler equations. It is from a multi-block code so I have the ability to easily generate restriction operators using geometric method. </div>


<div><br></div><div>I'm currently trying to setup a simple 2V cycle. For this I use:</div><div><br></div><div><div> ! Create the restriction operator between finest level and one below</div><div> call createRestrictOperator(RL1, 1)</div>


<div>     </div><div> call PCSetFromOptions(pc, PETScIerr)</div><div> nlevels=2</div><div> !call PCMGSetLevels(pc, nlevels, PETSC_NULL_OBJECT,PETScIerr) ! Doesn't work from fortran!!!</div><div> call PCMGSetGalerkin(pc, PETSC_TRUE, PETScIerr)</div>


<div> call PCMGSetRestriction(pc, 1, RL1, PETScierr)</div><div><br></div><div> call PCView(pc, PETScierr)</div></div><div><br></div><div>I run with the command line argument -pc_mg_levels 2 .</div><div><br></div><div>First issue: PCMGSetLevels does not work correctly from Fortran. The traceback I get when I run (with the commented-out line above and without the command line arguments) is:</div>


<div><br></div><div><div>[0]PETSC ERROR: ------------------------------------------------------------------------</div><div>[0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range</div>


<div>[0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger</div><div>[0]PETSC ERROR: or see <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC" target="_blank">http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC</a> ERROR: or try <a href="http://valgrind.org" target="_blank">http://valgrind.org</a> on GNU/linux and Apple Mac OS X to find memory corruption errors</div>


<div>[0]PETSC ERROR: likely location of problem given in stack below</div><div>[0]PETSC ERROR: ---------------------  Stack Frames ------------------------------------</div><div>[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,</div>


<div>[0]PETSC ERROR:       INSTEAD the line number of the start of the function</div><div>[0]PETSC ERROR:       is given.</div><div>[0]PETSC ERROR: [0] PetscCommDuplicate line 144 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/sys/objects/tagm.c</div>


<div>[0]PETSC ERROR: [0] PetscHeaderCreate_Private line 30 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/sys/objects/inherit.c</div><div>[0]PETSC ERROR: [0] KSPCreate line 557 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/ksp/ksp/interface/itcreate.c</div>


<div>[0]PETSC ERROR: [0] PCMGSetLevels line 180 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/ksp/pc/impls/mg/mg.c</div></div><div><br></div><div>I can get around that issue by using -pc_mg_levels 2 and not calling the PCMGSetLevels() command, but I would like to actually use this command and not have set it via manually setting an option. The result of the PCView is listed below which appears to be what I expect (given the default options)</div>


<div><br></div><div><div>PC Object: 1 MPI processes</div><div>  type: mg</div><div>    MG: type is MULTIPLICATIVE, levels=2 cycles=v</div><div>      Cycles per PCApply=1</div><div>      Using Galerkin computed coarse grid matrices</div>


<div>  Coarse grid solver -- level -------------------------------</div><div>    KSP Object:    (mg_coarse_)     1 MPI processes</div><div>      type: preonly</div><div>      maximum iterations=1, initial guess is zero</div>


<div>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>      left preconditioning</div><div>      using NONE norm type for convergence test</div><div>    PC Object:    (mg_coarse_)     1 MPI processes</div>


<div>      type: lu</div><div>        LU: out-of-place factorization</div><div>        tolerance for zero pivot 2.22045e-14</div><div>        matrix ordering: nd</div><div>  Down solver (pre-smoother) on level 1 -------------------------------</div>


<div>    KSP Object:    (mg_levels_1_)     1 MPI processes</div><div>      type: chebyshev</div><div>        Chebyshev: eigenvalue estimates:  min = 0, max = 0</div><div>      maximum iterations=2, initial guess is zero</div>


<div>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>      left preconditioning</div><div>      using NONE norm type for convergence test</div><div>    PC Object:    (mg_levels_1_)     1 MPI processes</div>


<div>      type: sor</div><div>        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1</div><div>  Up solver (post-smoother) same as down solver (pre-smoother)</div><div>  linear system matrix followed by preconditioner matrix:</div>


<div>  Matrix Object:   1 MPI processes</div><div>    type: seqbaij</div><div>    rows=92160, cols=92160, bs=5</div><div>    total: nonzeros=5582400, allocated nonzeros=5582400</div><div>    total number of mallocs used during MatSetValues calls =0</div>


<div>        block size is 5</div><div>  Matrix Object:   1 MPI processes</div><div>    type: seqbaij</div><div>    rows=92160, cols=92160, bs=5</div><div>    total: nonzeros=3089600, allocated nonzeros=3089600</div><div>


    total number of mallocs used during MatSetValues calls =0</div><div>        block size is 5</div></div><div><br></div><div><br></div><div>So that's ok. When it gets to actually solving my system with KSPSolve() I get a returned error code 73 which is "object in argument is in wrong state, e.g. unassembled mat"</div>
</blockquote><div><br></div><div>In general, please don't paraphrase error messages.</div><div><br></div><div>With these problems, try running in a debugger to get a real trace and with valgrind. The most common problem is stack corruption due to passing the wrong number or types of arguments.</div>
<div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">

<div><br></div><div>and the traceback is:</div><div><div>[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,</div><div>[0]PETSC ERROR:       INSTEAD the line number of the start of the function</div>


<div>[0]PETSC ERROR:       is given.</div><div>[0]PETSC ERROR: [0] MatGetVecs line 8142 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/mat/interface/matrix.c</div><div>[0]PETSC ERROR: [0] KSPGetVecs line 774 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/ksp/ksp/interface/iterativ.c</div>


<div>[0]PETSC ERROR: [0] PCSetUp_MG line 508 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/ksp/pc/impls/mg/mg.c</div><div>[0]PETSC ERROR: [0] PCSetUp line 810 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/ksp/pc/interface/precon.c</div>


<div>[0]PETSC ERROR: [0] KSPSetUp line 182 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/ksp/ksp/interface/itfunc.c</div><div>[0]PETSC ERROR: [0] KSPSolve line 351 /nfs/mica/home/kenway/Downloads/petsc-3.3/src/ksp/ksp/interface/itfunc.c</div>


</div><div><br></div><div>Which would indicate to me its an issue with one of the matrices involved, either the linear operator, A, the preconditioner P or my provided restriction matrix RL1. A and P i'm pretty sure are OK, since it works just fine when I use my usual BlockJacobi/Additive Schwartz with ILU(k) on the subdomains. Which leaves the Restriction operator as the issue. I added the following code before PCView()</div>

<div><br></div><div>     call MatGetVecs(RL1, vec1, vec2, ierr)</div><div>     call VecSet(vec1, 1.0, ierr)</div><div>     call MatMult(RL1, vec1, vec2, ierr)</div><div><br></div><div>which ran through just fine (I check each of the ierrs but those are not shown). </div>

<div><br></div><div>So I'm at a loss as to what is causing the error. </div><div>I'm using an up-to-date copy of petsc-3.3 from the repository. </div><div><br></div><div>Thanks,</div><div><br></div><div>Gaetan</div>


</blockquote></div><br></div></div>