[petsc-users] wrong no of nonzeros
    Barry Smith 
    bsmith at mcs.anl.gov
       
    Sat Nov 29 20:53:17 CST 2014
    
    
  
  Can you please rerun with the maint branch of PETSc (this is the release version with all fixes since the release). The line numbers you report do not match up with maint. From the code: (slightly different line numbers)
  } else {
    /*
       Destroy the blocks from the previous iteration
    */
    if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
      ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
      scall = MAT_INITIAL_MATRIX;
    }
  }
  /*
     Extract out the submatrices
  */
  ierr = MatGetSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
What is suppose to have happened is that the pc->flag would have been set to DIFFERENT_NONZERO_PATTERN and hence scall reset to to MAT_INITIAIL_MATRIX (hence you should never get to the state you are in).
Now it is possible we introduced a bug but we need to debug it off of maint; not some outdated 3.5.1
  Thanks
   Barry
BTW: Ignore Jed's gibberish about PCReset(), you've got no business messing with that.
 
> On Nov 29, 2014, at 8:37 PM, Derek Gaston <friedmud at gmail.com> wrote:
> 
> I am changing the stencil of my preconditioning matrix in the middle of a solve a solve... and I'm running into this:
> 
> [0]PETSC ERROR: Nonconforming object sizes
> [0]PETSC ERROR: Cannot reuse matrix. wrong no of nonzeros
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.5.1, Jul, 24, 2014 
> [0]PETSC ERROR: ../bison-opt on a arch-darwin-c-opt named dereksmacpro.home by gastdr Sat Nov 29 21:15:28 2014
> [0]PETSC ERROR: Configure options --prefix=/opt/moose/petsc/openmpi_petsc-3.5.1/clang-opt-superlu --with-hypre-dir=/opt/moose/hypre/openmpi_hypre-2.8.0b/clang-opt --with-debugging=no --with-pic=1 --with-shared-libraries=1 --with-cc=mpicc --with-cxx=mpicxx --with-fc=mpif90 --download-fblaslapack=1 --download-metis=1 --download-parmetis=1 --download-superlu_dist=1 CC=mpicc CXX=mpicxx FC=mpif90 F77=mpif77 F90=mpif90 CFLAGS="-fPIC -fopenmp" CXXFLAGS="-fPIC -fopenmp" FFLAGS="-fPIC -fopenmp" FCFLAGS="-fPIC -fopenmp" F90FLAGS="-fPIC -fopenmp" F77FLAGS="-fPIC -fopenmp" PETSC_DIR=/opt/moose/stack_src/petsc-3.5.1
> [0]PETSC ERROR: #1 MatGetSubMatrices_MPIAIJ_Local() line 1312 in /opt/moose/stack_src/petsc-3.5.1/src/mat/impls/aij/mpi/mpiov.c
> [0]PETSC ERROR: #2 MatGetSubMatrices_MPIAIJ() line 790 in /opt/moose/stack_src/petsc-3.5.1/src/mat/impls/aij/mpi/mpiov.c
> [0]PETSC ERROR: #3 MatGetSubMatrices() line 6352 in /opt/moose/stack_src/petsc-3.5.1/src/mat/interface/matrix.c
> [0]PETSC ERROR: #4 PCSetUp_ASM() line 369 in /opt/moose/stack_src/petsc-3.5.1/src/ksp/pc/impls/asm/asm.c
> [0]PETSC ERROR: #5 PCSetUp() line 902 in /opt/moose/stack_src/petsc-3.5.1/src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #6 KSPSetUp() line 305 in /opt/moose/stack_src/petsc-3.5.1/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #7 KSPSolve() line 417 in /opt/moose/stack_src/petsc-3.5.1/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #8 SNESSolve_NEWTONLS() line 232 in /opt/moose/stack_src/petsc-3.5.1/src/snes/impls/ls/ls.c
> [0]PETSC ERROR: #9 SNESSolve() line 3743 in /opt/moose/stack_src/petsc-3.5.1/src/snes/interface/snes.c
> [0]PETSC ERROR: #10 solve() line 559 in src/solvers/petsc_nonlinear_solver.C
> 
> 
> I know that in the past we used to have to specify SAME_NONZERO_PATTERN or DIFFERENT_NONZERO_PATTERN... but I thought those times were all past us ;-)
> 
> I'm using PETSc 3.5.1 on OSX ... underneath libMesh and MOOSE of course...
> 
> Any ideas where to start digging?
> 
> Thanks!
> Derek
    
    
More information about the petsc-users
mailing list