<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Thu, Mar 31, 2016 at 11:03 AM, Matthew Overholt <span dir="ltr"><<a href="mailto:overholt@capesim.com" target="_blank">overholt@capesim.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div lang="EN-US" link="blue" vlink="purple"><div><p class="MsoNormal">Hi,<u></u><u></u></p><p class="MsoNormal"><u></u> <u></u></p><p class="MsoNormal">I am just getting started with PETSc and am having difficulty with setting up my MATSBAIJ matrix.  I’m adapting the ksp/ex23.c example to my 3D FEM calculation; in my case the equation to solve is K*x = b.<u></u><u></u></p><p class="MsoNormal"><u></u> <u></u></p><p class="MsoNormal">For serial execution, the following works fine and gives the correct answer.<u></u><u></u></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatCreate(PETSC_COMM_WORLD,&K);CHKERRQ(ierr);<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetSizes(K,vlocal,vlocal,neqns,neqns);CHKERRQ(ierr);<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetType(K,MATSBAIJ);CHKERRQ(ierr);                // symmetric, block, sparse<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetOption(K,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);    // K is symmetric, positive-definite, sparse<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetOption(K,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);    // so only top tri is needed<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetUp(K);CHKERRQ(ierr);<u></u><u></u></span></p><p class="MsoNormal">where  <span style="font-size:9.0pt;font-family:"Courier New"">vlocal </span>is the result from the call to  <span style="font-size:9.0pt;font-family:"Courier New"">VecGetLocalSize(x,&vlocal) </span>for the solution vector (obviously the full size in the serial case).<u></u><u></u></p><p class="MsoNormal">However, for parallel execution the above crashes on a 11 SEGV Segmentation Violation Error on entry to the function <span style="font-size:9.0pt;font-family:"Courier New"">MatSetOption_MPISBAIJ()</span>, according to TotalView.<u></u><u></u></p><p class="MsoNormal"><u></u> <u></u></p><p class="MsoNormal">So instead I have been trying the following for parallel.<u></u><u></u></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   PetscInt blockSize = 1;         // use a block size of 1 since K is NOT block symmetric<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   PetscInt diagNZ = 5;            // # of non-zeros per row in upper diagonal portion of local submatrix<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   PetscInt offdiagNZ = 8;         // max # of non-zeros per row in off-diagonal portion of local submatrix<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatCreateSBAIJ(PETSC_COMM_WORLD,blockSize,vlocal,vlocal,neqns,neqns,diagNZ,NULL,offdiagNZ,NULL,&K);CHKERRQ(ierr);<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetOption(K,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);    // K is symmetric, positive-definite, sparse<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetOption(K,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);    // so only top tri is needed<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">   ierr = MatSetUp(K);CHKERRQ(ierr);<u></u><u></u></span></p><p class="MsoNormal">However, this fails during the process of setting matrix values (<span style="font-size:9.0pt;font-family:"Courier New"">MatSetValue()</span>) with the following error:<u></u><u></u></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">[0]PETSC ERROR: Argument out of range<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">[0]PETSC ERROR: new nonzero at (0,26) caused a malloc<u></u><u></u></span></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New"">…<u></u><u></u></span></p><p class="MsoNormal">(Note that in this case the matrix size is 96x96 split over 2 processors (vlocal = 48).)<u></u><u></u></p><p class="MsoNormal"><u></u> <u></u></p><p class="MsoNormal">If someone would please point me to the correct way to set up the (MATSBAIJ) matrix for a perfectly symmetric, positive-definite, sparse system, I’d appreciate it.</p></div></div></blockquote><div><br></div><div>I would start by making AIJ work in parallel. The switch to SBAIJ is then fairly easy.</div><div><br></div><div>The error appears to say that you have incorrectly allocated the number of nonzeros in the 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 lang="EN-US" link="blue" vlink="purple"><div><p class="MsoNormal">Thanks,<u></u><u></u></p><p class="MsoNormal">Matt Overholt<u></u><u></u></p><p class="MsoNormal">CapeSym, Inc.<u></u><u></u></p><p class="MsoNormal"><span style="font-size:9.0pt;font-family:"Courier New""><u></u> <u></u></span></p></div><br>
<table style="border-top:1px solid #aaabb6">
        <tbody><tr>
      <td style="width:55px;padding-top:18px"><a href="https://www.avast.com/sig-email?utm_medium=email&utm_source=link&utm_campaign=sig-email&utm_content=emailclient&utm_term=oa-2115-v2-b" target="_blank"><img src="https://ipmcdn.avast.com/images/2016/icons/icon-envelope-tick-round-orange-v1.png"></a></td>
                <td style="width:470px;padding-top:20px;color:#41424e;font-size:13px;font-family:Arial,Helvetica,sans-serif;line-height:18px">Virus-free. <a href="https://www.avast.com/sig-email?utm_medium=email&utm_source=link&utm_campaign=sig-email&utm_content=emailclient&utm_term=oa-2115-v2-b" style="color:#4453ea" target="_blank">www.avast.com</a>
                </td>
        </tr>
</tbody></table></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>