<div dir="ltr"><div dir="ltr">On Sat, Jul 18, 2020 at 6:11 AM Alex Fleeter <<a href="mailto:luis.saturday@gmail.com">luis.saturday@gmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">I guess using MatNest is still the most memory efficient one, am I correct? Does MatGetLocalSubMatrix demand extra memory for generating the sub matrices?</div></blockquote><div><br></div><div>No. LocalSubMatrix() is just a view into the backing matrix, whether it is AIJ or Nest, and it forwards the MatSetValues() calls.</div><div><br></div><div>   Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Thu, Jul 16, 2020 at 4:44 AM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr">On Thu, Jul 16, 2020 at 6:54 AM Yang Juntao <<a href="mailto:Y.Juntao@hotmail.com" target="_blank">Y.Juntao@hotmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">





<div lang="EN-US">
<div>
<p class="MsoNormal">Hi, Matt, <u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">In case, I want to put blocks of matrices(assembled independently) together as the system matrix just as in the example. What would be your recommended approach?<u></u><u></u></p>
<p class="MsoNormal">I first looked into MatCreateBlockMat, but it only supports sequential and then I was directed to MATNEST.  Do I have to assemble directly on “system matrix”.</p></div></div></blockquote><div><br></div><div>The idea is that you use</div><div><br></div><div>  <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatGetLocalSubMatrix.html" target="_blank">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatGetLocalSubMatrix.html</a></div><div><br></div><div>This does require the extra step of calling</div><div><br></div><div>  <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetLocalToGlobalMapping.html#MatSetLocalToGlobalMapping" target="_blank">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetLocalToGlobalMapping.html#MatSetLocalToGlobalMapping</a></div><div><br></div><div>but you must already have this in order to get the global indices.</div><div><br></div><div>If you use that, then you can assemble locally into AIJ or Nest with the same code.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div lang="EN-US"><div>
<p class="MsoNormal">“<u></u><u></u></p>
<p class="MsoNormal">Matrices of this type are nominally-sparse matrices in which each "entry" is a
<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat" target="_blank">
Mat</a> object. Each <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/Mat.html#Mat" target="_blank">
Mat</a> must have the same size and be sequential. The local and global sizes must be compatible with this decomposition. For matrices containing parallel submatrices and variable block sizes, see
<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATNEST.html#MATNEST" target="_blank">
MATNEST</a>. <u></u><u></u></p>
<p class="MsoNormal">“<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Regards<u></u><u></u></p>
<p class="MsoNormal">Juntao<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<div style="border-style:solid none none;border-top-width:1pt;border-top-color:rgb(225,225,225);padding:3pt 0in 0in">
<p class="MsoNormal"><b>From:</b> Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> <br>
<b>Sent:</b> Thursday, July 16, 2020 5:56 PM<br>
<b>To:</b> Karl Yang <<a href="mailto:y.juntao@hotmail.com" target="_blank">y.juntao@hotmail.com</a>><br>
<b>Cc:</b> Junchao Zhang <<a href="mailto:junchao.zhang@gmail.com" target="_blank">junchao.zhang@gmail.com</a>>; <a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a><br>
<b>Subject:</b> Re: [petsc-users] Errors encountered with KSPsolver with MATNEST and VECNEST<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<div>
<p class="MsoNormal">On Thu, Jul 16, 2020 at 5:43 AM Karl Yang <<a href="mailto:y.juntao@hotmail.com" target="_blank">y.juntao@hotmail.com</a>> wrote:<u></u><u></u></p>
</div>
<div>
<blockquote style="border-style:none none none solid;border-left-width:1pt;border-left-color:rgb(204,204,204);padding:0in 0in 0in 6pt;margin-left:4.8pt;margin-right:0in">
<div>
<p class="MsoNormal">Hi, Matt, <u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Sorry for the confusion. The error shows as below without SetIS().<u></u><u></u></p>
</div>
</blockquote>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">I see the problem now. It is about the interface for MatNest.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">MatNest is intended _only_ as an optimization. I think you should never have MatNest in user code. Here the<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">problem is that MatNest only supports extraction using the ISes it was created with, whereas detection is<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">finding a different IS. I recommend getting everything working with AIJ, and then use MatSetType() or<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">-mat_type at the end to try out Nest. The interfaces are the same.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">  Thanks,<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">    Matt<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<blockquote style="border-style:none none none solid;border-left-width:1pt;border-left-color:rgb(204,204,204);padding:0in 0in 0in 6pt;margin-left:4.8pt;margin-right:0in">
<div>
<p class="MsoNormal">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Arguments are incompatible<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Could not find index set<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020 <u></u>
<u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named ae3f0ab5a8c4 by Unknown Thu Jul 16 17:36:52 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #1 MatNestFindIS() line 365 in /usr/local/petsc/petsc-3.12.5/src/mat/impls/nest/matnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #2 MatNestFindSubMat() line 425 in /usr/local/petsc/petsc-3.12.5/src/mat/impls/nest/matnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #3 MatCreateSubMatrix_Nest() line 456 in /usr/local/petsc/petsc-3.12.5/src/mat/impls/nest/matnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #4 MatCreateSubMatrix() line 7859 in /usr/local/petsc/petsc-3.12.5/src/mat/interface/matrix.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #5 PCSetUp_FieldSplit() line 676 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #6 PCSetUp() line 894 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Corrupt argument: <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind</a><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Invalid Pointer to Object: Parameter # 1<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020 <u></u>
<u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named ae3f0ab5a8c4 by Unknown Thu Jul 16 17:36:52 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #7 MatDestroy() line 1266 in /usr/local/petsc/petsc-3.12.5/src/mat/interface/matrix.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #8 PCSetUp_FieldSplit() line 697 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #9 PCSetUp() line 894 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #10 KSPSetUp() line 376 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">KSP Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  type: fgmres<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    happy breakdown tolerance 1e-30<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  maximum iterations=10000, initial guess is zero<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  right preconditioning<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  using UNPRECONDITIONED norm type for convergence test<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">PC Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  type: fieldsplit<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  PC has not been set up so information may be incomplete<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    FieldSplit with Schur preconditioner, blocksize = 1, factorization FULL<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Preconditioner for the Schur complement formed from S itself<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Split info:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Split number 0 Defined by IS<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Split number 1 Defined by IS<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    KSP solver for A00 block<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      KSP Object: (fieldsplit_0_) 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        type: preonly<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        maximum iterations=10000, initial guess is zero<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        left preconditioning<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        using DEFAULT norm type for convergence test<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      PC Object: (fieldsplit_0_) 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        type not yet set<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        PC has not been set up so information may be incomplete<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    KSP solver for upper A00 in upper triangular factor <u></u>
<u></u></p>
</div>
<div>
<p class="MsoNormal">        not yet available<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    KSP solver for S = A11 - A10 inv(A00) A01 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        not yet available<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  linear system matrix = precond matrix:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  Mat Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    type: nest<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    rows=15, cols=15<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      Matrix object: <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        type=nest, rows=3, cols=3 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        MatNest structure: <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (0,0) : type=seqaij, rows=5, cols=5 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (0,1) : NULL <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (0,2) : type=seqaij, rows=5, cols=5 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (1,0) : NULL <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (1,1) : type=seqaij, rows=5, cols=5 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (1,2) : type=seqaij, rows=5, cols=5 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (2,0) : type=seqaij, rows=5, cols=5 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (2,1) : type=seqaij, rows=5, cols=5 <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">        (2,2) : NULL <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">Segmentation fault (core dumped)<u></u><u></u></p>
</div>
<p class="MsoNormal" style="margin-bottom:12pt"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Attached is the test code. <u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Regards<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">JT<u></u><u></u></p>
</div>
<p class="MsoNormal" style="margin-bottom:12pt"><br>
<br>
<br>
<u></u><u></u></p>
<div>
<p class="MsoNormal">On Jul 16 2020, at 5:30 pm, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:<u></u><u></u></p>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<div>
<div>
<p class="MsoNormal">On Thu, Jul 16, 2020 at 1:23 AM Karl Yang <<a href="https://link.getmailspring.com/link/1D2D2990-B215-4510-B196-F29BD06D4684@getmailspring.com/0?redirect=mailto%3Ay.juntao%40hotmail.com&recipient=a25lcGxleUBnbWFpbC5jb20%3D" title="mailto:y.juntao@hotmail.com" target="_blank">y.juntao@hotmail.com</a>>
 wrote:<u></u><u></u></p>
