[petsc-users] PCFieldSplit with MatNest

Lawrence Mitchell wencel at gmail.com
Wed Mar 13 09:15:27 CDT 2019



> On 13 Mar 2019, at 14:04, Matthew Knepley via petsc-users <petsc-users at mcs.anl.gov> wrote:
> 
> On Wed, Mar 13, 2019 at 9:44 AM Manuel Colera Rico via petsc-users <petsc-users at mcs.anl.gov> wrote:
> Yes:
> 
> [ 0]8416 bytes MatCreateSeqSBAIJWithArrays() line 2431 in 
> /opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c
> [ 0]8416 bytes MatCreateSeqSBAIJWithArrays() line 2431 in 
> /opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c
> [ 0]4544 bytes MatCreateSeqSBAIJWithArrays() line 2431 in 
> /opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c
> [ 0]4544 bytes MatCreateSeqSBAIJWithArrays() line 2431 in 
> /opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c
> 
> Junchao, do imax and ilen get missed in the Destroy with the user provides arrays?
> 
>   https://bitbucket.org/petsc/petsc/src/06a3e802b3873ffbfd04b71a0821522327dd9b04/src/mat/impls/sbaij/seq/sbaij.c#lines-2431

Looks like it. One probably needs something like:

diff --git a/src/mat/impls/sbaij/seq/sbaij.c b/src/mat/impls/sbaij/seq/sbaij.c
index 2b98394140..c39fc696d8 100644
--- a/src/mat/impls/sbaij/seq/sbaij.c
+++ b/src/mat/impls/sbaij/seq/sbaij.c
@@ -2442,6 +2442,7 @@ PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m
   sbaij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
   sbaij->free_a       = PETSC_FALSE;
   sbaij->free_ij      = PETSC_FALSE;
+  sbaij->free_imax_ilen = PETSC_TRUE;
 
   for (ii=0; ii<m; ii++) {
     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];

Lawrence


More information about the petsc-users mailing list