[petsc-users] Symmetric matrix: Setting entries below diagonal

Barry Smith bsmith at mcs.anl.gov
Fri Apr 7 11:30:33 CDT 2017


> On Apr 7, 2017, at 11:23 AM, Barletta, Ivano <ibarletta at inogs.it> wrote:
> 
> So, as far as I understand, the only benefit of PETSc with symmetric matrices
> is only when Matrix values are set, by reducing the overhead of MatSetValue calls?

    The benefits of using SBAIJ matrices are

1 You don't need to compute or set values below the diagonal

2) the matrix storage requires about 1/2 the memory since the lower diagonal part is not stored

   If you use AIJ or BAIJ matrices and MatSetOption() to indicate they are symmetric there is, of course no benefit of 1 or 2.

  Barry


> 
> Thanks,
> Ivano
> 
> 2017-04-07 17:19 GMT+02:00 Barry Smith <bsmith at mcs.anl.gov>:
> 
> > On Apr 7, 2017, at 6:40 AM, Florian Lindner <mailinglists at xgm.de> wrote:
> >
> > Hello,
> >
> > two questions about symmetric (MATSBAIJ) matrices.
> >
> > + Entries set with MatSetValue below the main diagonal are ignored. Is that by design?
> 
>    Yes
> 
> > I rather expected setting A_ij to
> > have the same effect as setting A_ji.
> 
>    You need to check the relationship between i and j and swap them if needed before the call.
> 
> >
> > + Has MatSetOption to MAT_SYMMETRIC and MAT_SYMMETRIC_ETERNAL any gain on MATSBAIJ matrices?
> 
>   No
> 
> >
> > Thanks,
> > Florian
> >
> > Test programm:
> >
> >
> > #include "petscmat.h"
> > #include "petscviewer.h"
> >
> > int main(int argc, char **argv)
> > {
> >  PetscInitialize(&argc, &argv, "", NULL);
> >  PetscErrorCode ierr = 0;
> >
> >  Mat A;
> >  ierr = MatCreate(PETSC_COMM_WORLD, &A); CHKERRQ(ierr);
> >  MatSetType(A, MATSBAIJ); CHKERRQ(ierr);
> >  ierr = MatSetSizes(A, 4, 4, PETSC_DECIDE, PETSC_DECIDE); CHKERRQ(ierr);
> >  ierr = MatSetUp(A); CHKERRQ(ierr);
> >  ierr = MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE); CHKERRQ(ierr);
> >  ierr = MatSetOption(A, MAT_SYMMETRY_ETERNAL, PETSC_TRUE); CHKERRQ(ierr);
> >
> >  // Stored
> >  ierr = MatSetValue(A, 1, 2, 21, INSERT_VALUES); CHKERRQ(ierr);
> >  ierr = MatSetValue(A, 1, 1, 11, INSERT_VALUES); CHKERRQ(ierr);
> >
> >  // Ignored
> >  ierr = MatSetValue(A, 2, 1, 22, INSERT_VALUES); CHKERRQ(ierr);
> >  ierr = MatSetValue(A, 3, 2, 32, INSERT_VALUES); CHKERRQ(ierr);
> >  ierr = MatSetValue(A, 3, 1, 31, INSERT_VALUES); CHKERRQ(ierr);
> >
> >  ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
> >  ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
> >
> >  PetscViewer viewer;
> >  ierr = PetscViewerCreate(PETSC_COMM_WORLD, &viewer); CHKERRQ(ierr);
> >  ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII); CHKERRQ(ierr);
> >  ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_DENSE); CHKERRQ(ierr);
> >  ierr = MatView(A, viewer); CHKERRQ(ierr);
> >  ierr = PetscViewerPopFormat(viewer); CHKERRQ(ierr);
> >  ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
> >
> >  PetscFinalize();
> >  return 0;
> > }
> 
> 



More information about the petsc-users mailing list