xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 8af6ec1cc8e49a30993fa2a7b177d7989c62ae61)
1*8af6ec1cSBarry 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)
51943db53SBarry Smith #define sfview_           PETSCSFVIEW
6*8af6ec1cSBarry Smith #define petscsfbcastbegin_ PETSCSFBCASTBEGIN
7*8af6ec1cSBarry Smith #define petscsfbcastend_   PETSCSFBCASTEND
8a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
91943db53SBarry Smith #define sfview_            petscsfview
10*8af6ec1cSBarry Smith #define petscsfbcastbegin_ petscsfbcastbegin
11*8af6ec1cSBarry Smith #define petscsfbcastend_   petscsfbcastend
12a6e9e4f7SMatthew G. Knepley #endif
13a6e9e4f7SMatthew G. Knepley 
141943db53SBarry Smith PETSC_EXTERN void PETSC_STDCALL petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
15a6e9e4f7SMatthew G. Knepley {
16a6e9e4f7SMatthew G. Knepley   PetscViewer v;
17a6e9e4f7SMatthew G. Knepley 
18a6e9e4f7SMatthew G. Knepley   PetscPatchDefaultViewers_Fortran(vin, v);
19a6e9e4f7SMatthew G. Knepley   *ierr = PetscSFView(*sf, v);
20a6e9e4f7SMatthew G. Knepley }
21*8af6ec1cSBarry Smith 
22*8af6ec1cSBarry 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))
23*8af6ec1cSBarry Smith {
24*8af6ec1cSBarry Smith   MPI_Datatype dtype;
25*8af6ec1cSBarry Smith   const void   *rootdata;
26*8af6ec1cSBarry Smith   void         *leafdata;
27*8af6ec1cSBarry Smith 
28*8af6ec1cSBarry Smith   dtype = MPI_Type_f2c(*unit);
29*8af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
30*8af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
31*8af6ec1cSBarry Smith   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata);
32*8af6ec1cSBarry Smith }
33*8af6ec1cSBarry Smith 
34*8af6ec1cSBarry 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))
35*8af6ec1cSBarry Smith {
36*8af6ec1cSBarry Smith   MPI_Datatype dtype;
37*8af6ec1cSBarry Smith   const void   *rootdata;
38*8af6ec1cSBarry Smith   void         *leafdata;
39*8af6ec1cSBarry Smith 
40*8af6ec1cSBarry Smith   dtype = MPI_Type_f2c(*unit);
41*8af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
42*8af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
43*8af6ec1cSBarry Smith   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata);
44*8af6ec1cSBarry Smith }
45