xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 5975b3b6e3931510e2a64a701673cbe1930c6f42)
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
99037d788SNicholas Arnold-Medabalimi   #define petscsfreducebegin_     PETSCSFREDUCEBEGIN
109037d788SNicholas Arnold-Medabalimi   #define petscsfreduceend_       PETSCSFREDUCEEND
117f139299SBarry Smith   #define f90arraysfnodecreate_   F90ARRAYSFNODECREATE
12fe2efc57SMark   #define petscsfviewfromoptions_ PETSCSFVIEWFROMOPTIONS
131fb7b255SJunchao Zhang   #define petscsfdestroy_         PETSCSFDESTROY
148dbb0df6SBarry Smith   #define petscsfsetgraph_        PETSCSFSETGRAPH
1594a885e8SJunchao Zhang   #define petscsfgetleafranks_    PETSCSFGETLEAFRANKS
1694a885e8SJunchao Zhang   #define petscsfgetrootranks_    PETSCSFGETROOTRANKS
17a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
185c87a30dSBarry Smith   #define petscsfgetgraph_        petscsfgetgraph
195c87a30dSBarry Smith   #define petscsfview_            petscsfview
208af6ec1cSBarry Smith   #define petscsfbcastbegin_      petscsfbcastbegin
218af6ec1cSBarry Smith   #define petscsfbcastend_        petscsfbcastend
229037d788SNicholas Arnold-Medabalimi   #define petscsfreducebegin_     petscsfreducebegin
239037d788SNicholas Arnold-Medabalimi   #define petscsfreduceend_       petscsfreduceend
247f139299SBarry Smith   #define f90arraysfnodecreate_   f90arraysfnodecreate
25fe2efc57SMark   #define petscsfviewfromoptions_ petscsfviewfromoptions
261fb7b255SJunchao Zhang   #define petscsfdestroy_         petscsfdestroy
278dbb0df6SBarry Smith   #define petscsfsetgraph_        petscsfsetgraph
2894a885e8SJunchao Zhang   #define petscsfgetleafranks_    petscsfgetleafranks
2994a885e8SJunchao Zhang   #define petscsfgetrootranks_    petscsfgetrootranks
30a6e9e4f7SMatthew G. Knepley #endif
31a6e9e4f7SMatthew G. Knepley 
3219caf8f3SSatish Balay PETSC_EXTERN void f90arraysfnodecreate_(const PetscInt *, PetscInt *, void *PETSC_F90_2PTR_PROTO_NOVAR);
337f139299SBarry Smith 
348dbb0df6SBarry Smith PETSC_EXTERN void petscsfsetgraph_(PetscSF *sf, PetscInt *nroots, PetscInt *nleaves, PetscInt *ilocal, PetscCopyMode *localmode, PetscSFNode *iremote, PetscCopyMode *remotemode, int *ierr)
358dbb0df6SBarry Smith {
368dbb0df6SBarry Smith   if (ilocal == PETSC_NULL_INTEGER_Fortran) ilocal = NULL;
378dbb0df6SBarry Smith   *ierr = PetscSFSetGraph(*sf, *nroots, *nleaves, ilocal, *localmode, iremote, *remotemode);
388dbb0df6SBarry Smith }
398dbb0df6SBarry Smith 
4019caf8f3SSatish Balay PETSC_EXTERN void petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
41a6e9e4f7SMatthew G. Knepley {
42a6e9e4f7SMatthew G. Knepley   PetscViewer v;
43a6e9e4f7SMatthew G. Knepley 
44a6e9e4f7SMatthew G. Knepley   PetscPatchDefaultViewers_Fortran(vin, v);
45a6e9e4f7SMatthew G. Knepley   *ierr = PetscSFView(*sf, v);
46a6e9e4f7SMatthew G. Knepley }
478af6ec1cSBarry Smith 
4861bf59e3SJunchao 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))
495c87a30dSBarry Smith {
505c87a30dSBarry Smith   const PetscInt    *ilocal;
515c87a30dSBarry Smith   const PetscSFNode *iremote;
528dbb0df6SBarry Smith   PetscInt           nl;
535c87a30dSBarry Smith 
54*5975b3b6SBarry Smith   *ierr = PetscSFGetGraph(*sf, nroots, nleaves, &ilocal, &iremote);
55*5975b3b6SBarry Smith   if (*ierr) return;
568dbb0df6SBarry Smith   nl = *nleaves;
578dbb0df6SBarry Smith   if (!ilocal) nl = 0;
588dbb0df6SBarry Smith   *ierr = F90Array1dCreate((void *)ilocal, MPIU_INT, 1, nl, ailocal PETSC_F90_2PTR_PARAM(pilocal));
597f139299SBarry Smith   /* this creates a memory leak */
608e383715SSatish Balay   f90arraysfnodecreate_((PetscInt *)iremote, nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
615c87a30dSBarry Smith }
625c87a30dSBarry Smith 
6394a885e8SJunchao 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))
6494a885e8SJunchao Zhang {
6594a885e8SJunchao Zhang   const PetscMPIInt *iranks   = NULL;
6694a885e8SJunchao Zhang   const PetscInt    *ioffset  = NULL;
6794a885e8SJunchao Zhang   const PetscInt    *irootloc = NULL;
6894a885e8SJunchao Zhang 
69*5975b3b6SBarry Smith   *ierr = PetscSFGetLeafRanks(*sf, niranks, &iranks, &ioffset, &irootloc);
70*5975b3b6SBarry Smith   if (*ierr) return;
71*5975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)irootloc, MPIU_INT, 1, ioffset[*niranks], airootloc PETSC_F90_2PTR_PARAM(pirootloc));
72*5975b3b6SBarry Smith   if (*ierr) return;
73*5975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)iranks, MPI_INT, 1, *niranks, airanks PETSC_F90_2PTR_PARAM(piranks));
74*5975b3b6SBarry Smith   if (*ierr) return;
75*5975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)ioffset, MPIU_INT, 1, *niranks + 1, aioffset PETSC_F90_2PTR_PARAM(pioffset));
76*5975b3b6SBarry Smith   if (*ierr) return;
7794a885e8SJunchao Zhang }
7894a885e8SJunchao Zhang 
79*5975b3b6SBarry Smith PETSC_EXTERN void petscsfgetrootranks_(PetscSF *sf, PetscInt *nranks, F90Array1d *aranks, F90Array1d *aroffset, F90Array1d *armine, F90Array1d *arremote, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(pranks) PETSC_F90_2PTR_PROTO(proffset) PETSC_F90_2PTR_PROTO(prmine) PETSC_F90_2PTR_PROTO(prremote))
8094a885e8SJunchao Zhang {
8194a885e8SJunchao Zhang   const PetscMPIInt *ranks   = NULL;
8294a885e8SJunchao Zhang   const PetscInt    *roffset = NULL;
8394a885e8SJunchao Zhang   const PetscInt    *rmine   = NULL;
8494a885e8SJunchao Zhang   const PetscInt    *rremote = NULL;
8594a885e8SJunchao Zhang 
86*5975b3b6SBarry Smith   *ierr = PetscSFGetRootRanks(*sf, nranks, &ranks, &roffset, &rmine, &rremote);
87*5975b3b6SBarry Smith   if (*ierr) return;
88*5975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)ranks, MPI_INT, 1, *nranks, aranks PETSC_F90_2PTR_PARAM(pranks));
89*5975b3b6SBarry Smith   if (*ierr) return;
90*5975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)roffset, MPIU_INT, 1, *nranks + 1, aroffset PETSC_F90_2PTR_PARAM(proffset));
91*5975b3b6SBarry Smith   if (*ierr) return;
92*5975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)rmine, MPIU_INT, 1, roffset[*nranks], armine PETSC_F90_2PTR_PARAM(prmine));
93*5975b3b6SBarry Smith   if (*ierr) return;
94*5975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)rremote, MPIU_INT, 1, roffset[*nranks], arremote PETSC_F90_2PTR_PARAM(prremote));
95*5975b3b6SBarry Smith   if (*ierr) return;
9694a885e8SJunchao Zhang }
9794a885e8SJunchao Zhang 
986f7e44deSSatish Balay #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
99ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
1006f7e44deSSatish Balay {
1016f7e44deSSatish Balay   MPI_Datatype dtype;
102ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1036f7e44deSSatish Balay 
104*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
105*5975b3b6SBarry Smith   if (*ierr) return;
106ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr, cop);
1076f7e44deSSatish Balay }
1086f7e44deSSatish Balay 
109ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
1106f7e44deSSatish Balay {
1116f7e44deSSatish Balay   MPI_Datatype dtype;
112ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1136f7e44deSSatish Balay 
114*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
115*5975b3b6SBarry Smith   if (*ierr) return;
116ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr, cop);
1176f7e44deSSatish Balay }
1186f7e44deSSatish Balay 
1199037d788SNicholas Arnold-Medabalimi PETSC_EXTERN void petscsfreducebegin_(PetscSF *sf, MPI_Fint *unit, const void *lptr, void *rptr, MPI_Fint *op, PetscErrorCode *ierr)
1209037d788SNicholas Arnold-Medabalimi {
1219037d788SNicholas Arnold-Medabalimi   MPI_Datatype dtype;
1229037d788SNicholas Arnold-Medabalimi   MPI_Op       cop = MPI_Op_f2c(*op);
1239037d788SNicholas Arnold-Medabalimi 
124*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
125*5975b3b6SBarry Smith   if (*ierr) return;
1269037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceBegin(*sf, dtype, lptr, rptr, cop);
1279037d788SNicholas Arnold-Medabalimi }
1289037d788SNicholas Arnold-Medabalimi 
1299037d788SNicholas Arnold-Medabalimi PETSC_EXTERN void petscsfreduceend_(PetscSF *sf, MPI_Fint *unit, const void *lptr, void *rptr, MPI_Fint *op, PetscErrorCode *ierr)
1309037d788SNicholas Arnold-Medabalimi {
1319037d788SNicholas Arnold-Medabalimi   MPI_Datatype dtype;
1329037d788SNicholas Arnold-Medabalimi   MPI_Op       cop = MPI_Op_f2c(*op);
1339037d788SNicholas Arnold-Medabalimi 
134*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
135*5975b3b6SBarry Smith   if (*ierr) return;
1369037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceEnd(*sf, dtype, lptr, rptr, cop);
1379037d788SNicholas Arnold-Medabalimi }
1389037d788SNicholas Arnold-Medabalimi 
1396f7e44deSSatish Balay #else
1406f7e44deSSatish Balay 
141ad227feaSJunchao 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))
1428af6ec1cSBarry Smith {
1438af6ec1cSBarry Smith   MPI_Datatype dtype;
1448af6ec1cSBarry Smith   const void  *rootdata;
1458af6ec1cSBarry Smith   void        *leafdata;
146ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1478af6ec1cSBarry Smith 
148*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
149*5975b3b6SBarry Smith   if (*ierr) return;
150*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
151*5975b3b6SBarry Smith   if (*ierr) return;
152*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
153*5975b3b6SBarry Smith   if (*ierr) return;
154ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata, cop);
1558af6ec1cSBarry Smith }
1568af6ec1cSBarry Smith 
157ad227feaSJunchao 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))
1588af6ec1cSBarry Smith {
1598af6ec1cSBarry Smith   MPI_Datatype dtype;
1608af6ec1cSBarry Smith   const void  *rootdata;
1618af6ec1cSBarry Smith   void        *leafdata;
162ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1638af6ec1cSBarry Smith 
164*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
165*5975b3b6SBarry Smith   if (*ierr) return;
166*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
167*5975b3b6SBarry Smith   if (*ierr) return;
168*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
169*5975b3b6SBarry Smith   if (*ierr) return;
170ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata, cop);
1718af6ec1cSBarry Smith }
1729037d788SNicholas Arnold-Medabalimi 
1739037d788SNicholas Arnold-Medabalimi PETSC_EXTERN void petscsfreducebegin_(PetscSF *sf, MPI_Fint *unit, F90Array1d *lptr, F90Array1d *rptr, MPI_Fint *op, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(lptrd) PETSC_F90_2PTR_PROTO(rptrd))
1749037d788SNicholas Arnold-Medabalimi {
1759037d788SNicholas Arnold-Medabalimi   MPI_Datatype dtype;
1769037d788SNicholas Arnold-Medabalimi   const void  *leafdata;
1779037d788SNicholas Arnold-Medabalimi   void        *rootdata;
1789037d788SNicholas Arnold-Medabalimi   MPI_Op       cop = MPI_Op_f2c(*op);
1799037d788SNicholas Arnold-Medabalimi 
180*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
181*5975b3b6SBarry Smith   if (*ierr) return;
182*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
183*5975b3b6SBarry Smith   if (*ierr) return;
184*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
185*5975b3b6SBarry Smith   if (*ierr) return;
1869037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceBegin(*sf, dtype, leafdata, rootdata, cop);
1879037d788SNicholas Arnold-Medabalimi }
1889037d788SNicholas Arnold-Medabalimi 
1899037d788SNicholas Arnold-Medabalimi PETSC_EXTERN void petscsfreduceend_(PetscSF *sf, MPI_Fint *unit, F90Array1d *lptr, F90Array1d *rptr, MPI_Fint *op, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(lptrd) PETSC_F90_2PTR_PROTO(rptrd))
1909037d788SNicholas Arnold-Medabalimi {
1919037d788SNicholas Arnold-Medabalimi   MPI_Datatype dtype;
1929037d788SNicholas Arnold-Medabalimi   const void  *leafdata;
1939037d788SNicholas Arnold-Medabalimi   void        *rootdata;
1949037d788SNicholas Arnold-Medabalimi   MPI_Op       cop = MPI_Op_f2c(*op);
1959037d788SNicholas Arnold-Medabalimi 
196*5975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
197*5975b3b6SBarry Smith   if (*ierr) return;
198*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
199*5975b3b6SBarry Smith   if (*ierr) return;
200*5975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
201*5975b3b6SBarry Smith   if (*ierr) return;
2029037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceEnd(*sf, dtype, leafdata, rootdata, cop);
2039037d788SNicholas Arnold-Medabalimi }
2049037d788SNicholas Arnold-Medabalimi 
20519caf8f3SSatish Balay PETSC_EXTERN void petscsfviewfromoptions_(PetscSF *ao, PetscObject obj, char *type, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
206fe2efc57SMark {
207fe2efc57SMark   char *t;
208fe2efc57SMark 
209fe2efc57SMark   FIXCHAR(type, len, t);
210b14c0cbaSBlaise Bourdin   CHKFORTRANNULLOBJECT(obj);
211*5975b3b6SBarry Smith   *ierr = PetscSFViewFromOptions(*ao, obj, t);
212*5975b3b6SBarry Smith   if (*ierr) return;
213fe2efc57SMark   FREECHAR(type, t);
214fe2efc57SMark }
2156f7e44deSSatish Balay 
2161fb7b255SJunchao Zhang PETSC_EXTERN void petscsfdestroy_(PetscSF *x, int *ierr)
2171fb7b255SJunchao Zhang {
2181fb7b255SJunchao Zhang   PETSC_FORTRAN_OBJECT_F_DESTROYED_TO_C_NULL(x);
219*5975b3b6SBarry Smith   *ierr = PetscSFDestroy(x);
220*5975b3b6SBarry Smith   if (*ierr) return;
2211fb7b255SJunchao Zhang   PETSC_FORTRAN_OBJECT_C_NULL_TO_F_DESTROYED(x);
2221fb7b255SJunchao Zhang }
2231fb7b255SJunchao Zhang 
2246f7e44deSSatish Balay #endif
225