xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 8d1b7334283fc0fb9ef522d6119108a52625ee8c)
18af6ec1cSBarry Smith #include <petsc/private/f90impl.h>
2a6e9e4f7SMatthew G. Knepley #include <petsc/private/sfimpl.h>
3a6e9e4f7SMatthew G. Knepley 
4a6e9e4f7SMatthew G. Knepley #if defined(PETSC_HAVE_FORTRAN_CAPS)
55c87a30dSBarry Smith #define petscsfview_       PETSCSFVIEW
6*8d1b7334SBarry Smith #define petscsfgetgraph_   PETSCSFGETGRAPH
78af6ec1cSBarry Smith #define petscsfbcastbegin_ PETSCSFBCASTBEGIN
88af6ec1cSBarry Smith #define petscsfbcastend_   PETSCSFBCASTEND
9a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
105c87a30dSBarry Smith #define petscsfgetgraph_   petscsfgetgraph
115c87a30dSBarry Smith #define petscsfview_       petscsfview
128af6ec1cSBarry Smith #define petscsfbcastbegin_ petscsfbcastbegin
138af6ec1cSBarry Smith #define petscsfbcastend_   petscsfbcastend
14a6e9e4f7SMatthew G. Knepley #endif
15a6e9e4f7SMatthew G. Knepley 
161943db53SBarry Smith PETSC_EXTERN void PETSC_STDCALL petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
17a6e9e4f7SMatthew G. Knepley {
18a6e9e4f7SMatthew G. Knepley   PetscViewer v;
19a6e9e4f7SMatthew G. Knepley 
20a6e9e4f7SMatthew G. Knepley   PetscPatchDefaultViewers_Fortran(vin, v);
21a6e9e4f7SMatthew G. Knepley   *ierr = PetscSFView(*sf, v);
22a6e9e4f7SMatthew G. Knepley }
238af6ec1cSBarry Smith 
245c87a30dSBarry Smith 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))
255c87a30dSBarry Smith {
265c87a30dSBarry Smith   const PetscInt    *ilocal;
275c87a30dSBarry Smith   const PetscSFNode *iremote;
285c87a30dSBarry Smith 
295c87a30dSBarry Smith   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
306a43dea8SBarry Smith   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal));
315c87a30dSBarry Smith   /* 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 */
325c87a30dSBarry Smith   *ierr = F90Array1dCreate((void*)iremote,MPIU_INT,1,2*(*nleaves), airemote PETSC_F90_2PTR_PARAM(pilocal));
335c87a30dSBarry Smith }
345c87a30dSBarry Smith 
358af6ec1cSBarry Smith 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))
368af6ec1cSBarry Smith {
378af6ec1cSBarry Smith   MPI_Datatype dtype;
388af6ec1cSBarry Smith   const void   *rootdata;
398af6ec1cSBarry Smith   void         *leafdata;
408af6ec1cSBarry Smith 
415c87a30dSBarry Smith 
425c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
438af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
448af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
458af6ec1cSBarry Smith   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata);
468af6ec1cSBarry Smith }
478af6ec1cSBarry Smith 
488af6ec1cSBarry Smith 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))
498af6ec1cSBarry Smith {
508af6ec1cSBarry Smith   MPI_Datatype dtype;
518af6ec1cSBarry Smith   const void   *rootdata;
528af6ec1cSBarry Smith   void         *leafdata;
538af6ec1cSBarry Smith 
545c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
558af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
568af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
578af6ec1cSBarry Smith   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata);
588af6ec1cSBarry Smith }
59