</div>
</div>
<p class="MsoNormal"><img border="0" width="1" height="1" style="width: 0.0104in; height: 0.0104in;" id="gmail-m_-8310778005465918319gmail-m_-7740916875831783563gmail-m_8175375066635148063_x0000_i1025" src="https://link.getmailspring.com/open/1D2D2990-B215-4510-B196-F29BD06D4684@getmailspring.com?me=05842e25&recipient=a25lcGxleUBnbWFpbC5jb20%3D" alt="Sent from Mailspring"><u></u><u></u></p>
<div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal">Hi, Junchao,<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Thank you. It works. :)<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">May I ask some follow up question regarding the same example about fieldsplit. I am trying to learn about matnest with fieldsplit preconditioner on manual page 94. But I get a bit lost trying to understand fields, ISs once it is not eactly
 the same on the manual.<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">In this same example, I tried to use PCFIELDSPLIT and  PCFieldSplitSetDetectSaddlePoint as follows.<u></u><u></u></p>
</div>
</blockquote>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">You cannot mix DetectSaddlePoint() with SetIS() since the detection tries to make the proper ISes. Maybe we should check that.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">  Thanks,<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">     Matt<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal">IS rows[3];<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">IS cols[3];<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">MatNestGetISs(testMesh->system, rows, cols);<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">KSP ksp;<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">PC  pc;<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">KSPCreate(PETSC_COMM_WORLD, &ksp);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">KSPSetType(ksp, KSPFGMRES);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">KSPSetOperators(ksp, testMesh->system, testMesh->system);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">KSPGetPC(ksp, &pc);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">PCFieldSplitSetIS(pc, "u", rows[0]);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">PCFieldSplitSetIS(pc, "v", rows[1]);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">PCFieldSplitSetIS(pc, "p", rows[2]);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">PCSetType(pc, PCFIELDSPLIT);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">PCFieldSplitSetDetectSaddlePoint(pc, PETSC_TRUE);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">KSPSetUp(ksp);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">KSPSolve(ksp, testMesh->rhs, testMesh->solution);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> <u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">But I will get the following error,<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Arguments are incompatible<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Could not find index set<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" title="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named f601294ac804 by Unknown Thu Jul 16 13:08:39 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #1 MatNestFindIS() line 365 in /usr/local/petsc/petsc-3.12.5/src/mat/impls/nest/matnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #2 MatNestFindSubMat() line 425 in /usr/local/petsc/petsc-3.12.5/src/mat/impls/nest/matnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #3 MatCreateSubMatrix_Nest() line 456 in /usr/local/petsc/petsc-3.12.5/src/mat/impls/nest/matnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #4 MatCreateSubMatrix() line 7859 in /usr/local/petsc/petsc-3.12.5/src/mat/interface/matrix.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #5 PCSetUp_FieldSplit() line 676 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #6 PCSetUp() line 894 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #7 KSPSetUp() line 376 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Corrupt argument: <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind" title="https://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind</a><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Invalid Pointer to Object: Parameter # 1<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" title="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named f601294ac804 by Unknown Thu Jul 16 13:08:39 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #8 MatDestroy() line 1266 in /usr/local/petsc/petsc-3.12.5/src/mat/interface/matrix.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #9 PCSetUp_FieldSplit() line 697 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #10 PCSetUp() line 894 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #11 KSPSetUp() line 376 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #12 KSPSolve() line 703 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">I read about the examples on the manual and mostly are based on 2<em><span style="font-family:Calibri,sans-serif">2 blocks but still cannot figure out how pcfieldsplit should be set in this 3</span></em> *<em><span style="font-family:Calibri,sans-serif">3
 blocks, and how IS should be set accordingly. I read about -fieldsplit_1_pc_type none option suggested on the documentation, but I am not sure is it still valid in a 33 blocks.</span></em><u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal"><em><span style="font-family:Calibri,sans-serif">Regards</span></em><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><em><span style="font-family:Calibri,sans-serif">JT</span></em><u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">On Jul 16 2020, at 6:57 am, Junchao Zhang <<a href="mailto:junchao.zhang@gmail.com" title="mailto:junchao.zhang@gmail.com" target="_blank">junchao.zhang@gmail.com</a>> wrote:<u></u><u></u></p>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<div>
