hi Barry,

On 23/10/17 15:04, Barry Smith wrote:
   Adrian,

    Sorry for the delay in generating the bindings for you. I have put them in 
the branch

barry/fortran-PetscSFBcastBegin

the example src/vec/sf/examples/tutorials/ex1f.F90 attempts to test them. 
Please let us know if this works for you and if you need anything else. They 
will only work for MPIU_INT and MPIU_REAL and MPIU_SCALAR because that is the 
only Fortran pointer arrays we support currently.

Oh, thanks very much for that. I really only need support for integers in my application at present, so what you've done should be fine for now.

I checked out your new branch but it is giving a bunch of warnings and errors when I try to make it.

The warnings are:

In file included from /home/acro018/software/PETSc/code/include/petscsys.h:1631:0, from /home/acro018/software/PETSc/code/include/petsc/private/petscimpl.h:8, from /home/acro018/software/PETSc/code/include/petsc/private/fortranimpl.h:6, from /home/acro018/software/PETSc/code/include/petsc/private/f90impl.h:4, from /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:1: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array1dCreate’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:60:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array1dAccess’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:78:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array1dDestroy’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:93:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array2dCreate’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:151:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array2dAccess’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:166:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array2dDestroy’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:181:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array3dCreate’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:239:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array3dAccess’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:254:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array3dDestroy’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:269:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array4dCreate’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:320:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array4dAccess’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:335:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)
^
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c: In function ‘F90Array4dDestroy’: /home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c:344:80: warning: cast from pointer to integer of different size [-Wpointer-to-int-cast] } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported MPI_Datatype: %d",(PetscInt)type);
^
/home/acro018/software/PETSc/code/include/petscerror.h:126:122: note: in definition of macro ‘SETERRQ1’ #define SETERRQ1(comm,ierr,s,a1) return PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,ierr,PETSC_ERROR_INITIAL,s,a1)



The errors look like they're from some unrelated DMPlex box mesh changes:

/home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c: In function ‘dmplexcreateboxmesh_’: /home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c:105:11: error: incompatible type for argument 4 of ‘DMPlexCreateBoxMesh’
 *__ierr = DMPlexCreateBoxMesh(
           ^
In file included from /home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c:30:0: /home/acro018/software/PETSc/code/include/petscdmplex.h:124:29: note: expected ‘const PetscInt *’ but argument is of type ‘PetscBool’ PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
                             ^
/home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c:106:52: warning: passing argument 5 of ‘DMPlexCreateBoxMesh’ from incompatible pointer type
  MPI_Comm_f2c(*(comm)),*dim,*numFaces,*interpolate,dm);
                                                    ^
In file included from /home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c:30:0: /home/acro018/software/PETSc/code/include/petscdmplex.h:124:29: note: expected ‘const PetscReal *’ but argument is of type ‘struct _p_DM **’ PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
                             ^
/home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c:105:11: error: too few arguments to function ‘DMPlexCreateBoxMesh’
 *__ierr = DMPlexCreateBoxMesh(
           ^
In file included from /home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c:30:0: /home/acro018/software/PETSc/code/include/petscdmplex.h:124:29: note: declared here PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
                             ^
/home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c: In function ‘dmplexcreatehexboxmesh_’: /home/acro018/software/PETSc/code/src/dm/impls/plex/ftn-auto/plexcreatef.c:109:1: warning: implicit declaration of function ‘DMPlexCreateHexBoxMesh’ [-Wimplicit-function-declaration]
 *__ierr = DMPlexCreateHexBoxMesh(
 ^
gmakefile:227: recipe for target 'linux-gnu-c-opt/obj/src/dm/impls/plex/ftn-auto/plexcreatef.o' failed make[2]: *** [linux-gnu-c-opt/obj/src/dm/impls/plex/ftn-auto/plexcreatef.o] Error 1


I was shocked, shocked to discover that it looks like PetscSF has essentially 
no bindings for Fortran. We can add more relatively easily if you need them. I 
don't understand how you could previously test your code given that 
PetscSFCreate() didn't even have a Fortran stub generated?

My tests are based on DMPlex and I was just using DMGetPointSF() and DMPlexDistribute() to get some SFs to test.


I have also written some Fortran bindings for PetscSFGetGraph() and PetscSFSetGraph(). However it looks like you have made the PetscSFNode type accessible from Fortran in your branch? In which case my bindings could probably be modified to make use of that. Currently they look like this (using a 2-by-n array of int instead of an array of PetscSFNode for the root data array):

+      Interface
+         Subroutine PetscSFGetGraph(sf,nroots,nleaves,
+     &      larray,rarray,ierr)
+          use petscisdef
+          PetscSF :: sf
+          PetscInt :: nroots, nleaves
+          PetscInt, pointer :: larray(:)
+          PetscInt, pointer :: rarray(:,:)
+          PetscErrorCode :: ierr
+        End Subroutine
+      End Interface


+PETSC_EXTERN void PETSC_STDCALL petscsfgetgraph_(PetscSF *sf, PetscInt *nroots, PetscInt *nleaves, F90Array1d *lptr, F90Array2d *rptr, int *ierr PETSC_F90_2PTR_PROTO(lptrd) PETSC_F90_2PTR_PROTO(rptrd))
+{
+  const PetscInt *ilocal;
+  const PetscSFNode *iremote;
+
+ *ierr = PetscSFGetGraph(*sf, nroots, nleaves, &ilocal, &iremote); if (*ierr) return;
+
+ *ierr = F90Array1dCreate((void*) ilocal, PETSC_INT, 1, *nleaves, lptr PETSC_F90_2PTR_PARAM(lptrd)); if (*ierr) return; + *ierr = F90Array2dCreate((void*) iremote, PETSC_INT, 1, 2, 1, *nleaves, rptr PETSC_F90_2PTR_PARAM(rptrd));
+}

- Adrian

--
Dr Adrian Croucher
Senior Research Fellow
Department of Engineering Science
University of Auckland, New Zealand
email: [email protected]
tel: +64 (0)9 923 4611

Reply via email to