<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, Apr 7, 2017 at 11:23 AM, Barletta, Ivano <span dir="ltr"><<a href="mailto:ibarletta@inogs.it" target="_blank">ibarletta@inogs.it</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div><div><div>So, as far as I understand, the only benefit of PETSc with symmetric matrices<br></div>is only when Matrix values are set, by reducing the overhead of MatSetValue calls?<br></div></div></div></blockquote><div><br></div><div>It halves the storage. There is a slight advantage from not having to load the lower triangle.</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 dir="ltr"><div><div></div>Thanks,<br></div>Ivano<br></div><div class="gmail_extra"><br><div class="gmail_quote">2017-04-07 17:19 GMT+02:00 Barry Smith <span dir="ltr"><<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>></span>:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><span><br>
> On Apr 7, 2017, at 6:40 AM, Florian Lindner <<a href="mailto:mailinglists@xgm.de" target="_blank">mailinglists@xgm.de</a>> wrote:<br>
><br>
> Hello,<br>
><br>
> two questions about symmetric (MATSBAIJ) matrices.<br>
><br>
> + Entries set with MatSetValue below the main diagonal are ignored. Is that by design?<br>
<br>
</span> Yes<br>
<span><br>
> I rather expected setting A_ij to<br>
> have the same effect as setting A_ji.<br>
<br>
</span> You need to check the relationship between i and j and swap them if needed before the call.<br>
<span><br>
><br>
> + Has MatSetOption to MAT_SYMMETRIC and MAT_SYMMETRIC_ETERNAL any gain on MATSBAIJ matrices?<br>
<br>
</span> No<br>
<div class="m_-5085421560058038031HOEnZb"><div class="m_-5085421560058038031h5"><br>
><br>
> Thanks,<br>
> Florian<br>
><br>
> Test programm:<br>
><br>
><br>
> #include "petscmat.h"<br>
> #include "petscviewer.h"<br>
><br>
> int main(int argc, char **argv)<br>
> {<br>
> PetscInitialize(&argc, &argv, "", NULL);<br>
> PetscErrorCode ierr = 0;<br>
><br>
> Mat A;<br>
> ierr = MatCreate(PETSC_COMM_WORLD, &A); CHKERRQ(ierr);<br>
> MatSetType(A, MATSBAIJ); CHKERRQ(ierr);<br>
> ierr = MatSetSizes(A, 4, 4, PETSC_DECIDE, PETSC_DECIDE); CHKERRQ(ierr);<br>
> ierr = MatSetUp(A); CHKERRQ(ierr);<br>
> ierr = MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE); CHKERRQ(ierr);<br>
> ierr = MatSetOption(A, MAT_SYMMETRY_ETERNAL, PETSC_TRUE); CHKERRQ(ierr);<br>
><br>
> // Stored<br>
> ierr = MatSetValue(A, 1, 2, 21, INSERT_VALUES); CHKERRQ(ierr);<br>
> ierr = MatSetValue(A, 1, 1, 11, INSERT_VALUES); CHKERRQ(ierr);<br>
><br>
> // Ignored<br>
> ierr = MatSetValue(A, 2, 1, 22, INSERT_VALUES); CHKERRQ(ierr);<br>
> ierr = MatSetValue(A, 3, 2, 32, INSERT_VALUES); CHKERRQ(ierr);<br>
> ierr = MatSetValue(A, 3, 1, 31, INSERT_VALUES); CHKERRQ(ierr);<br>
><br>
> ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);<br>
> ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);<br>
><br>
> PetscViewer viewer;<br>
> ierr = PetscViewerCreate(PETSC_COMM_W<wbr>ORLD, &viewer); CHKERRQ(ierr);<br>
> ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII); CHKERRQ(ierr);<br>
> ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_DENSE); CHKERRQ(ierr);<br>
> ierr = MatView(A, viewer); CHKERRQ(ierr);<br>
> ierr = PetscViewerPopFormat(viewer); CHKERRQ(ierr);<br>
> ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);<br>
><br>
> PetscFinalize();<br>
> return 0;<br>
> }<br>
<br>
</div></div></blockquote></div><br></div>
</blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature" data-smartmail="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>