<p class="MsoNormal">Adding the MatNestSetVecType() line fixed the problem.  But I don't know why petsc had this weird design originally.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">      MatCreateNest(PETSC_COMM_WORLD, 3, NULL, 3, NULL, mblock, &testMesh->system);<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"> +   MatNestSetVecType(testMesh->system,VECNEST);<u></u><u></u></p>
</div>
<div>
<div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">--Junchao Zhang<u></u><u></u></p>
</div>
</div>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<div>
<div>
<p class="MsoNormal">On Wed, Jul 15, 2020 at 2:16 PM Karl Yang <<a href="mailto:y.juntao@hotmail.com" title="mailto:y.juntao@hotmail.com" target="_blank">y.juntao@hotmail.com</a>> wrote:<u></u><u></u></p>
</div>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal">Hi, Barry,<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">I created a minimum example to reproduce the error. The example code is attached.<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Mat Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  type: nest<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  Matrix object:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    type=nest, rows=3, cols=3<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    MatNest structure:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (0,0) : type=seqaij, rows=5, cols=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (0,1) : NULL<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (0,2) : type=seqaij, rows=5, cols=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (1,0) : NULL<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (1,1) : type=seqaij, rows=5, cols=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (1,2) : type=seqaij, rows=5, cols=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (2,0) : type=seqaij, rows=5, cols=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (2,1) : type=seqaij, rows=5, cols=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (2,2) : NULL<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">Vec Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  type: nest<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  VecNest, rows=3,  structure:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  (0) : type=seq, rows=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Vec Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      type: seq<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  (1) : type=seq, rows=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Vec Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      type: seq<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  (2) : type=seq, rows=5<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Vec Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      type: seq<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    1.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Invalid argument<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Nest vector argument 2 not setup.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" title="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named f601294ac804 by Unknown Thu Jul 16 03:11:02 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #1 VecCopy_Nest() line 68 in /usr/local/petsc/petsc-3.12.5/src/vec/vec/impls/nest/vecnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #2 VecCopy() line 1577 in /usr/local/petsc/petsc-3.12.5/src/vec/vec/interface/vector.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #3 KSPInitialResidual() line 61 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itres.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #4 KSPSolve_GMRES() line 236 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/gmres.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #5 KSPSolve() line 760 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">On Jul 16 2020, at 2:39 am, Barry Smith <<a href="mailto:bsmith@petsc.dev" title="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> wrote:<u></u><u></u></p>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">  Could be a bug in PETSc, please send a sample program that produces the problem.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">  Barry<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal">On Jul 15, 2020, at 1:28 PM, Karl Yang <<a href="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/0?redirect=mailto%3Ay.juntao%40hotmail.com&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" title="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/0?redirect=mailto%3Ay.juntao%40hotmail.com&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" target="_blank">y.juntao@hotmail.com</a>>
 wrote:<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<div>
