xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 8dbb0df66754dee4fb72dff2ad56e76db1e6b7c7)
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
10fe2efc57SMark #define petscsfviewfromoptions_ PETSCSFVIEWFROMOPTIONS
111fb7b255SJunchao Zhang #define petscsfdestroy_         PETSCSFDESTROY
12*8dbb0df6SBarry Smith #define petscsfsetgraph_        PETSCSFSETGRAPH
13a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
145c87a30dSBarry Smith #define petscsfgetgraph_        petscsfgetgraph
155c87a30dSBarry Smith #define petscsfview_            petscsfview
168af6ec1cSBarry Smith #define petscsfbcastbegin_      petscsfbcastbegin
178af6ec1cSBarry Smith #define petscsfbcastend_        petscsfbcastend
187f139299SBarry Smith #define f90arraysfnodecreate_   f90arraysfnodecreate
19fe2efc57SMark #define petscsfviewfromoptions_ petscsfviewfromoptions
201fb7b255SJunchao Zhang #define petscsfdestroy_         petscsfdestroy
21*8dbb0df6SBarry Smith #define petscsfsetgraph_        petscsfsetgraph
22a6e9e4f7SMatthew G. Knepley #endif
23a6e9e4f7SMatthew G. Knepley 
2419caf8f3SSatish Balay PETSC_EXTERN void f90arraysfnodecreate_(const PetscInt *,PetscInt *,void * PETSC_F90_2PTR_PROTO_NOVAR);
257f139299SBarry Smith 
26*8dbb0df6SBarry Smith PETSC_EXTERN void  petscsfsetgraph_(PetscSF *sf,PetscInt *nroots,PetscInt *nleaves, PetscInt *ilocal,PetscCopyMode *localmode, PetscSFNode *iremote,PetscCopyMode *remotemode, int *ierr)
27*8dbb0df6SBarry Smith {
28*8dbb0df6SBarry Smith   if (ilocal == PETSC_NULL_INTEGER_Fortran) ilocal = NULL;
29*8dbb0df6SBarry Smith   *ierr = PetscSFSetGraph(*sf,*nroots,*nleaves,ilocal,*localmode,iremote,*remotemode);
30*8dbb0df6SBarry Smith }
31*8dbb0df6SBarry Smith 
3219caf8f3SSatish Balay PETSC_EXTERN void petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
33a6e9e4f7SMatthew G. Knepley {
34a6e9e4f7SMatthew G. Knepley   PetscViewer v;
35a6e9e4f7SMatthew G. Knepley 
36a6e9e4f7SMatthew G. Knepley   PetscPatchDefaultViewers_Fortran(vin, v);
37a6e9e4f7SMatthew G. Knepley   *ierr = PetscSFView(*sf, v);
38a6e9e4f7SMatthew G. Knepley }
398af6ec1cSBarry Smith 
4061bf59e3SJunchao Zhang PETSC_EXTERN void  petscsfgetgraph_(PetscSF *sf,PetscInt *nroots,PetscInt *nleaves, F90Array1d  *ailocal, F90Array1d  *airemote, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(pilocal) PETSC_F90_2PTR_PROTO(piremote))
415c87a30dSBarry Smith {
425c87a30dSBarry Smith   const PetscInt    *ilocal;
435c87a30dSBarry Smith   const PetscSFNode *iremote;
44*8dbb0df6SBarry Smith   PetscInt          nl;
455c87a30dSBarry Smith 
465c87a30dSBarry Smith   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
47*8dbb0df6SBarry Smith   nl = *nleaves;
48*8dbb0df6SBarry Smith   if (!ilocal) nl = 0;
49*8dbb0df6SBarry Smith   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,nl, ailocal PETSC_F90_2PTR_PARAM(pilocal));
507f139299SBarry Smith   /* this creates a memory leak */
518e383715SSatish Balay   f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
525c87a30dSBarry Smith }
535c87a30dSBarry Smith 
546f7e44deSSatish Balay #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
55ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
566f7e44deSSatish Balay {
576f7e44deSSatish Balay   MPI_Datatype dtype;
58ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
596f7e44deSSatish Balay 
606f7e44deSSatish Balay   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
61ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr, cop);
626f7e44deSSatish Balay }
636f7e44deSSatish Balay 
64ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
656f7e44deSSatish Balay {
666f7e44deSSatish Balay   MPI_Datatype dtype;
67ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
686f7e44deSSatish Balay 
696f7e44deSSatish Balay   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
70ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr, cop);
716f7e44deSSatish Balay }
726f7e44deSSatish Balay 
736f7e44deSSatish Balay #else
746f7e44deSSatish Balay 
75ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, MPI_Fint *op, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd))
768af6ec1cSBarry Smith {
778af6ec1cSBarry Smith   MPI_Datatype dtype;
788af6ec1cSBarry Smith   const void   *rootdata;
798af6ec1cSBarry Smith   void         *leafdata;
80ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
818af6ec1cSBarry Smith 
825c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
838af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
848af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
85ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata, cop);
868af6ec1cSBarry Smith }
878af6ec1cSBarry Smith 
88ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, MPI_Fint *op, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd))
898af6ec1cSBarry Smith {
908af6ec1cSBarry Smith   MPI_Datatype dtype;
918af6ec1cSBarry Smith   const void   *rootdata;
928af6ec1cSBarry Smith   void         *leafdata;
93ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
948af6ec1cSBarry Smith 
955c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
968af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
978af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
98ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata, cop);
998af6ec1cSBarry Smith }
10019caf8f3SSatish Balay PETSC_EXTERN void petscsfviewfromoptions_(PetscSF *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
101fe2efc57SMark {
102fe2efc57SMark   char *t;
103fe2efc57SMark 
104fe2efc57SMark   FIXCHAR(type,len,t);
105b14c0cbaSBlaise Bourdin   CHKFORTRANNULLOBJECT(obj);
106fe2efc57SMark   *ierr = PetscSFViewFromOptions(*ao,obj,t);if (*ierr) return;
107fe2efc57SMark   FREECHAR(type,t);
108fe2efc57SMark }
1096f7e44deSSatish Balay 
1101fb7b255SJunchao Zhang PETSC_EXTERN void petscsfdestroy_(PetscSF *x,int *ierr)
1111fb7b255SJunchao Zhang {
1121fb7b255SJunchao Zhang   PETSC_FORTRAN_OBJECT_F_DESTROYED_TO_C_NULL(x);
1131fb7b255SJunchao Zhang   *ierr = PetscSFDestroy(x); if (*ierr) return;
1141fb7b255SJunchao Zhang   PETSC_FORTRAN_OBJECT_C_NULL_TO_F_DESTROYED(x);
1151fb7b255SJunchao Zhang }
1161fb7b255SJunchao Zhang 
1176f7e44deSSatish Balay #endif
118