xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 6f7e44deaea2893df81d7a817cb6f425adb425ac)
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
68d1b7334SBarry Smith #define petscsfgetgraph_      PETSCSFGETGRAPH
78af6ec1cSBarry Smith #define petscsfbcastbegin_    PETSCSFBCASTBEGIN
88af6ec1cSBarry Smith #define petscsfbcastend_      PETSCSFBCASTEND
97f139299SBarry Smith #define f90arraysfnodecreate_ F90ARRAYSFNODECREATE
10a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
115c87a30dSBarry Smith #define petscsfgetgraph_      petscsfgetgraph
125c87a30dSBarry Smith #define petscsfview_          petscsfview
138af6ec1cSBarry Smith #define petscsfbcastbegin_    petscsfbcastbegin
148af6ec1cSBarry Smith #define petscsfbcastend_      petscsfbcastend
157f139299SBarry Smith #define f90arraysfnodecreate_ f90arraysfnodecreate
16a6e9e4f7SMatthew G. Knepley #endif
17a6e9e4f7SMatthew G. Knepley 
187f139299SBarry Smith PETSC_EXTERN void PETSC_STDCALL f90arraysfnodecreate_(const PetscInt *,PetscInt *,void * PETSC_F90_2PTR_PROTO_NOVAR);
197f139299SBarry Smith 
201943db53SBarry Smith PETSC_EXTERN void PETSC_STDCALL petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
21a6e9e4f7SMatthew G. Knepley {
22a6e9e4f7SMatthew G. Knepley   PetscViewer v;
23a6e9e4f7SMatthew G. Knepley 
24a6e9e4f7SMatthew G. Knepley   PetscPatchDefaultViewers_Fortran(vin, v);
25a6e9e4f7SMatthew G. Knepley   *ierr = PetscSFView(*sf, v);
26a6e9e4f7SMatthew G. Knepley }
278af6ec1cSBarry Smith 
287f139299SBarry Smith 
295c87a30dSBarry 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))
305c87a30dSBarry Smith {
315c87a30dSBarry Smith   const PetscInt    *ilocal;
325c87a30dSBarry Smith   const PetscSFNode *iremote;
335c87a30dSBarry Smith 
345c87a30dSBarry Smith   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
356a43dea8SBarry Smith   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal));
367f139299SBarry Smith   /* this creates a memory leak */
378e383715SSatish Balay   f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
385c87a30dSBarry Smith }
395c87a30dSBarry Smith 
40*6f7e44deSSatish Balay #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
41*6f7e44deSSatish Balay PETSC_EXTERN void PETSC_STDCALL petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, int *ierr)
42*6f7e44deSSatish Balay {
43*6f7e44deSSatish Balay   MPI_Datatype dtype;
44*6f7e44deSSatish Balay 
45*6f7e44deSSatish Balay   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
46*6f7e44deSSatish Balay   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr);
47*6f7e44deSSatish Balay }
48*6f7e44deSSatish Balay 
49*6f7e44deSSatish Balay 
50*6f7e44deSSatish Balay PETSC_EXTERN void PETSC_STDCALL petscsfbcastend_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, int *ierr)
51*6f7e44deSSatish Balay {
52*6f7e44deSSatish Balay   MPI_Datatype dtype;
53*6f7e44deSSatish Balay 
54*6f7e44deSSatish Balay   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
55*6f7e44deSSatish Balay   *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr);
56*6f7e44deSSatish Balay }
57*6f7e44deSSatish Balay 
58*6f7e44deSSatish Balay #else
59*6f7e44deSSatish Balay 
608af6ec1cSBarry 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))
618af6ec1cSBarry Smith {
628af6ec1cSBarry Smith   MPI_Datatype dtype;
638af6ec1cSBarry Smith   const void   *rootdata;
648af6ec1cSBarry Smith   void         *leafdata;
658af6ec1cSBarry Smith 
665c87a30dSBarry Smith 
675c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
688af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
698af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
708af6ec1cSBarry Smith   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata);
718af6ec1cSBarry Smith }
728af6ec1cSBarry Smith 
738af6ec1cSBarry 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))
748af6ec1cSBarry Smith {
758af6ec1cSBarry Smith   MPI_Datatype dtype;
768af6ec1cSBarry Smith   const void   *rootdata;
778af6ec1cSBarry Smith   void         *leafdata;
788af6ec1cSBarry Smith 
795c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
808af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
818af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
828af6ec1cSBarry Smith   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata);
838af6ec1cSBarry Smith }
84*6f7e44deSSatish Balay 
85*6f7e44deSSatish Balay #endif
86