<p class="MsoNormal">Hi, Barry,<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Thanks for your reply. I output the matrix and vec with matview and vecview. And I think I do have the values set as it has been shown.<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Mat Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  type: nest<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  Matrix object:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    type=nest, rows=3, cols=3<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    MatNest structure:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (0,0) : type=seqaij, rows=25, cols=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (0,1) : NULL<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (0,2) : type=seqaij, rows=25, cols=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (1,0) : NULL<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (1,1) : type=seqaij, rows=25, cols=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (1,2) : type=seqaij, rows=25, cols=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (2,0) : type=seqaij, rows=25, cols=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (2,1) : type=seqaij, rows=25, cols=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    (2,2) : NULL<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">Vec Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  type: nest<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  VecNest, rows=3,  structure:<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  (0) : type=seq, rows=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Vec Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      type: seq<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -51.0776<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -90.8404<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    58.7011<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    90.8404<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    51.0776<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  (1) : name="Vec_0x84000000_0", type=seq, rows=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Vec Object: Vec_0x84000000_0 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      type: seq<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -51.0776<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -90.8404<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    58.7011<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    90.8404<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    -75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    51.0776<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    75.8279<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">  (2) : type=seq, rows=25<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    Vec Object: 1 MPI processes<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">      type: seq<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">    0.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Invalid argument<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Nest vector argument 4 not setup.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: See <a href="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/1?redirect=https%3A%2F%2Fwww.mcs.anl.gov%2Fpetsc%2Fdocumentation%2Ffaq.html&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" title="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/1?redirect=https%3A%2F%2Fwww.mcs.anl.gov%2Fpetsc%2Fdocumentation%2Ffaq.html&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named f601294ac804 by Unknown Thu Jul 16 02:15:26 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #1 VecPointwiseMult_Nest() line 245 in /usr/local/petsc/petsc-3.12.5/src/vec/vec/impls/nest/vecnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #2 VecPointwiseMult() line 1107 in /usr/local/petsc/petsc-3.12.5/src/vec/vec/interface/vector.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #3 PCApply_Jacobi() line 272 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/jacobi/jacobi.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #4 PCApply() line 444 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #5 KSP_PCApply() line 281 in /usr/local/petsc/petsc-3.12.5/include/petsc/private/kspimpl.h<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #6 KSPInitialResidual() line 65 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itres.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #7 KSPSolve_GMRES() line 236 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/gmres.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #8 KSPSolve() line 760 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">On Jul 16 2020, at 2:13 am, Barry Smith <<a href="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/2?redirect=mailto%3Absmith%40petsc.dev&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" title="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/2?redirect=mailto%3Absmith%40petsc.dev&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" target="_blank">bsmith@petsc.dev</a>>
 wrote:<u></u><u></u></p>
