<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, Sep 12, 2014 at 2:56 AM, Florian Lindner <span dir="ltr"><<a href="mailto:mailinglists@xgm.de" target="_blank">mailinglists@xgm.de</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Hello,<br>
<br>
I have a matrix that have the option set MAT_SYMMETRY_ETERNAL and set some values in the upper triangular. When reading values I was expecting that Petsc makes it a symmetric matrix, but the lower triangular is empty like it was initialized.<br></blockquote><div><br></div><div>This is not the point of the functions. The mechanics of symmetry are handled by the matrix type, SBAIJ. This</div><div>option is an optimization hint. It says that the matrix will respond that it is symmetric if asked, which may be</div><div>much much cheaper than checking whether it is actually symmetric.</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">
Thanks,<br>
Florian<br>
<br>
Example code:<br>
<br>
<br>
#include <iostream><br>
#include "petscmat.h"<br>
#include "petscviewer.h"<br>
<br>
using namespace std;<br>
<br>
// Compiling with: mpic++ -g3 -Wall -I ~/software/petsc/include -I ~/software/petsc/arch-linux2-c-debug/include -L ~/software/petsc/arch-linux2-c-debug/lib -lpetsc test.cpp<br>
<br>
int main(int argc, char **args)<br>
{<br>
PetscInitialize(&argc, &args, "", NULL);<br>
<br>
PetscErrorCode ierr = 0;<br>
int N = 4;<br>
Mat matrix;<br>
<br>
// Create dense matrix, but code should work for sparse, too (I hope)<br>
// dense is more convenient to MatView.<br>
ierr = MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, N, N, NULL, &matrix); CHKERRQ(ierr);<br>
ierr = MatSetUp(matrix); CHKERRQ(ierr);<br>
ierr = MatSetOption(matrix, MAT_SYMMETRY_ETERNAL, PETSC_TRUE); CHKERRQ(ierr);<br>
<br>
MatSetValue(matrix, 1, 1, 1, INSERT_VALUES);<br>
MatSetValue(matrix, 1, 2, 2, INSERT_VALUES);<br>
MatSetValue(matrix, 1, 3, 3, INSERT_VALUES);<br>
MatSetValue(matrix, 2, 3, 4, INSERT_VALUES);<br>
<br>
ierr = MatAssemblyBegin(matrix, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);<br>
ierr = MatAssemblyEnd(matrix, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);<br>
<br>
const PetscScalar *vals;<br>
ierr = MatGetRow(matrix, 2, NULL, NULL, &vals);<br>
cout << "Vals = " << vals[0] << " " << vals[1] << " " << vals[2] << " " << vals[3] << endl;<br>
// prints: Vals = 0 0 0 4<br>
// excepted: Vals = 0 2 0 4<br>
ierr = MatRestoreRow(matrix, 2, NULL, NULL, &vals);<br>
<br>
ierr = MatView(matrix, PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(ierr);<br>
<br>
PetscFinalize();<br>
return 0;<br>
}<br>
</blockquote></div><br><br clear="all"><div><br></div>-- <br>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>