xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 94a885e832907c4332a272c2ada6a4d063ea1de0)
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
128dbb0df6SBarry Smith #define petscsfsetgraph_        PETSCSFSETGRAPH
13*94a885e8SJunchao Zhang #define petscsfgetleafranks_    PETSCSFGETLEAFRANKS
14*94a885e8SJunchao Zhang #define petscsfgetrootranks_    PETSCSFGETROOTRANKS
15a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
165c87a30dSBarry Smith #define petscsfgetgraph_        petscsfgetgraph
175c87a30dSBarry Smith #define petscsfview_            petscsfview
188af6ec1cSBarry Smith #define petscsfbcastbegin_      petscsfbcastbegin
198af6ec1cSBarry Smith #define petscsfbcastend_        petscsfbcastend
207f139299SBarry Smith #define f90arraysfnodecreate_   f90arraysfnodecreate
21fe2efc57SMark #define petscsfviewfromoptions_ petscsfviewfromoptions
221fb7b255SJunchao Zhang #define petscsfdestroy_         petscsfdestroy
238dbb0df6SBarry Smith #define petscsfsetgraph_        petscsfsetgraph
24*94a885e8SJunchao Zhang #define petscsfgetleafranks_    petscsfgetleafranks
25*94a885e8SJunchao Zhang #define petscsfgetrootranks_    petscsfgetrootranks
26a6e9e4f7SMatthew G. Knepley #endif
27a6e9e4f7SMatthew G. Knepley 
2819caf8f3SSatish Balay PETSC_EXTERN void f90arraysfnodecreate_(const PetscInt *,PetscInt *,void * PETSC_F90_2PTR_PROTO_NOVAR);
297f139299SBarry Smith 
308dbb0df6SBarry Smith PETSC_EXTERN void  petscsfsetgraph_(PetscSF *sf,PetscInt *nroots,PetscInt *nleaves, PetscInt *ilocal,PetscCopyMode *localmode, PetscSFNode *iremote,PetscCopyMode *remotemode, int *ierr)
318dbb0df6SBarry Smith {
328dbb0df6SBarry Smith   if (ilocal == PETSC_NULL_INTEGER_Fortran) ilocal = NULL;
338dbb0df6SBarry Smith   *ierr = PetscSFSetGraph(*sf,*nroots,*nleaves,ilocal,*localmode,iremote,*remotemode);
348dbb0df6SBarry Smith }
358dbb0df6SBarry Smith 
3619caf8f3SSatish Balay PETSC_EXTERN void petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
37a6e9e4f7SMatthew G. Knepley {
38a6e9e4f7SMatthew G. Knepley   PetscViewer v;
39a6e9e4f7SMatthew G. Knepley 
40a6e9e4f7SMatthew G. Knepley   PetscPatchDefaultViewers_Fortran(vin, v);
41a6e9e4f7SMatthew G. Knepley   *ierr = PetscSFView(*sf, v);
42a6e9e4f7SMatthew G. Knepley }
438af6ec1cSBarry Smith 
4461bf59e3SJunchao 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))
455c87a30dSBarry Smith {
465c87a30dSBarry Smith   const PetscInt    *ilocal;
475c87a30dSBarry Smith   const PetscSFNode *iremote;
488dbb0df6SBarry Smith   PetscInt          nl;
495c87a30dSBarry Smith 
505c87a30dSBarry Smith   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
518dbb0df6SBarry Smith   nl = *nleaves;
528dbb0df6SBarry Smith   if (!ilocal) nl = 0;
538dbb0df6SBarry Smith   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,nl, ailocal PETSC_F90_2PTR_PARAM(pilocal));
547f139299SBarry Smith   /* this creates a memory leak */
558e383715SSatish Balay   f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
565c87a30dSBarry Smith }
575c87a30dSBarry Smith 
58*94a885e8SJunchao Zhang PETSC_EXTERN void petscsfgetleafranks_(PetscSF *sf, PetscInt *niranks, F90Array1d *airanks, F90Array1d *aioffset, F90Array1d *airootloc, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(piranks) PETSC_F90_2PTR_PROTO(pioffset) PETSC_F90_2PTR_PROTO(pirootloc))
59*94a885e8SJunchao Zhang {
60*94a885e8SJunchao Zhang   const PetscMPIInt *iranks = NULL;
61*94a885e8SJunchao Zhang   const PetscInt    *ioffset = NULL;
62*94a885e8SJunchao Zhang   const PetscInt    *irootloc = NULL;
63*94a885e8SJunchao Zhang 
64*94a885e8SJunchao Zhang   *ierr = PetscSFGetLeafRanks(*sf, niranks, &iranks, &ioffset, &irootloc);if (*ierr) return;
65*94a885e8SJunchao Zhang   *ierr = F90Array1dCreate((void*)irootloc, MPIU_INT, 1, ioffset[*niranks], airootloc PETSC_F90_2PTR_PARAM(pirootloc));if (*ierr) return;
66*94a885e8SJunchao Zhang   *ierr = F90Array1dCreate((void*)iranks, MPI_INT, 1, *niranks, airanks PETSC_F90_2PTR_PARAM(piranks));if (*ierr) return;
67*94a885e8SJunchao Zhang   *ierr = F90Array1dCreate((void*)ioffset, MPIU_INT, 1, *niranks+1, aioffset PETSC_F90_2PTR_PARAM(pioffset));if (*ierr) return;
68*94a885e8SJunchao Zhang }
69*94a885e8SJunchao Zhang 
70*94a885e8SJunchao Zhang PETSC_EXTERN void petscsfgetrootranks_(PetscSF *sf, PetscInt *nranks, F90Array1d *aranks, F90Array1d *aroffset, F90Array1d *armine, F90Array1d *arremote,
71*94a885e8SJunchao Zhang  PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(pranks) PETSC_F90_2PTR_PROTO(proffset) PETSC_F90_2PTR_PROTO(prmine) PETSC_F90_2PTR_PROTO(prremote))
72*94a885e8SJunchao Zhang {
73*94a885e8SJunchao Zhang   const PetscMPIInt *ranks = NULL;
74*94a885e8SJunchao Zhang   const PetscInt    *roffset = NULL;
75*94a885e8SJunchao Zhang   const PetscInt    *rmine = NULL;
76*94a885e8SJunchao Zhang   const PetscInt    *rremote = NULL;
77*94a885e8SJunchao Zhang 
78*94a885e8SJunchao Zhang   *ierr = PetscSFGetRootRanks(*sf, nranks, &ranks, &roffset, &rmine, &rremote);if (*ierr) return;
79*94a885e8SJunchao Zhang   *ierr = F90Array1dCreate((void*)ranks, MPI_INT, 1, *nranks, aranks PETSC_F90_2PTR_PARAM(pranks));if (*ierr) return;
80*94a885e8SJunchao Zhang   *ierr = F90Array1dCreate((void*)roffset, MPIU_INT, 1, *nranks+1, aroffset PETSC_F90_2PTR_PARAM(proffset));if (*ierr) return;
81*94a885e8SJunchao Zhang   *ierr = F90Array1dCreate((void*)rmine, MPIU_INT, 1, roffset[*nranks], armine PETSC_F90_2PTR_PARAM(prmine));if (*ierr) return;
82*94a885e8SJunchao Zhang   *ierr = F90Array1dCreate((void*)rremote, MPIU_INT, 1, roffset[*nranks], arremote PETSC_F90_2PTR_PARAM(prremote));if (*ierr) return;
83*94a885e8SJunchao Zhang }
84*94a885e8SJunchao Zhang 
856f7e44deSSatish Balay #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
86ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
876f7e44deSSatish Balay {
886f7e44deSSatish Balay   MPI_Datatype dtype;
89ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
906f7e44deSSatish Balay 
916f7e44deSSatish Balay   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
92ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr, cop);
936f7e44deSSatish Balay }
946f7e44deSSatish Balay 
95ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
966f7e44deSSatish Balay {
976f7e44deSSatish Balay   MPI_Datatype dtype;
98ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
996f7e44deSSatish Balay 
1006f7e44deSSatish Balay   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
101ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr, cop);
1026f7e44deSSatish Balay }
1036f7e44deSSatish Balay 
1046f7e44deSSatish Balay #else
1056f7e44deSSatish Balay 
106ad227feaSJunchao 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))
1078af6ec1cSBarry Smith {
1088af6ec1cSBarry Smith   MPI_Datatype dtype;
1098af6ec1cSBarry Smith   const void   *rootdata;
1108af6ec1cSBarry Smith   void         *leafdata;
111ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1128af6ec1cSBarry Smith 
1135c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
1148af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
1158af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
116ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata, cop);
1178af6ec1cSBarry Smith }
1188af6ec1cSBarry Smith 
119ad227feaSJunchao 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))
1208af6ec1cSBarry Smith {
1218af6ec1cSBarry Smith   MPI_Datatype dtype;
1228af6ec1cSBarry Smith   const void   *rootdata;
1238af6ec1cSBarry Smith   void         *leafdata;
124ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1258af6ec1cSBarry Smith 
1265c87a30dSBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
1278af6ec1cSBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
1288af6ec1cSBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
129ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata, cop);
1308af6ec1cSBarry Smith }
13119caf8f3SSatish Balay PETSC_EXTERN void petscsfviewfromoptions_(PetscSF *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
132fe2efc57SMark {
133fe2efc57SMark   char *t;
134fe2efc57SMark 
135fe2efc57SMark   FIXCHAR(type,len,t);
136b14c0cbaSBlaise Bourdin   CHKFORTRANNULLOBJECT(obj);
137fe2efc57SMark   *ierr = PetscSFViewFromOptions(*ao,obj,t);if (*ierr) return;
138fe2efc57SMark   FREECHAR(type,t);
139fe2efc57SMark }
1406f7e44deSSatish Balay 
1411fb7b255SJunchao Zhang PETSC_EXTERN void petscsfdestroy_(PetscSF *x,int *ierr)
1421fb7b255SJunchao Zhang {
1431fb7b255SJunchao Zhang   PETSC_FORTRAN_OBJECT_F_DESTROYED_TO_C_NULL(x);
1441fb7b255SJunchao Zhang   *ierr = PetscSFDestroy(x); if (*ierr) return;
1451fb7b255SJunchao Zhang   PETSC_FORTRAN_OBJECT_C_NULL_TO_F_DESTROYED(x);
1461fb7b255SJunchao Zhang }
1471fb7b255SJunchao Zhang 
1486f7e44deSSatish Balay #endif
149