</div>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">  From the error message the second vector you are using to form the nest vector has not been completely created/set up. Make sure you called either VecCreateSeq(), VecCreateMPI() or did VecCreate(), VecSetType(), VecSetUp() and if it is
 the right hand side make sure you put numerical values in it.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">   Barry<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
<blockquote style="margin-top:5pt;margin-bottom:5pt">
<div>
<p class="MsoNormal">On Jul 15, 2020, at 1:05 PM, Karl Yang <<a href="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/3?redirect=https%3A%2F%2Flink.getmailspring.com%2Flink%2F540F752A-DB59-4878-9C98-B8F4B4D9FD0D%40getmailspring.com%2F0%3Fredirect%3Dmailto%253Ay.juntao%2540hotmail.com%26recipient%3DYnNtaXRoQHBldHNjLmRldg%253D%253D&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" title="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/3?redirect=https%3A%2F%2Flink.getmailspring.com%2Flink%2F540F752A-DB59-4878-9C98-B8F4B4D9FD0D%40getmailspring.com%2F0%3Fredirect%3Dmailto%253Ay.juntao%2540hotmail.c" target="_blank">y.juntao@hotmail.com</a>>
 wrote:<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<div>
<p class="MsoNormal">Hello,<u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal">I'm trying to solve a 2*2 block matrix system.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">I managed to form a 2*2 MatNest and a 21 VecNest. However when I try to  make use of ksp solve, I encounter the following error. There is only one or two examples on MatNest and VecNest. I could not figure out what could be the trouble.<u></u><u></u></p>
</div>
<p class="MsoNormal" style="margin-bottom:12pt"><u></u> <u></u></p>
<div>
<p class="MsoNormal">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Invalid argument<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Nest vector argument 2 not setup.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: See <a href="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/4?redirect=https%3A%2F%2Flink.getmailspring.com%2Flink%2F540F752A-DB59-4878-9C98-B8F4B4D9FD0D%40getmailspring.com%2F1%3Fredirect%3Dhttps%253A%252F%252Fwww.mcs.anl.gov%252Fpetsc%252Fdocumentation%252Ffaq.html%26recipient%3DYnNtaXRoQHBldHNjLmRldg%253D%253D&recipient=cGV0c2MtdXNlcnNAbWNzLmFubC5nb3Y%3D" title="https://link.getmailspring.com/link/060AB2C3-F3D7-408C-99CF-C4D1B27DD59F@getmailspring.com/4?redirect=https%3A%2F%2Flink.getmailspring.com%2Flink%2F540F752A-DB59-4878-9C98-B8F4B4D9FD0D%40getmailspring.com%2F1%3Fredirect%3Dhttps%253A%252F%252Fwww.mcs.anl.g" target="_blank">
https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named f601294ac804 by Unknown Thu Jul 16 01:56:11 2020<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #1 VecCopy_Nest() line 68 in /usr/local/petsc/petsc-3.12.5/src/vec/vec/impls/nest/vecnest.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #2 VecCopy() line 1577 in /usr/local/petsc/petsc-3.12.5/src/vec/vec/interface/vector.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #3 KSPInitialResidual() line 61 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itres.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #4 KSPSolve_GMRES() line 236 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/gmres.c<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">[0]PETSC ERROR: #5 KSPSolve() line 760 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c<u></u><u></u></p>
</div>
<p class="MsoNormal" style="margin-bottom:12pt"><u></u> <u></u></p>
<div>
<p class="MsoNormal">Regards<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">JT<u></u><u></u></p>
</div>
</div>
</blockquote>
</div>
</div>
</blockquote>
</div>
</blockquote>
</div>
</div>
</blockquote>
</blockquote>
</div>
</blockquote>
</blockquote>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal">--<u></u><u></u></p>
</div>
<div>
<div>
<div>
<div>
<div>
<div>
<div>
<div>
<p class="MsoNormal">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal">-- Norbert Wiener<u></u><u></u></p>
</div>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<div>
<p class="MsoNormal"><a href="https://link.getmailspring.com/link/1D2D2990-B215-4510-B196-F29BD06D4684@getmailspring.com/1?redirect=http%3A%2F%2Fwww.cse.buffalo.edu%2F~knepley%2F&recipient=a25lcGxleUBnbWFpbC5jb20%3D" title="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><u></u><u></u></p>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</blockquote>
</div>
<p class="MsoNormal"><br clear="all">
<u></u><u></u></p>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<p class="MsoNormal">-- <u></u><u></u></p>
<div>
<div>
<div>
<div>
<div>
<div>
<div>
<p class="MsoNormal">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<u></u><u></u></p>
</div>
<div>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
<div>
<p class="MsoNormal"><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><u></u><u></u></p>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>

</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>