<html dir="ltr" adlesseunifierdata="["{\"w\":false,\"id\":\"com.kwizzu.fastesttube\",\"name\":\"FastestTube\",\"isComponentMode\":true}"]">
<head>
<meta http-equiv="Content-Type" content="text/html; charset=iso-8859-1">
<style type="text/css" id="owaParaStyle"></style><style type="text/css"></style><style type="text/css"></style>
</head>
<body fpstyle="1" ocsi="0">
<div style="direction: ltr;font-family: Tahoma;color: #000000;font-size: 10pt;">Dear Matt,
<div><br>
</div>
<div>thanks for your prompt response.</div>
<div><br>
</div>
<div>You are were right (KSPView output below): I was trying to precondition the off-diagonal blocks (I messed up the ordering of the blocks, stored as ((00), (01), (11, 10)) where I should have simply used row-major ordering).</div>
<div><br>
</div>
<div>One question: is PCFieldSplitSetIS is most effective way to go about this? Or is there a way to pass the DMDAs directly rather than extracting the ISs from A?</div>
<div><br>
</div>
<div>Thanks a lot,</div>
<div>Vincent</div>
<div><br>
</div>
<div>
<div>KSP Object: 8 MPI processes</div>
<div>  type: cg</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 PRECONDITIONED norm type for convergence test</div>
<div>PC Object: 8 MPI processes</div>
<div>  type: fieldsplit</div>
<div>    FieldSplit with MULTIPLICATIVE composition: total splits = 3</div>
<div>    Solver info for each split is in the following KSP objects:</div>
<div>    Split number 0 Defined by IS</div>
<div>    KSP Object:    (fieldsplit_0_)     8 MPI processes</div>
<div>      type: cg</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:    (fieldsplit_0_)     8 MPI processes</div>
<div>      type: jacobi</div>
<div>      PC has not been set up so information may be incomplete</div>
<div>      linear system matrix = precond matrix:</div>
<div>      Mat Object:      (fieldsplit_0_)       8 MPI processes</div>
<div>        type: mpiaij</div>
<div>        rows=262144, cols=262144</div>
<div>        total: nonzeros=7.07789e+06, allocated nonzeros=7.07789e+06</div>
<div>        total number of mallocs used during MatSetValues calls =0</div>
<div>    Split number 1 Defined by IS</div>
<div>    KSP Object:    (fieldsplit_1_)     8 MPI processes</div>
<div>      type: cg</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:    (fieldsplit_1_)     8 MPI processes</div>
<div>      type: jacobi</div>
<div>      PC has not been set up so information may be incomplete</div>
<div>      linear system matrix = precond matrix:</div>
<div>      Mat Object:      (fieldsplit_1_)       8 MPI processes</div>
<div>        type: shell</div>
<div>        rows=262144, cols=262144</div>
<div>    Split number 2 Defined by IS</div>
<div>    KSP Object:    (fieldsplit_2_)     8 MPI processes</div>
<div>      type: cg</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:    (fieldsplit_2_)     8 MPI processes</div>
<div>      type: jacobi</div>
<div>      PC has not been set up so information may be incomplete</div>
<div>      linear system matrix = precond matrix:</div>
<div>      Mat Object:      (fieldsplit_2_)       8 MPI processes</div>
<div>        type: shell</div>
<div>        rows=262144, cols=262144</div>
<div>  linear system matrix = precond matrix:</div>
<div>  Mat Object:   8 MPI processes</div>
<div>    type: nest</div>
<div>    rows=786432, cols=786432</div>
<div>      Matrix object: </div>
<div>        type=nest, rows=3, cols=3 </div>
<div>        MatNest structure: </div>
<div>        (0,0) : prefix="fieldsplit_0_", type=mpiaij, rows=262144, cols=262144 </div>
<div>        (0,1) : type=shell, rows=262144, cols=262144 </div>
<div>        (0,2) : type=shell, rows=262144, cols=262144 </div>
<div>        (1,0) : type=mpiaij, rows=262144, cols=262144 </div>
<div>        (1,1) : prefix="fieldsplit_1_", type=shell, rows=262144, cols=262144 </div>
<div>        (1,2) : type=shell, rows=262144, cols=262144 </div>
<div>        (2,0) : type=mpiaij, rows=262144, cols=262144 </div>
<div>        (2,1) : type=shell, rows=262144, cols=262144 </div>
<div>        (2,2) : prefix="fieldsplit_2_", type=shell, rows=262144, cols=262144 </div>
</div>
<div><br>
</div>
<div><br>
</div>
<div>
<div style="font-family: Times New Roman; color: #000000; font-size: 16px">
<hr tabindex="-1">
<div id="divRpF562472" style="direction: ltr;"><font face="Tahoma" size="2" color="#000000"><b>From:</b> Matthew Knepley [knepley@gmail.com]<br>
<b>Sent:</b> Friday, April 29, 2016 7:42 AM<br>
<b>To:</b> Le Chenadec, Vincent<br>
<b>Cc:</b> petsc-users@mcs.anl.gov<br>
<b>Subject:</b> Re: [petsc-users] PCFIELDSPLIT on a MATNEST with MATSHELL off-diagonal blocks<br>
</font><br>
</div>
<div></div>
<div>
<div dir="ltr">
<div class="gmail_extra">
<div class="gmail_quote">On Fri, Apr 29, 2016 at 7:34 AM, Le Chenadec, Vincent <span dir="ltr">
<<a href="mailto:vlc@illinois.edu" target="_blank">vlc@illinois.edu</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex; border-left:1px #ccc solid; padding-left:1ex">
<div>
<div style="direction:ltr; font-family:Tahoma; color:#000000; font-size:10pt">
<div>Dear all,</div>
<div><br>
</div>
<div>I'm trying to use PCFIELDSPLIT as a preconditioner on a block matrix A,</div>
<div><br>
</div>
[A_{00}   A_{01}]
<div>[A_{10}   A_{11}]</div>
<div><br>
</div>
<div>where the diagonal blocks (A_{00} and A_{11}) are diagonally dominant. I would like to use either one of the block preconditioners (not even based on the Schur complement, simply either one of the Jacobi/Gauss-Seidel ones).</div>
</div>
</div>
</blockquote>
<div><br>
</div>
<div>For any questions about solvers, we need to see the output of -ksp_view.</div>
<div><br>
</div>
<div>Below, it seems clear that either a) One of your diagonal blocks is a MatShell, or b) you are trying to invert</div>
<div>an off-diagonal block.</div>
<div><br>
</div>
<div>   Matt</div>
<div> </div>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex; border-left:1px #ccc solid; padding-left:1ex">
<div>
<div style="direction:ltr; font-family:Tahoma; color:#000000; font-size:10pt">
<div>The issue I'm running into concerns the formats I would like to use, so I'll start with that:</div>
<div>
<ol>
<li><span style="font-size:10pt">A_{00} and A_{11} are both in MATMPIAIJ format (created using DMCreateMatrix on 2D DMDAs of different sizes),</span></li><li><span style="font-size:10pt">A_{10} and A_{01} are both in MATSHELL format (created using DMCreateShell). T</span><span style="font-size:10pt">he reason why I used this format is that it seemed convenient, since the 2D DMDAs are of different sizes and I
 assumed I had to provide the MATOP_MULT operation only,</span></li><li><span style="font-size:10pt">Finally the blocks are stored in MATNEST format (MatCreateNest).</span></li></ol>
