[petsc-users] snessetjacobian question

Jed Brown jedbrown at mcs.anl.gov
Thu Sep 22 10:36:31 CDT 2011


On Thu, Sep 22, 2011 at 17:21, Milan Mitrovic <milan.v.mitrovic at gmail.com>wrote:

> > Don't destroy the matrix, just call the preallocation routines when the
> > nonzero pattern changes.
> > You can actually destroy the matrices from C, but I don't think there is
> a
> > pure Fortran way to do it (it needs another pointer indirection) and you
> > should never need to anyway.
>
> How should I create a matrix without supplying the preallocation? just
> set everything to zero and then call the prealloc routines later?
>

Just create the matrix and set the sizes, but don't preallocate. The wrapper
you are calling is just a convenience function:

PetscErrorCode  MatCreateMPIAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt
M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const
PetscInt o_nnz[],Mat *A)
{
  PetscErrorCode ierr;
  PetscMPIInt    size;

  PetscFunctionBegin;
  ierr = MatCreate(comm,A);CHKERRQ(ierr);
  ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
  ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
  if (size > 1) {
    ierr = MatSetType(*A,MATMPIAIJ);CHKERRQ(ierr);
    ierr =
MatMPIAIJSetPreallocation(*A,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr);
  } else {
    ierr = MatSetType(*A,MATSEQAIJ);CHKERRQ(ierr);
    ierr = MatSeqAIJSetPreallocation(*A,d_nz,d_nnz);CHKERRQ(ierr);
  }
  PetscFunctionReturn(0);
}




>
> > What do you mean "complains"? You can certainly use a different
> > preconditioning matrix.
>
> I dont want to use a different matrix... I want to use the same... I
> just dont want to get the error that I pasted in the first message...
> and I have no idea why I get it in the first place...
>

A hanging reference was created so the same object was destroyed multiple
times. I'm kinda surprised there wasn't a check earlier to verify that a
valid matrix was present. In any case, you shouldn't create a new matrix
inside the Jacobian evaluation function, at least not when calling from
Fortran.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20110922/29682bc2/attachment.htm>


More information about the petsc-users mailing list