xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision ffeef943c8ee50edff320d8a3135bb0c94853e4c)
1 #include <petsc/private/f90impl.h>
2 #include <petsc/private/sfimpl.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5   #define petscsfgetgraph_      PETSCSFGETGRAPH
6   #define petscsfbcastbegin_    PETSCSFBCASTBEGIN
7   #define petscsfbcastend_      PETSCSFBCASTEND
8   #define petscsfreducebegin_   PETSCSFREDUCEBEGIN
9   #define petscsfreduceend_     PETSCSFREDUCEEND
10   #define f90arraysfnodecreate_ F90ARRAYSFNODECREATE
11   #define petscsfdestroy_       PETSCSFDESTROY
12   #define petscsfsetgraph_      PETSCSFSETGRAPH
13   #define petscsfgetleafranks_  PETSCSFGETLEAFRANKS
14   #define petscsfgetrootranks_  PETSCSFGETROOTRANKS
15 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
16   #define petscsfgetgraph_      petscsfgetgraph
17   #define petscsfbcastbegin_    petscsfbcastbegin
18   #define petscsfbcastend_      petscsfbcastend
19   #define petscsfreducebegin_   petscsfreducebegin
20   #define petscsfreduceend_     petscsfreduceend
21   #define f90arraysfnodecreate_ f90arraysfnodecreate
22   #define petscsfdestroy_       petscsfdestroy
23   #define petscsfsetgraph_      petscsfsetgraph
24   #define petscsfgetleafranks_  petscsfgetleafranks
25   #define petscsfgetrootranks_  petscsfgetrootranks
26 #endif
27 
28 PETSC_EXTERN void f90arraysfnodecreate_(const PetscInt *, PetscInt *, void *PETSC_F90_2PTR_PROTO_NOVAR);
29 
30 PETSC_EXTERN void petscsfsetgraph_(PetscSF *sf, PetscInt *nroots, PetscInt *nleaves, PetscInt *ilocal, PetscCopyMode *localmode, PetscSFNode *iremote, PetscCopyMode *remotemode, int *ierr)
31 {
32   if (ilocal == PETSC_NULL_INTEGER_Fortran) ilocal = NULL;
33   *ierr = PetscSFSetGraph(*sf, *nroots, *nleaves, ilocal, *localmode, iremote, *remotemode);
34 }
35 
36 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))
37 {
38   const PetscInt    *ilocal;
39   const PetscSFNode *iremote;
40   PetscInt           nl;
41 
42   *ierr = PetscSFGetGraph(*sf, nroots, nleaves, &ilocal, &iremote);
43   if (*ierr) return;
44   nl = *nleaves;
45   if (!ilocal) nl = 0;
46   *ierr = F90Array1dCreate((void *)ilocal, MPIU_INT, 1, nl, ailocal PETSC_F90_2PTR_PARAM(pilocal));
47   /* this creates a memory leak */
48   f90arraysfnodecreate_((PetscInt *)iremote, nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
49 }
50 
51 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))
52 {
53   const PetscMPIInt *iranks   = NULL;
54   const PetscInt    *ioffset  = NULL;
55   const PetscInt    *irootloc = NULL;
56 
57   *ierr = PetscSFGetLeafRanks(*sf, niranks, &iranks, &ioffset, &irootloc);
58   if (*ierr) return;
59   *ierr = F90Array1dCreate((void *)irootloc, MPIU_INT, 1, ioffset[*niranks], airootloc PETSC_F90_2PTR_PARAM(pirootloc));
60   if (*ierr) return;
61   *ierr = F90Array1dCreate((void *)iranks, MPI_INT, 1, *niranks, airanks PETSC_F90_2PTR_PARAM(piranks));
62   if (*ierr) return;
63   *ierr = F90Array1dCreate((void *)ioffset, MPIU_INT, 1, *niranks + 1, aioffset PETSC_F90_2PTR_PARAM(pioffset));
64   if (*ierr) return;
65 }
66 
67 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))
68 {
69   const PetscMPIInt *ranks   = NULL;
70   const PetscInt    *roffset = NULL;
71   const PetscInt    *rmine   = NULL;
72   const PetscInt    *rremote = NULL;
73 
74   *ierr = PetscSFGetRootRanks(*sf, nranks, &ranks, &roffset, &rmine, &rremote);
75   if (*ierr) return;
76   *ierr = F90Array1dCreate((void *)ranks, MPI_INT, 1, *nranks, aranks PETSC_F90_2PTR_PARAM(pranks));
77   if (*ierr) return;
78   *ierr = F90Array1dCreate((void *)roffset, MPIU_INT, 1, *nranks + 1, aroffset PETSC_F90_2PTR_PARAM(proffset));
79   if (*ierr) return;
80   *ierr = F90Array1dCreate((void *)rmine, MPIU_INT, 1, roffset[*nranks], armine PETSC_F90_2PTR_PARAM(prmine));
81   if (*ierr) return;
82   *ierr = F90Array1dCreate((void *)rremote, MPIU_INT, 1, roffset[*nranks], arremote PETSC_F90_2PTR_PARAM(prremote));
83   if (*ierr) return;
84 }
85 
86 #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
87 PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
88 {
89   MPI_Datatype dtype;
90   MPI_Op       cop = MPI_Op_f2c(*op);
91 
92   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
93   if (*ierr) return;
94   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr, cop);
95 }
96 
97 PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
98 {
99   MPI_Datatype dtype;
100   MPI_Op       cop = MPI_Op_f2c(*op);
101 
102   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
103   if (*ierr) return;
104   *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr, cop);
105 }
106 
107 PETSC_EXTERN void petscsfreducebegin_(PetscSF *sf, MPI_Fint *unit, const void *lptr, void *rptr, MPI_Fint *op, PetscErrorCode *ierr)
108 {
109   MPI_Datatype dtype;
110   MPI_Op       cop = MPI_Op_f2c(*op);
111 
112   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
113   if (*ierr) return;
114   *ierr = PetscSFReduceBegin(*sf, dtype, lptr, rptr, cop);
115 }
116 
117 PETSC_EXTERN void petscsfreduceend_(PetscSF *sf, MPI_Fint *unit, const void *lptr, void *rptr, MPI_Fint *op, PetscErrorCode *ierr)
118 {
119   MPI_Datatype dtype;
120   MPI_Op       cop = MPI_Op_f2c(*op);
121 
122   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
123   if (*ierr) return;
124   *ierr = PetscSFReduceEnd(*sf, dtype, lptr, rptr, cop);
125 }
126 
127 #else
128 
129 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))
130 {
131   MPI_Datatype dtype;
132   const void  *rootdata;
133   void        *leafdata;
134   MPI_Op       cop = MPI_Op_f2c(*op);
135 
136   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
137   if (*ierr) return;
138   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
139   if (*ierr) return;
140   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
141   if (*ierr) return;
142   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata, cop);
143 }
144 
145 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))
146 {
147   MPI_Datatype dtype;
148   const void  *rootdata;
149   void        *leafdata;
150   MPI_Op       cop = MPI_Op_f2c(*op);
151 
152   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
153   if (*ierr) return;
154   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
155   if (*ierr) return;
156   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
157   if (*ierr) return;
158   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata, cop);
159 }
160 
161 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))
162 {
163   MPI_Datatype dtype;
164   const void  *leafdata;
165   void        *rootdata;
166   MPI_Op       cop = MPI_Op_f2c(*op);
167 
168   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
169   if (*ierr) return;
170   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
171   if (*ierr) return;
172   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
173   if (*ierr) return;
174   *ierr = PetscSFReduceBegin(*sf, dtype, leafdata, rootdata, cop);
175 }
176 
177 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))
178 {
179   MPI_Datatype dtype;
180   const void  *leafdata;
181   void        *rootdata;
182   MPI_Op       cop = MPI_Op_f2c(*op);
183 
184   *ierr = PetscMPIFortranDatatypeToC(*unit, &dtype);
185   if (*ierr) return;
186   *ierr = F90Array1dAccess(rptr, dtype, (void **)&rootdata PETSC_F90_2PTR_PARAM(rptrd));
187   if (*ierr) return;
188   *ierr = F90Array1dAccess(lptr, dtype, (void **)&leafdata PETSC_F90_2PTR_PARAM(lptrd));
189   if (*ierr) return;
190   *ierr = PetscSFReduceEnd(*sf, dtype, leafdata, rootdata, cop);
191 }
192 
193 PETSC_EXTERN void petscsfdestroy_(PetscSF *x, int *ierr)
194 {
195   PETSC_FORTRAN_OBJECT_F_DESTROYED_TO_C_NULL(x);
196   *ierr = PetscSFDestroy(x);
197   if (*ierr) return;
198   PETSC_FORTRAN_OBJECT_C_NULL_TO_F_DESTROYED(x);
199 }
200 
201 #endif
202