</div>
<div><span style="font-size:10pt">With this choice, matrix-vector multiplications work fine (the DMDAs being packed in a DMComposite), which makes me think that so far the code is ok. Also, w</span><span style="font-size:10pt">hen I attempt to invert the matrix
 without preconditioner, I get the correct result.</span></div>
<div><span style="font-size:10pt"><br>
</span></div>
<div><span style="font-size:10pt">Now here is where things do not work: </span><span style="font-size:10pt">if I try to use a preconditioner, I get the error pasted below (for a simple PCJACOBI)</span></div>
<div><br>
</div>
<div>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: No support for this operation for this object type</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: Mat type shell</p>
</div>
<div>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #1 MatGetDiagonal() line 4306 in petsc/3.6.3_dbg/src/mat/interface/matrix.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #2 PCSetUp_Jacobi() line 180 in petsc/3.6.3_dbg/src/ksp/pc/impls/jacobi/jacobi.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #3 PCSetUp_Jacobi_NonSymmetric() line 261 in petsc/3.6.3_dbg/src/ksp/pc/impls/jacobi/jacobi.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #4 PCApply_Jacobi() line 286 in petsc/3.6.3_dbg/src/ksp/pc/impls/jacobi/jacobi.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #5 PCApply() line 483 in petsc/3.6.3_dbg/src/ksp/pc/interface/precon.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #6 KSP_PCApply() line 242 in petsc/3.6.3_dbg/include/petsc/private/kspimpl.h</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #7 KSPSolve_CG() line 139 in petsc/3.6.3_dbg/src/ksp/ksp/impls/cg/cg.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #8 KSPSolve() line 604 in petsc/3.6.3_dbg/src/ksp/ksp/interface/itfunc.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #9 PCApply_FieldSplit() line 988 in petsc/3.6.3_dbg/src/ksp/pc/impls/fieldsplit/fieldsplit.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #10 PCApply() line 483 in petsc/3.6.3_dbg/src/ksp/pc/interface/precon.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #11 KSP_PCApply() line 242 in petsc/3.6.3_dbg/include/petsc/private/kspimpl.h</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #12 KSPSolve_CG() line 139 in petsc/3.6.3_dbg/src/ksp/ksp/impls/cg/cg.c</p>
<p style="margin-right:0px; margin-left:0px; font-size:12px; line-height:normal; font-family:Menlo">
[0]PETSC ERROR: #13 KSPSolve() line 604 in petsc/3.6.3_dbg/src/ksp/ksp/interface/itfunc.c</p>
</div>
<div><br>
</div>
<div>I may not be calling the appropriate routines, h<span style="font-size:10pt">ere is what</span><span style="font-size:10pt"> I'm basically doing:</span></div>
<div>
<div>
<ol>
<li><span style="font-size:10pt">KSP type set to KSPCG (A is symmetric),</span><span style="font-size:10pt"> PC type to PCFIELDSPLIT,</span></li><li><span style="font-size:10pt">I then use MatNestGetISs on A to supply the ISs to PCFieldSplitSetIS,</span></li><li><span style="font-size:10pt">Finally I use PCFieldSplitGetSubKSP to use set the subKSP for each diagonal block.</span></li></ol>
</div>
<div><span style="font-size:10pt">As I mentioned, as long as the subKSPs don't have any preconditioning (PCNONE), KSPSolve does not complain, only when I try to using one (even PCJACOBI) do I get similar errors). That </span><span style="font-size:10pt">makes
 me suspect that I'm doing something wrong, </span><span style="font-size:10pt">but I could not find an example that was doing the same thing...</span></div>
<div><span style="font-size:10pt"><br>
</span></div>
<div><span style="font-size:10pt">Any pointers would be much appreciated!</span></div>
</div>
<div><br>
</div>
<div>Thanks,</div>
<div>Vincent</div>
</div>
</div>
</blockquote>
</div>
<br>
<br clear="all">
<div><br>
</div>
-- <br>
<div class="gmail_signature">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
-- Norbert Wiener</div>
</div>
</div>
</div>
</div>
</div>
</div>
</body>
</html>