xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision fe2efc57b7777594bce7568e90822861480cbdc8)
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
10*fe2efc57SMark #define petscsfviewfromoptions_ PETSCSFVIEWFROMOPTIONS
11a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
125c87a30dSBarry Smith #define petscsfgetgraph_      petscsfgetgraph
135c87a30dSBarry Smith #define petscsfview_          petscsfview
148af6ec1cSBarry Smith #define petscsfbcastbegin_    petscsfbcastbegin
158af6ec1cSBarry Smith #define petscsfbcastend_      petscsfbcastend
167f139299SBarry Smith #define f90arraysfnodecreate_ f90arraysfnodecreate
17*fe2efc57SMark #define petscsfviewfromoptions_ petscsfviewfromoptions
18a6e9e4f7SMatthew G. Knepley #endif
19a6e9e4f7SMatthew G. Knepley 
207f139299SBarry Smith PETSC_EXTERN void PETSC_STDCALL f90arraysfnodecreate_(const PetscInt *,PetscInt *,void * PETSC_F90_2PTR_PROTO_NOVAR);
217f139299SBarry Smith 
221943db53SBarry Smith PETSC_EXTERN void PETSC_STDCALL petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
23a6e9e4f7SMatthew G. Knepley {
24a6e9e4f7SMatthew G. Knepley   PetscViewer v;
25a6e9e4f7SMatthew G. Knepley 
26a6e9e4f7SMatthew G. Knepley   PetscPatchDefaultViewers_Fortran(vin, v);
27a6e9e4f7SMatthew G. Knepley   *ierr = PetscSFView(*sf, v);
28a6e9e4f7SMatthew G. Knepley }
298af6ec1cSBarry Smith 
307f139299SBarry Smith 
315c87a30dSBarry 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))
325c87a30dSBarry Smith {
335c87a30dSBarry Smith   const PetscInt    *ilocal;
345c87a30dSBarry Smith   const PetscSFNode *iremote;
355c87a30dSBarry Smith 
365c87a30dSBarry Smith   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
376a43dea8SBarry Smith   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal));
387f139299SBarry Smith   /* this creates a memory leak */
398e383715SSatish Balay   f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
405c87a30dSBarry Smith }
415c87a30dSBarry Smith 
428af6ec1cSBarry 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))
438af6ec1cSBarry Smith {
448af6ec1cSBarry Smith   MPI_Datatype dtype;
458af6ec1cSBarry Smith   const void   *rootdata;
468af6ec1cSBarry Smith   void         *leafdata;
478af6ec1cSBarry Smith 
485c87a30dSBarry Smith 
495c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
508af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
518af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
528af6ec1cSBarry Smith   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata);
538af6ec1cSBarry Smith }
548af6ec1cSBarry Smith 
558af6ec1cSBarry 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))
568af6ec1cSBarry Smith {
578af6ec1cSBarry Smith   MPI_Datatype dtype;
588af6ec1cSBarry Smith   const void   *rootdata;
598af6ec1cSBarry Smith   void         *leafdata;
608af6ec1cSBarry Smith 
615c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
628af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
638af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
648af6ec1cSBarry Smith   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata);
658af6ec1cSBarry Smith }
66*fe2efc57SMark PETSC_EXTERN void PETSC_STDCALL petscsfviewfromoptions_(PetscSF *ao,PetscObject obj,char* type PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
67*fe2efc57SMark {
68*fe2efc57SMark   char *t;
69*fe2efc57SMark 
70*fe2efc57SMark   FIXCHAR(type,len,t);
71*fe2efc57SMark   *ierr = PetscSFViewFromOptions(*ao,obj,t);if (*ierr) return;
72*fe2efc57SMark   FREECHAR(type,t);
73*fe2efc57SMark }
74