[petsc-users] How to check NULL pointer in Fortran

Barry Smith bsmith at petsc.dev
Fri Nov 19 07:20:39 CST 2021


  it is very possible there is a bug in the Fortran interface for this function. It looks like whoever wrote the Fortran interface did not think through the special case of contiguous entries.

PETSC_EXTERN void  petscsfgetgraph_(PetscSF *sf,PetscInt *nroots,PetscInt *nleaves, F90Array1d  *ailocal, F90Array1d  *airemote, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(pilocal) PETSC_F90_2PTR_PROTO(piremote))
{
  const PetscInt    *ilocal;
  const PetscSFNode *iremote;

  *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
  *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal));
  /* this creates a memory leak */
  f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
}

PetscErrorCode F90Array1dCreate(void *array,MPI_Datatype type,PetscInt start,PetscInt len,F90Array1d *ptr PETSC_F90_2PTR_PROTO(ptrd))
{
  PetscFunctionBegin;
  if (type == MPIU_SCALAR) {
    if (!len) array = PETSC_NULL_SCALAR_Fortran;
    f90array1dcreatescalar_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
  } else if (type == MPIU_REAL) {
    if (!len) array = PETSC_NULL_REAL_Fortran;
    f90array1dcreatereal_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
  } else if (type == MPIU_INT) {
    if (!len) array = PETSC_NULL_INTEGER_Fortran;
    f90array1dcreateint_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
  } else if (type == MPIU_FORTRANADDR) {
    f90array1dcreatefortranaddr_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
  } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
  PetscFunctionReturn(0);
}

It blindly calls F90Array1dCreate() with a len = *nleaves > 0 but a NULL array ilocal then f90array1dcreateint_ must create an invalid object.

Barry



> On Nov 19, 2021, at 12:34 AM, 袁煕 <yuanxi at advancesoft.jp> wrote:
> 
> Dear PETSc-team,
> 
> I am using function PetscSFGetGraph in my program like
> -----------------------
> call PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr)
> -----------------------
> 
> In some cases, it works well. But in some cases, I encountered following error
> ------------------
> PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range 
> ------------------
> I found it is due to the action of reading gmine, which is a fortran pointer point to an array. It is reasonable because PETSc manual tells me "if returned value is NULL, it means leaves are in contiguous storage". The problem is that I cannot find if gmine is a null pointer. 
> 1)  I cannot use "if (gmine==PETSC_NULL_INTEGER)" because my intel compiler would return a compile error "A scalar-valued expression is required in this context"
> 2)  When using standard style of checking a null pointer in fortran, "associated(gmine)", it returns "T". Even in cases such action of "print *, gmine(1)" would give rise to above Segmentation Violation error.
> Is there any means to check the NULL pointer in Fortran in above cases.
> Many thanks,
> Yuan



More information about the petsc-users mailing list