xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision b0dcfd164860a975c76f90dabf1036901aab1c4e)
1*6dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
2a6e9e4f7SMatthew G. Knepley #include <petsc/private/sfimpl.h>
3a6e9e4f7SMatthew G. Knepley 
4a6e9e4f7SMatthew G. Knepley #if defined(PETSC_HAVE_FORTRAN_CAPS)
5ce78bad3SBarry Smith   #define petscsfrestoregraph_     PETSCSFRESTOREGRAPH
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
1194a885e8SJunchao Zhang   #define petscsfgetleafranks_     PETSCSFGETLEAFRANKS
1294a885e8SJunchao Zhang   #define petscsfgetrootranks_     PETSCSFGETROOTRANKS
13ce78bad3SBarry Smith   #define f90array1dcreatesfnode_  F90ARRAY1DCREATESFNODE
14ce78bad3SBarry Smith   #define f90array1ddestroysfnode_ F90ARRAY1DDESTROYSFNODE
15a6e9e4f7SMatthew G. Knepley #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
16ce78bad3SBarry Smith   #define petscsfrestoregraph_     petscsfrestoregraph
175c87a30dSBarry Smith   #define petscsfgetgraph_         petscsfgetgraph
188af6ec1cSBarry Smith   #define petscsfbcastbegin_       petscsfbcastbegin
198af6ec1cSBarry Smith   #define petscsfbcastend_         petscsfbcastend
209037d788SNicholas Arnold-Medabalimi   #define petscsfreducebegin_      petscsfreducebegin
219037d788SNicholas Arnold-Medabalimi   #define petscsfreduceend_        petscsfreduceend
2294a885e8SJunchao Zhang   #define petscsfgetleafranks_     petscsfgetleafranks
2394a885e8SJunchao Zhang   #define petscsfgetrootranks_     petscsfgetrootranks
2465960ac8SPierre Jolivet   #define f90array1dcreatesfnode_  f90array1dcreatesfnode
2565960ac8SPierre Jolivet   #define f90array1ddestroysfnode_ f90array1ddestroysfnode
26a6e9e4f7SMatthew G. Knepley #endif
27a6e9e4f7SMatthew G. Knepley 
28ce78bad3SBarry Smith PETSC_EXTERN void f90array1dcreatesfnode_(const PetscSFNode *, PetscInt *, PetscInt *, void *PETSC_F90_2PTR_PROTO_NOVAR);
29ce78bad3SBarry Smith PETSC_EXTERN void f90array1ddestroysfnode_(void *PETSC_F90_2PTR_PROTO_NOVAR);
307f139299SBarry Smith 
petscsfgetgraph_(PetscSF * sf,PetscInt * nroots,PetscInt * nleaves,F90Array1d * ailocal,F90Array1d * airemote,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (pilocal)PETSC_F90_2PTR_PROTO (piremote))3161bf59e3SJunchao 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))
325c87a30dSBarry Smith {
335c87a30dSBarry Smith   const PetscInt    *ilocal;
345c87a30dSBarry Smith   const PetscSFNode *iremote;
35ce78bad3SBarry Smith   PetscInt           nl, one = 1;
365c87a30dSBarry Smith 
375975b3b6SBarry Smith   *ierr = PetscSFGetGraph(*sf, nroots, nleaves, &ilocal, &iremote);
385975b3b6SBarry Smith   if (*ierr) return;
398dbb0df6SBarry Smith   nl = *nleaves;
408dbb0df6SBarry Smith   if (!ilocal) nl = 0;
418dbb0df6SBarry Smith   *ierr = F90Array1dCreate((void *)ilocal, MPIU_INT, 1, nl, ailocal PETSC_F90_2PTR_PARAM(pilocal));
427f139299SBarry Smith   /* this creates a memory leak */
43ce78bad3SBarry Smith   f90array1dcreatesfnode_(iremote, &one, nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
44ce78bad3SBarry Smith }
45ce78bad3SBarry Smith 
petscsfrestoregraph_(PetscSF * sf,PetscInt * nroots,PetscInt * nleaves,F90Array1d * ailocal,F90Array1d * airemote,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (pilocal)PETSC_F90_2PTR_PROTO (piremote))46ce78bad3SBarry Smith PETSC_EXTERN void petscsfrestoregraph_(PetscSF *sf, PetscInt *nroots, PetscInt *nleaves, F90Array1d *ailocal, F90Array1d *airemote, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(pilocal) PETSC_F90_2PTR_PROTO(piremote))
47ce78bad3SBarry Smith {
48ce78bad3SBarry Smith   *ierr = F90Array1dDestroy(ailocal, MPIU_INT PETSC_F90_2PTR_PARAM(pilocal));
49ce78bad3SBarry Smith   if (*ierr) return;
50ce78bad3SBarry Smith   f90array1ddestroysfnode_(airemote PETSC_F90_2PTR_PARAM(piremote));
515c87a30dSBarry Smith }
525c87a30dSBarry Smith 
petscsfgetleafranks_(PetscSF * sf,PetscMPIInt * niranks,F90Array1d * airanks,F90Array1d * aioffset,F90Array1d * airootloc,PetscErrorCode * ierr PETSC_F90_2PTR_PROTO (piranks)PETSC_F90_2PTR_PROTO (pioffset)PETSC_F90_2PTR_PROTO (pirootloc))536497c311SBarry Smith PETSC_EXTERN void petscsfgetleafranks_(PetscSF *sf, PetscMPIInt *niranks, F90Array1d *airanks, F90Array1d *aioffset, F90Array1d *airootloc, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(piranks) PETSC_F90_2PTR_PROTO(pioffset) PETSC_F90_2PTR_PROTO(pirootloc))
5494a885e8SJunchao Zhang {
5594a885e8SJunchao Zhang   const PetscMPIInt *iranks   = NULL;
5694a885e8SJunchao Zhang   const PetscInt    *ioffset  = NULL;
5794a885e8SJunchao Zhang   const PetscInt    *irootloc = NULL;
5894a885e8SJunchao Zhang 
595975b3b6SBarry Smith   *ierr = PetscSFGetLeafRanks(*sf, niranks, &iranks, &ioffset, &irootloc);
605975b3b6SBarry Smith   if (*ierr) return;
615975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)irootloc, MPIU_INT, 1, ioffset[*niranks], airootloc PETSC_F90_2PTR_PARAM(pirootloc));
625975b3b6SBarry Smith   if (*ierr) return;
635975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)iranks, MPI_INT, 1, *niranks, airanks PETSC_F90_2PTR_PARAM(piranks));
645975b3b6SBarry Smith   if (*ierr) return;
655975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)ioffset, MPIU_INT, 1, *niranks + 1, aioffset PETSC_F90_2PTR_PARAM(pioffset));
665975b3b6SBarry Smith   if (*ierr) return;
6794a885e8SJunchao Zhang }
6894a885e8SJunchao Zhang 
petscsfgetrootranks_(PetscSF * sf,PetscMPIInt * 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))696497c311SBarry Smith PETSC_EXTERN void petscsfgetrootranks_(PetscSF *sf, PetscMPIInt *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))
7094a885e8SJunchao Zhang {
7194a885e8SJunchao Zhang   const PetscMPIInt *ranks   = NULL;
7294a885e8SJunchao Zhang   const PetscInt    *roffset = NULL;
7394a885e8SJunchao Zhang   const PetscInt    *rmine   = NULL;
7494a885e8SJunchao Zhang   const PetscInt    *rremote = NULL;
7594a885e8SJunchao Zhang 
765975b3b6SBarry Smith   *ierr = PetscSFGetRootRanks(*sf, nranks, &ranks, &roffset, &rmine, &rremote);
775975b3b6SBarry Smith   if (*ierr) return;
785975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)ranks, MPI_INT, 1, *nranks, aranks PETSC_F90_2PTR_PARAM(pranks));
795975b3b6SBarry Smith   if (*ierr) return;
805975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)roffset, MPIU_INT, 1, *nranks + 1, aroffset PETSC_F90_2PTR_PARAM(proffset));
815975b3b6SBarry Smith   if (*ierr) return;
825975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)rmine, MPIU_INT, 1, roffset[*nranks], armine PETSC_F90_2PTR_PARAM(prmine));
835975b3b6SBarry Smith   if (*ierr) return;
845975b3b6SBarry Smith   *ierr = F90Array1dCreate((void *)rremote, MPIU_INT, 1, roffset[*nranks], arremote PETSC_F90_2PTR_PARAM(prremote));
855975b3b6SBarry Smith   if (*ierr) return;
8694a885e8SJunchao Zhang }
8794a885e8SJunchao Zhang 
886f7e44deSSatish Balay #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
petscsfbcastbegin_(PetscSF * sf,MPI_Fint * unit,const void * rptr,void * lptr,MPI_Fint * op,PetscErrorCode * ierr)89ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
906f7e44deSSatish Balay {
916f7e44deSSatish Balay   MPI_Datatype dtype;
92ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
936f7e44deSSatish Balay 
945975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
955975b3b6SBarry Smith   if (*ierr) return;
96ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr, cop);
976f7e44deSSatish Balay }
986f7e44deSSatish Balay 
petscsfbcastend_(PetscSF * sf,MPI_Fint * unit,const void * rptr,void * lptr,MPI_Fint * op,PetscErrorCode * ierr)99ad227feaSJunchao Zhang PETSC_EXTERN void petscsfbcastend_(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 
1045975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
1055975b3b6SBarry Smith   if (*ierr) return;
106ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr, cop);
1076f7e44deSSatish Balay }
1086f7e44deSSatish Balay 
petscsfreducebegin_(PetscSF * sf,MPI_Fint * unit,const void * lptr,void * rptr,MPI_Fint * op,PetscErrorCode * ierr)1099037d788SNicholas Arnold-Medabalimi PETSC_EXTERN void petscsfreducebegin_(PetscSF *sf, MPI_Fint *unit, const void *lptr, void *rptr, MPI_Fint *op, PetscErrorCode *ierr)
1109037d788SNicholas Arnold-Medabalimi {
1119037d788SNicholas Arnold-Medabalimi   MPI_Datatype dtype;
1129037d788SNicholas Arnold-Medabalimi   MPI_Op       cop = MPI_Op_f2c(*op);
1139037d788SNicholas Arnold-Medabalimi 
1145975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
1155975b3b6SBarry Smith   if (*ierr) return;
1169037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceBegin(*sf, dtype, lptr, rptr, cop);
1179037d788SNicholas Arnold-Medabalimi }
1189037d788SNicholas Arnold-Medabalimi 
petscsfreduceend_(PetscSF * sf,MPI_Fint * unit,const void * lptr,void * rptr,MPI_Fint * op,PetscErrorCode * ierr)1199037d788SNicholas Arnold-Medabalimi PETSC_EXTERN void petscsfreduceend_(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 
1245975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
1255975b3b6SBarry Smith   if (*ierr) return;
1269037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceEnd(*sf, dtype, lptr, rptr, cop);
1279037d788SNicholas Arnold-Medabalimi }
1289037d788SNicholas Arnold-Medabalimi 
1296f7e44deSSatish Balay #else
1306f7e44deSSatish Balay 
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))131ad227feaSJunchao 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))
1328af6ec1cSBarry Smith {
1338af6ec1cSBarry Smith   MPI_Datatype dtype;
1348af6ec1cSBarry Smith   const void  *rootdata;
1358af6ec1cSBarry Smith   void        *leafdata;
136ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1378af6ec1cSBarry Smith 
1385975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
1395975b3b6SBarry Smith   if (*ierr) return;
1405975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
1415975b3b6SBarry Smith   if (*ierr) return;
1425975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
1435975b3b6SBarry Smith   if (*ierr) return;
144ad227feaSJunchao Zhang   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata, cop);
1458af6ec1cSBarry Smith }
1468af6ec1cSBarry Smith 
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))147ad227feaSJunchao 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))
1488af6ec1cSBarry Smith {
1498af6ec1cSBarry Smith   MPI_Datatype dtype;
1508af6ec1cSBarry Smith   const void  *rootdata;
1518af6ec1cSBarry Smith   void        *leafdata;
152ad227feaSJunchao Zhang   MPI_Op       cop = MPI_Op_f2c(*op);
1538af6ec1cSBarry Smith 
1545975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
1555975b3b6SBarry Smith   if (*ierr) return;
1565975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
1575975b3b6SBarry Smith   if (*ierr) return;
1585975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
1595975b3b6SBarry Smith   if (*ierr) return;
160ad227feaSJunchao Zhang   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata, cop);
1618af6ec1cSBarry Smith }
1629037d788SNicholas Arnold-Medabalimi 
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))1639037d788SNicholas 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))
1649037d788SNicholas Arnold-Medabalimi {
1659037d788SNicholas Arnold-Medabalimi   MPI_Datatype dtype;
1669037d788SNicholas Arnold-Medabalimi   const void  *leafdata;
1679037d788SNicholas Arnold-Medabalimi   void        *rootdata;
1689037d788SNicholas Arnold-Medabalimi   MPI_Op       cop = MPI_Op_f2c(*op);
1699037d788SNicholas Arnold-Medabalimi 
1705975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
1715975b3b6SBarry Smith   if (*ierr) return;
1725975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
1735975b3b6SBarry Smith   if (*ierr) return;
1745975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
1755975b3b6SBarry Smith   if (*ierr) return;
1769037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceBegin(*sf, dtype, leafdata, rootdata, cop);
1779037d788SNicholas Arnold-Medabalimi }
1789037d788SNicholas Arnold-Medabalimi 
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))1799037d788SNicholas 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))
1809037d788SNicholas Arnold-Medabalimi {
1819037d788SNicholas Arnold-Medabalimi   MPI_Datatype dtype;
1829037d788SNicholas Arnold-Medabalimi   const void  *leafdata;
1839037d788SNicholas Arnold-Medabalimi   void        *rootdata;
1849037d788SNicholas Arnold-Medabalimi   MPI_Op       cop = MPI_Op_f2c(*op);
1859037d788SNicholas Arnold-Medabalimi 
1865975b3b6SBarry Smith   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
1875975b3b6SBarry Smith   if (*ierr) return;
1885975b3b6SBarry Smith   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
1895975b3b6SBarry Smith   if (*ierr) return;
1905975b3b6SBarry Smith   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
1915975b3b6SBarry Smith   if (*ierr) return;
1929037d788SNicholas Arnold-Medabalimi   *ierr = PetscSFReduceEnd(*sf, dtype, leafdata, rootdata, cop);
1939037d788SNicholas Arnold-Medabalimi }
1946f7e44deSSatish Balay #endif
195