xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 18ebb3a5280c02fe10d3ba9cebbb96a615bebc44)
1 #include <petsc/private/f90impl.h>
2 #include <petsc/private/sfimpl.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5 #define petscsfview_       PETSCSFVIEW
6 #define petscsfgetgraph_   PETSCSFEGTGRAPH
7 #define petscsfbcastbegin_ PETSCSFBCASTBEGIN
8 #define petscsfbcastend_   PETSCSFBCASTEND
9 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
10 #define petscsfgetgraph_   petscsfgetgraph
11 #define petscsfview_       petscsfview
12 #define petscsfbcastbegin_ petscsfbcastbegin
13 #define petscsfbcastend_   petscsfbcastend
14 #endif
15 
16 PETSC_EXTERN void PETSC_STDCALL petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
17 {
18   PetscViewer v;
19 
20   PetscPatchDefaultViewers_Fortran(vin, v);
21   *ierr = PetscSFView(*sf, v);
22 }
23 
24 PETSC_EXTERN void PETSC_STDCALL  petscsfgetgraph_(PetscSF *sf,PetscInt *nroots,PetscInt *nleaves, F90Array1d  *ailocal, F90Array1d  *airemote, int *ierr PETSC_F90_2PTR_PROTO(pilocal) PETSC_F90_2PTR_PROTO(piremote))
25 {
26   const PetscInt    *ilocal;
27   const PetscSFNode *iremote;
28 
29   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
30   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal));
31   /* this is iffy, since airemote is actually an array of type(PetscSFNode) not an array of PetscInt; works with gfortran, needs testing with all other compilers */
32   *ierr = F90Array1dCreate((void*)iremote,MPIU_INT,1,2*(*nleaves), airemote PETSC_F90_2PTR_PARAM(pilocal));
33 }
34 
35 PETSC_EXTERN void PETSC_STDCALL petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, int *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd))
36 {
37   MPI_Datatype dtype;
38   const void   *rootdata;
39   void         *leafdata;
40 
41 
42   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
43   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
44   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
45   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata);
46 }
47 
48 PETSC_EXTERN void PETSC_STDCALL petscsfbcastend_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, int *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd))
49 {
50   MPI_Datatype dtype;
51   const void   *rootdata;
52   void         *leafdata;
53 
54   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
55   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
56   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
57   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata);
58 }
59