16dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
26dd63270SBarry Smith #include <petscdmplex.h>
36dd63270SBarry Smith
46dd63270SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
56dd63270SBarry Smith #define dmplexgetcone_ DMPLEXGETCONE
66dd63270SBarry Smith #define dmplexrestorecone_ DMPLEXRESTORECONE
76dd63270SBarry Smith #define dmplexgetconeorientation_ DMPLEXGETCONEORIENTATION
86dd63270SBarry Smith #define dmplexrestoreconeorientation_ DMPLEXRESTORECONEORIENTATION
96dd63270SBarry Smith #define dmplexgetsupport_ DMPLEXGETSUPPORT
106dd63270SBarry Smith #define dmplexrestoresupport_ DMPLEXRESTORESUPPORT
116dd63270SBarry Smith #define dmplexgettransitiveclosure_ DMPLEXGETTRANSITIVECLOSURE
126dd63270SBarry Smith #define dmplexrestoretransitiveclosure_ DMPLEXRESTORETRANSITIVECLOSURE
136dd63270SBarry Smith #define dmplexvecgetclosure_ DMPLEXVECGETCLOSURE
146dd63270SBarry Smith #define dmplexvecrestoreclosure_ DMPLEXVECRESTORECLOSURE
156dd63270SBarry Smith #define dmplexvecsetclosure_ DMPLEXVECSETCLOSURE
166dd63270SBarry Smith #define dmplexmatsetclosure_ DMPLEXMATSETCLOSURE
17*c3871b17SNoam T #define dmplexgetclosureindices_ DMPLEXGETCLOSUREINDICES
18*c3871b17SNoam T #define dmplexrestoreclosureindices_ DMPLEXRESTORECLOSUREINDICES
196dd63270SBarry Smith #define dmplexgetjoin_ DMPLEXGETJOIN
206dd63270SBarry Smith #define dmplexgetfulljoin_ DMPLEXGETFULLJOIN
216dd63270SBarry Smith #define dmplexrestorejoin_ DMPLEXRESTOREJOIN
226dd63270SBarry Smith #define dmplexgetmeet_ DMPLEXGETMEET
236dd63270SBarry Smith #define dmplexgetfullmeet_ DMPLEXGETFULLMEET
246dd63270SBarry Smith #define dmplexrestoremeet_ DMPLEXRESTOREMEET
256dd63270SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
266dd63270SBarry Smith #define dmplexgetcone_ dmplexgetcone
276dd63270SBarry Smith #define dmplexrestorecone_ dmplexrestorecone
286dd63270SBarry Smith #define dmplexgetconeorientation_ dmplexgetconeorientation
296dd63270SBarry Smith #define dmplexrestoreconeorientation_ dmplexrestoreconeorientation
306dd63270SBarry Smith #define dmplexgetsupport_ dmplexgetsupport
316dd63270SBarry Smith #define dmplexrestoresupport_ dmplexrestoresupport
326dd63270SBarry Smith #define dmplexgettransitiveclosure_ dmplexgettransitiveclosure
336dd63270SBarry Smith #define dmplexrestoretransitiveclosure_ dmplexrestoretransitiveclosure
346dd63270SBarry Smith #define dmplexvecgetclosure_ dmplexvecgetclosure
356dd63270SBarry Smith #define dmplexvecrestoreclosure_ dmplexvecrestoreclosure
366dd63270SBarry Smith #define dmplexvecsetclosure_ dmplexvecsetclosure
376dd63270SBarry Smith #define dmplexmatsetclosure_ dmplexmatsetclosure
38*c3871b17SNoam T #define dmplexgetclosureindices_ dmplexgetclosureindices
39*c3871b17SNoam T #define dmplexrestoreclosureindices_ dmplexrestoreclosureindices
406dd63270SBarry Smith #define dmplexgetjoin_ dmplexgetjoin
416dd63270SBarry Smith #define dmplexgetfulljoin_ dmplexgetfulljoin
426dd63270SBarry Smith #define dmplexrestorejoin_ dmplexrestorejoin
436dd63270SBarry Smith #define dmplexgetmeet_ dmplexgetmeet
446dd63270SBarry Smith #define dmplexgetfullmeet_ dmplexgetfullmeet
456dd63270SBarry Smith #define dmplexrestoremeet_ dmplexrestoremeet
466dd63270SBarry Smith #endif
476dd63270SBarry Smith
dmplexgetcone_(DM * dm,PetscInt * p,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))486dd63270SBarry Smith PETSC_EXTERN void dmplexgetcone_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
496dd63270SBarry Smith {
506dd63270SBarry Smith const PetscInt *v;
516dd63270SBarry Smith PetscInt n;
526dd63270SBarry Smith
536dd63270SBarry Smith *ierr = DMPlexGetConeSize(*dm, *p, &n);
546dd63270SBarry Smith if (*ierr) return;
556dd63270SBarry Smith *ierr = DMPlexGetCone(*dm, *p, &v);
566dd63270SBarry Smith if (*ierr) return;
576dd63270SBarry Smith *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
586dd63270SBarry Smith }
596dd63270SBarry Smith
dmplexrestorecone_(DM * dm,PetscInt * p,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))606dd63270SBarry Smith PETSC_EXTERN void dmplexrestorecone_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
616dd63270SBarry Smith {
626dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
636dd63270SBarry Smith if (*ierr) return;
646dd63270SBarry Smith }
656dd63270SBarry Smith
dmplexgetconeorientation_(DM * dm,PetscInt * p,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))666dd63270SBarry Smith PETSC_EXTERN void dmplexgetconeorientation_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
676dd63270SBarry Smith {
686dd63270SBarry Smith const PetscInt *v;
696dd63270SBarry Smith PetscInt n;
706dd63270SBarry Smith
716dd63270SBarry Smith *ierr = DMPlexGetConeSize(*dm, *p, &n);
726dd63270SBarry Smith if (*ierr) return;
736dd63270SBarry Smith *ierr = DMPlexGetConeOrientation(*dm, *p, &v);
746dd63270SBarry Smith if (*ierr) return;
756dd63270SBarry Smith *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
766dd63270SBarry Smith }
776dd63270SBarry Smith
dmplexrestoreconeorientation_(DM * dm,PetscInt * p,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))786dd63270SBarry Smith PETSC_EXTERN void dmplexrestoreconeorientation_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
796dd63270SBarry Smith {
806dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
816dd63270SBarry Smith if (*ierr) return;
826dd63270SBarry Smith }
836dd63270SBarry Smith
dmplexgetsupport_(DM * dm,PetscInt * p,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))846dd63270SBarry Smith PETSC_EXTERN void dmplexgetsupport_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
856dd63270SBarry Smith {
866dd63270SBarry Smith const PetscInt *v;
876dd63270SBarry Smith PetscInt n;
886dd63270SBarry Smith
896dd63270SBarry Smith *ierr = DMPlexGetSupportSize(*dm, *p, &n);
906dd63270SBarry Smith if (*ierr) return;
916dd63270SBarry Smith *ierr = DMPlexGetSupport(*dm, *p, &v);
926dd63270SBarry Smith if (*ierr) return;
936dd63270SBarry Smith *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
946dd63270SBarry Smith }
956dd63270SBarry Smith
dmplexrestoresupport_(DM * dm,PetscInt * p,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))966dd63270SBarry Smith PETSC_EXTERN void dmplexrestoresupport_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
976dd63270SBarry Smith {
986dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
996dd63270SBarry Smith if (*ierr) return;
1006dd63270SBarry Smith }
1016dd63270SBarry Smith
dmplexgettransitiveclosure_(DM * dm,PetscInt * p,PetscBool * useCone,PetscInt * N,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))1026dd63270SBarry Smith PETSC_EXTERN void dmplexgettransitiveclosure_(DM *dm, PetscInt *p, PetscBool *useCone, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
1036dd63270SBarry Smith {
1046dd63270SBarry Smith PetscInt *v = NULL;
1056dd63270SBarry Smith PetscInt n;
1066dd63270SBarry Smith
1076dd63270SBarry Smith CHKFORTRANNULL(N);
1086dd63270SBarry Smith *ierr = DMPlexGetTransitiveClosure(*dm, *p, *useCone, &n, &v);
1096dd63270SBarry Smith if (*ierr) return;
1106dd63270SBarry Smith *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n * 2, ptr PETSC_F90_2PTR_PARAM(ptrd));
1116dd63270SBarry Smith if (N) *N = n;
1126dd63270SBarry Smith }
1136dd63270SBarry Smith
dmplexrestoretransitiveclosure_(DM * dm,PetscInt * p,PetscBool * useCone,PetscInt * N,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))1146dd63270SBarry Smith PETSC_EXTERN void dmplexrestoretransitiveclosure_(DM *dm, PetscInt *p, PetscBool *useCone, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
1156dd63270SBarry Smith {
1166dd63270SBarry Smith PetscInt *array;
1176dd63270SBarry Smith
1186dd63270SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_INT, (void **)&array PETSC_F90_2PTR_PARAM(ptrd));
1196dd63270SBarry Smith if (*ierr) return;
1206dd63270SBarry Smith *ierr = DMPlexRestoreTransitiveClosure(*dm, *p, *useCone, NULL, &array);
1216dd63270SBarry Smith if (*ierr) return;
1226dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
1236dd63270SBarry Smith if (*ierr) return;
1246dd63270SBarry Smith }
1256dd63270SBarry Smith
dmplexvecgetclosure_(DM * dm,PetscSection * section,Vec * x,PetscInt * point,PetscInt * N,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))1266dd63270SBarry Smith PETSC_EXTERN void dmplexvecgetclosure_(DM *dm, PetscSection *section, Vec *x, PetscInt *point, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
1276dd63270SBarry Smith {
1286dd63270SBarry Smith PetscScalar *v = NULL;
1296dd63270SBarry Smith PetscInt n;
1306dd63270SBarry Smith
1316dd63270SBarry Smith CHKFORTRANNULL(N);
1326dd63270SBarry Smith *ierr = DMPlexVecGetClosure(*dm, *section, *x, *point, &n, &v);
1336dd63270SBarry Smith if (*ierr) return;
1346dd63270SBarry Smith *ierr = F90Array1dCreate((void *)v, MPIU_SCALAR, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
1356dd63270SBarry Smith if (N) *N = n;
1366dd63270SBarry Smith }
1376dd63270SBarry Smith
dmplexvecrestoreclosure_(DM * dm,PetscSection * section,Vec * v,PetscInt * point,PetscInt * N,F90Array1d * ptr,int * ierr PETSC_F90_2PTR_PROTO (ptrd))1386dd63270SBarry Smith PETSC_EXTERN void dmplexvecrestoreclosure_(DM *dm, PetscSection *section, Vec *v, PetscInt *point, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
1396dd63270SBarry Smith {
1406dd63270SBarry Smith PetscScalar *array;
1416dd63270SBarry Smith
1426dd63270SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&array PETSC_F90_2PTR_PARAM(ptrd));
1436dd63270SBarry Smith if (*ierr) return;
1446dd63270SBarry Smith *ierr = DMPlexVecRestoreClosure(*dm, *section, *v, *point, NULL, &array);
1456dd63270SBarry Smith if (*ierr) return;
1466dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
1476dd63270SBarry Smith if (*ierr) return;
1486dd63270SBarry Smith }
1496dd63270SBarry Smith
dmplexgetclosureindices_(DM * dm,PetscSection * section,PetscSection * idxSection,PetscInt * point,PetscBool * useClPerm,PetscInt * numIndices,F90Array1d * idxPtr,PetscInt * outOffsets,F90Array1d * valPtr,int * ierr PETSC_F90_2PTR_PROTO (idxPtrd)PETSC_F90_2PTR_PROTO (valPtrd))150*c3871b17SNoam T PETSC_EXTERN void dmplexgetclosureindices_(DM *dm, PetscSection *section, PetscSection *idxSection, PetscInt *point, PetscBool *useClPerm, PetscInt *numIndices, F90Array1d *idxPtr, PetscInt *outOffsets, F90Array1d *valPtr, int *ierr PETSC_F90_2PTR_PROTO(idxPtrd) PETSC_F90_2PTR_PROTO(valPtrd))
151*c3871b17SNoam T {
152*c3871b17SNoam T PetscInt *indices;
153*c3871b17SNoam T PetscScalar *values;
154*c3871b17SNoam T
155*c3871b17SNoam T CHKFORTRANNULL(outOffsets);
156*c3871b17SNoam T if (FORTRANNULLSCALARPOINTER(valPtr)) *ierr = DMPlexGetClosureIndices(*dm, *section, *idxSection, *point, *useClPerm, numIndices, &indices, outOffsets, NULL);
157*c3871b17SNoam T else *ierr = DMPlexGetClosureIndices(*dm, *section, *idxSection, *point, *useClPerm, numIndices, &indices, outOffsets, &values);
158*c3871b17SNoam T if (*ierr) return;
159*c3871b17SNoam T *ierr = F90Array1dCreate((void *)indices, MPIU_INT, 1, *numIndices, idxPtr PETSC_F90_2PTR_PARAM(idxPtrd));
160*c3871b17SNoam T if (*ierr) return;
161*c3871b17SNoam T if (FORTRANNULLSCALARPOINTER(valPtr)) *ierr = F90Array1dCreate((void *)values, MPIU_SCALAR, 1, *numIndices, valPtr PETSC_F90_2PTR_PARAM(valPtrd));
162*c3871b17SNoam T }
163*c3871b17SNoam T
dmplexrestoreclosureindices_(DM * dm,PetscSection * section,PetscSection * idxSection,PetscInt * point,PetscBool * useClPerm,PetscInt * numIndices,F90Array1d * idxPtr,PetscInt * outOffsets,F90Array1d * valPtr,int * ierr PETSC_F90_2PTR_PROTO (idxPtrd)PETSC_F90_2PTR_PROTO (valPtrd))164*c3871b17SNoam T PETSC_EXTERN void dmplexrestoreclosureindices_(DM *dm, PetscSection *section, PetscSection *idxSection, PetscInt *point, PetscBool *useClPerm, PetscInt *numIndices, F90Array1d *idxPtr, PetscInt *outOffsets, F90Array1d *valPtr, int *ierr PETSC_F90_2PTR_PROTO(idxPtrd) PETSC_F90_2PTR_PROTO(valPtrd))
165*c3871b17SNoam T {
166*c3871b17SNoam T PetscInt *indices;
167*c3871b17SNoam T PetscScalar *values = NULL;
168*c3871b17SNoam T
169*c3871b17SNoam T CHKFORTRANNULL(outOffsets);
170*c3871b17SNoam T *ierr = F90Array1dAccess(idxPtr, MPIU_INT, (void **)&indices PETSC_F90_2PTR_PARAM(idxPtrd));
171*c3871b17SNoam T if (*ierr) return;
172*c3871b17SNoam T if (!FORTRANNULLSCALARPOINTER(valPtr)) {
173*c3871b17SNoam T *ierr = F90Array1dAccess(valPtr, MPIU_SCALAR, (void **)&values PETSC_F90_2PTR_PARAM(valPtrd));
174*c3871b17SNoam T if (*ierr) return;
175*c3871b17SNoam T *ierr = DMPlexRestoreClosureIndices(*dm, *section, *idxSection, *point, *useClPerm, numIndices, &indices, outOffsets, &values);
176*c3871b17SNoam T } else *ierr = DMPlexRestoreClosureIndices(*dm, *section, *idxSection, *point, *useClPerm, numIndices, &indices, outOffsets, NULL);
177*c3871b17SNoam T if (*ierr) return;
178*c3871b17SNoam T *ierr = F90Array1dDestroy(idxPtr, MPIU_INT PETSC_F90_2PTR_PARAM(idxPtrd));
179*c3871b17SNoam T if (*ierr) return;
180*c3871b17SNoam T if (!FORTRANNULLSCALARPOINTER(valPtr)) *ierr = F90Array1dDestroy(valPtr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(valPtrd));
181*c3871b17SNoam T }
182*c3871b17SNoam T
dmplexgetjoin_(DM * dm,PetscInt * numPoints,PetscInt * points,PetscInt * N,F90Array1d * cptr,int * ierr PETSC_F90_2PTR_PROTO (cptrd))1836dd63270SBarry Smith PETSC_EXTERN void dmplexgetjoin_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
1846dd63270SBarry Smith {
1856dd63270SBarry Smith const PetscInt *coveredPoints;
1866dd63270SBarry Smith PetscInt n;
1876dd63270SBarry Smith
1886dd63270SBarry Smith CHKFORTRANNULL(N);
1896dd63270SBarry Smith *ierr = DMPlexGetJoin(*dm, *numPoints, points, &n, &coveredPoints);
1906dd63270SBarry Smith if (*ierr) return;
1916dd63270SBarry Smith *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
1926dd63270SBarry Smith if (N) *N = n;
1936dd63270SBarry Smith }
1946dd63270SBarry Smith
dmplexgetfulljoin_(DM * dm,PetscInt * numPoints,PetscInt * points,PetscInt * N,F90Array1d * cptr,int * ierr PETSC_F90_2PTR_PROTO (cptrd))1956dd63270SBarry Smith PETSC_EXTERN void dmplexgetfulljoin_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
1966dd63270SBarry Smith {
1976dd63270SBarry Smith const PetscInt *coveredPoints;
1986dd63270SBarry Smith PetscInt n;
1996dd63270SBarry Smith
2006dd63270SBarry Smith CHKFORTRANNULL(N);
2016dd63270SBarry Smith *ierr = DMPlexGetFullJoin(*dm, *numPoints, points, &n, &coveredPoints);
2026dd63270SBarry Smith if (*ierr) return;
2036dd63270SBarry Smith *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
2046dd63270SBarry Smith if (N) *N = n;
2056dd63270SBarry Smith }
2066dd63270SBarry Smith
dmplexrestorejoin_(DM * dm,PetscInt * numPoints,PetscInt * points,PetscInt * N,F90Array1d * cptr,int * ierr PETSC_F90_2PTR_PROTO (cptrd))2076dd63270SBarry Smith PETSC_EXTERN void dmplexrestorejoin_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
2086dd63270SBarry Smith {
2096dd63270SBarry Smith PetscInt *coveredPoints;
2106dd63270SBarry Smith
2116dd63270SBarry Smith *ierr = F90Array1dAccess(cptr, MPIU_INT, (void **)&coveredPoints PETSC_F90_2PTR_PARAM(cptrd));
2126dd63270SBarry Smith if (*ierr) return;
2136dd63270SBarry Smith *ierr = DMPlexRestoreJoin(*dm, 0, NULL, NULL, (const PetscInt **)&coveredPoints);
2146dd63270SBarry Smith if (*ierr) return;
2156dd63270SBarry Smith *ierr = F90Array1dDestroy(cptr, MPIU_INT PETSC_F90_2PTR_PARAM(cptrd));
2166dd63270SBarry Smith if (*ierr) return;
2176dd63270SBarry Smith }
2186dd63270SBarry Smith
dmplexgetmeet_(DM * dm,PetscInt * numPoints,PetscInt * points,PetscInt * N,F90Array1d * cptr,int * ierr PETSC_F90_2PTR_PROTO (cptrd))2196dd63270SBarry Smith PETSC_EXTERN void dmplexgetmeet_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
2206dd63270SBarry Smith {
2216dd63270SBarry Smith const PetscInt *coveredPoints;
2226dd63270SBarry Smith PetscInt n;
2236dd63270SBarry Smith
2246dd63270SBarry Smith CHKFORTRANNULL(N);
2256dd63270SBarry Smith *ierr = DMPlexGetMeet(*dm, *numPoints, points, &n, &coveredPoints);
2266dd63270SBarry Smith if (*ierr) return;
2276dd63270SBarry Smith *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
2286dd63270SBarry Smith if (N) *N = n;
2296dd63270SBarry Smith }
2306dd63270SBarry Smith
dmplexgetfullmeet_(DM * dm,PetscInt * numPoints,PetscInt * points,PetscInt * N,F90Array1d * cptr,int * ierr PETSC_F90_2PTR_PROTO (cptrd))2316dd63270SBarry Smith PETSC_EXTERN void dmplexgetfullmeet_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
2326dd63270SBarry Smith {
2336dd63270SBarry Smith const PetscInt *coveredPoints;
2346dd63270SBarry Smith PetscInt n;
2356dd63270SBarry Smith
2366dd63270SBarry Smith CHKFORTRANNULL(N);
2376dd63270SBarry Smith if (*ierr) return;
2386dd63270SBarry Smith *ierr = DMPlexGetFullMeet(*dm, *numPoints, points, &n, &coveredPoints);
2396dd63270SBarry Smith if (*ierr) return;
2406dd63270SBarry Smith *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
2416dd63270SBarry Smith if (N) *N = n;
2426dd63270SBarry Smith }
2436dd63270SBarry Smith
dmplexrestoremeet_(DM * dm,PetscInt * numPoints,PetscInt * points,PetscInt * N,F90Array1d * cptr,int * ierr PETSC_F90_2PTR_PROTO (cptrd))2446dd63270SBarry Smith PETSC_EXTERN void dmplexrestoremeet_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
2456dd63270SBarry Smith {
2466dd63270SBarry Smith PetscInt *coveredPoints;
2476dd63270SBarry Smith
2486dd63270SBarry Smith *ierr = F90Array1dAccess(cptr, MPIU_INT, (void **)&coveredPoints PETSC_F90_2PTR_PARAM(cptrd));
2496dd63270SBarry Smith if (*ierr) return;
2506dd63270SBarry Smith *ierr = DMPlexRestoreMeet(*dm, 0, NULL, NULL, (const PetscInt **)&coveredPoints);
2516dd63270SBarry Smith if (*ierr) return;
2526dd63270SBarry Smith *ierr = F90Array1dDestroy(cptr, MPIU_INT PETSC_F90_2PTR_PARAM(cptrd));
2536dd63270SBarry Smith if (*ierr) return;
2546dd63270SBarry Smith }
255