[petsc-dev] error on MatCreateTranspose with MatSchur

Stephan Kramer s.kramer at imperial.ac.uk
Fri Feb 28 11:09:14 CST 2014


Dear petsc devs,

Building our code against today's 'master', I came across the following error that worked fine with petsc 3.4.

In our code we create a MatSchur matrix for which the A01 is created with MatCreateTranspose (as it's the transpose of A10). So we have:

! assembly of A00, Ap00, A10 and A11 in the normal way
call MatCreateTranspose(A10, A01, ierr)
call MatCreateSchurComplement(A00, Ap00, A01, A10, A11, N, ierr)

Then in the next KSPSolve with the schur complement matrix N, I get the following error:


[0]PETSC ERROR: --------------------- Error Message ------------------------------------
[0]PETSC ERROR: Nonconforming object sizes!
[0]PETSC ERROR: Mat mat,Vec x: global dim -1 5!
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Development GIT revision: v3.4.3-3444-g2480567  GIT Date: 2014-02-28 11:14:27 +0100
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: ../../bin/fluidity on a arch-linux2-c-debug named stommel by skramer Fri Feb 28 16:43:06 2014
[0]PETSC ERROR: Libraries linked from /data/stephan/git/petsc/arch-linux2-c-debug/lib
[0]PETSC ERROR: Configure run at Fri Feb 28 15:20:21 2014
[0]PETSC ERROR: Configure options --download-f-blas-lapack=1 --download-blacs=1 --download-scalapack=1 --download-ptscotch=1 --download-mumps=1 --download-hypre=1 --download-umfpack=1
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: MatMult() line 2235 in /data/stephan/git/petsc/src/mat/interface/matrix.c
[0]PETSC ERROR: MatMult_SchurComplement() line 82 in /data/stephan/git/petsc/src/ksp/ksp/utils/schurm.c
[0]PETSC ERROR: MatMult() line 2244 in /data/stephan/git/petsc/src/mat/interface/matrix.c
[0]PETSC ERROR: KSP_MatMult() line 204 in /data/stephan/git/petsc/include/petsc-private/kspimpl.h
[0]PETSC ERROR: KSPSolve_CG() line 218 in /data/stephan/git/petsc/src/ksp/ksp/impls/cg/cg.c
[0]PETSC ERROR: KSPSolve() line 458 in /data/stephan/git/petsc/src/ksp/ksp/interface/itfunc.c

This seems to be caused by the fact that the global sizes of the tranpose A01 matrix have been set with PETSC_DECIDE (-1) and they never get filled in. I don't really know how this worked before. I 
tried adding a 'call MatSetup(A01, ierr)' between the MatCreateTranspose call and the MatCreateSchurComplement call, but that didn't help.

...

Oh actually, digging a little deeper this may be caused by the following change:

$ git diff -r 33d57670fcdbf57d9203d482728f549b81403a0e~ src/mat/impls/transpose/transm.c
diff --git a/src/mat/impls/transpose/transm.c b/src/mat/impls/transpose/transm.c
index b935782..91cd6ee 100644
--- a/src/mat/impls/transpose/transm.c
+++ b/src/mat/impls/transpose/transm.c
@@ -112,10 +112,7 @@ PetscErrorCode  MatCreateTranspose(Mat A,Mat *N)
    (*N)->ops->multtransposeadd = MatMultTransposeAdd_Transpose;
    (*N)->assembled             = PETSC_TRUE;

-  ierr = PetscLayoutSetBlockSize((*N)->rmap,A->cmap->bs);CHKERRQ(ierr);
-  ierr = PetscLayoutSetBlockSize((*N)->cmap,A->rmap->bs);CHKERRQ(ierr);
-  ierr = PetscLayoutSetUp((*N)->rmap);CHKERRQ(ierr);
-  ierr = PetscLayoutSetUp((*N)->cmap);CHKERRQ(ierr);
+  ierr = MatSetBlockSizes(*N,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs));CHKERRQ(ierr);
    ierr = MatSetUp(*N);CHKERRQ(ierr);
    PetscFunctionReturn(0);
  }

So previously the global sizes were set in these PetscLayoutSetup() calls.

Cheers
Stephan



More information about the petsc-dev mailing list