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