xref: /petsc/src/dm/impls/plex/ftn-custom/zplexf90.c (revision 6dd63270497ad23dcf16ae500a87ff2b2a0b7474)
1*6dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
2*6dd63270SBarry Smith #include <petscdmplex.h>
3*6dd63270SBarry Smith 
4*6dd63270SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
5*6dd63270SBarry Smith   #define dmplexgetcone_                  DMPLEXGETCONE
6*6dd63270SBarry Smith   #define dmplexrestorecone_              DMPLEXRESTORECONE
7*6dd63270SBarry Smith   #define dmplexgetconeorientation_       DMPLEXGETCONEORIENTATION
8*6dd63270SBarry Smith   #define dmplexrestoreconeorientation_   DMPLEXRESTORECONEORIENTATION
9*6dd63270SBarry Smith   #define dmplexgetsupport_               DMPLEXGETSUPPORT
10*6dd63270SBarry Smith   #define dmplexrestoresupport_           DMPLEXRESTORESUPPORT
11*6dd63270SBarry Smith   #define dmplexgettransitiveclosure_     DMPLEXGETTRANSITIVECLOSURE
12*6dd63270SBarry Smith   #define dmplexrestoretransitiveclosure_ DMPLEXRESTORETRANSITIVECLOSURE
13*6dd63270SBarry Smith   #define dmplexvecgetclosure_            DMPLEXVECGETCLOSURE
14*6dd63270SBarry Smith   #define dmplexvecrestoreclosure_        DMPLEXVECRESTORECLOSURE
15*6dd63270SBarry Smith   #define dmplexvecsetclosure_            DMPLEXVECSETCLOSURE
16*6dd63270SBarry Smith   #define dmplexmatsetclosure_            DMPLEXMATSETCLOSURE
17*6dd63270SBarry Smith   #define dmplexgetjoin_                  DMPLEXGETJOIN
18*6dd63270SBarry Smith   #define dmplexgetfulljoin_              DMPLEXGETFULLJOIN
19*6dd63270SBarry Smith   #define dmplexrestorejoin_              DMPLEXRESTOREJOIN
20*6dd63270SBarry Smith   #define dmplexgetmeet_                  DMPLEXGETMEET
21*6dd63270SBarry Smith   #define dmplexgetfullmeet_              DMPLEXGETFULLMEET
22*6dd63270SBarry Smith   #define dmplexrestoremeet_              DMPLEXRESTOREMEET
23*6dd63270SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
24*6dd63270SBarry Smith   #define dmplexgetcone_                  dmplexgetcone
25*6dd63270SBarry Smith   #define dmplexrestorecone_              dmplexrestorecone
26*6dd63270SBarry Smith   #define dmplexgetconeorientation_       dmplexgetconeorientation
27*6dd63270SBarry Smith   #define dmplexrestoreconeorientation_   dmplexrestoreconeorientation
28*6dd63270SBarry Smith   #define dmplexgetsupport_               dmplexgetsupport
29*6dd63270SBarry Smith   #define dmplexrestoresupport_           dmplexrestoresupport
30*6dd63270SBarry Smith   #define dmplexgettransitiveclosure_     dmplexgettransitiveclosure
31*6dd63270SBarry Smith   #define dmplexrestoretransitiveclosure_ dmplexrestoretransitiveclosure
32*6dd63270SBarry Smith   #define dmplexvecgetclosure_            dmplexvecgetclosure
33*6dd63270SBarry Smith   #define dmplexvecrestoreclosure_        dmplexvecrestoreclosure
34*6dd63270SBarry Smith   #define dmplexvecsetclosure_            dmplexvecsetclosure
35*6dd63270SBarry Smith   #define dmplexmatsetclosure_            dmplexmatsetclosure
36*6dd63270SBarry Smith   #define dmplexgetjoin_                  dmplexgetjoin
37*6dd63270SBarry Smith   #define dmplexgetfulljoin_              dmplexgetfulljoin
38*6dd63270SBarry Smith   #define dmplexrestorejoin_              dmplexrestorejoin
39*6dd63270SBarry Smith   #define dmplexgetmeet_                  dmplexgetmeet
40*6dd63270SBarry Smith   #define dmplexgetfullmeet_              dmplexgetfullmeet
41*6dd63270SBarry Smith   #define dmplexrestoremeet_              dmplexrestoremeet
42*6dd63270SBarry Smith #endif
43*6dd63270SBarry Smith 
44*6dd63270SBarry Smith PETSC_EXTERN void dmplexgetcone_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
45*6dd63270SBarry Smith {
46*6dd63270SBarry Smith   const PetscInt *v;
47*6dd63270SBarry Smith   PetscInt        n;
48*6dd63270SBarry Smith 
49*6dd63270SBarry Smith   *ierr = DMPlexGetConeSize(*dm, *p, &n);
50*6dd63270SBarry Smith   if (*ierr) return;
51*6dd63270SBarry Smith   *ierr = DMPlexGetCone(*dm, *p, &v);
52*6dd63270SBarry Smith   if (*ierr) return;
53*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
54*6dd63270SBarry Smith }
55*6dd63270SBarry Smith 
56*6dd63270SBarry Smith PETSC_EXTERN void dmplexrestorecone_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
57*6dd63270SBarry Smith {
58*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
59*6dd63270SBarry Smith   if (*ierr) return;
60*6dd63270SBarry Smith }
61*6dd63270SBarry Smith 
62*6dd63270SBarry Smith PETSC_EXTERN void dmplexgetconeorientation_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
63*6dd63270SBarry Smith {
64*6dd63270SBarry Smith   const PetscInt *v;
65*6dd63270SBarry Smith   PetscInt        n;
66*6dd63270SBarry Smith 
67*6dd63270SBarry Smith   *ierr = DMPlexGetConeSize(*dm, *p, &n);
68*6dd63270SBarry Smith   if (*ierr) return;
69*6dd63270SBarry Smith   *ierr = DMPlexGetConeOrientation(*dm, *p, &v);
70*6dd63270SBarry Smith   if (*ierr) return;
71*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
72*6dd63270SBarry Smith }
73*6dd63270SBarry Smith 
74*6dd63270SBarry Smith PETSC_EXTERN void dmplexrestoreconeorientation_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
75*6dd63270SBarry Smith {
76*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
77*6dd63270SBarry Smith   if (*ierr) return;
78*6dd63270SBarry Smith }
79*6dd63270SBarry Smith 
80*6dd63270SBarry Smith PETSC_EXTERN void dmplexgetsupport_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
81*6dd63270SBarry Smith {
82*6dd63270SBarry Smith   const PetscInt *v;
83*6dd63270SBarry Smith   PetscInt        n;
84*6dd63270SBarry Smith 
85*6dd63270SBarry Smith   *ierr = DMPlexGetSupportSize(*dm, *p, &n);
86*6dd63270SBarry Smith   if (*ierr) return;
87*6dd63270SBarry Smith   *ierr = DMPlexGetSupport(*dm, *p, &v);
88*6dd63270SBarry Smith   if (*ierr) return;
89*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
90*6dd63270SBarry Smith }
91*6dd63270SBarry Smith 
92*6dd63270SBarry Smith PETSC_EXTERN void dmplexrestoresupport_(DM *dm, PetscInt *p, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
93*6dd63270SBarry Smith {
94*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
95*6dd63270SBarry Smith   if (*ierr) return;
96*6dd63270SBarry Smith }
97*6dd63270SBarry Smith 
98*6dd63270SBarry Smith PETSC_EXTERN void dmplexgettransitiveclosure_(DM *dm, PetscInt *p, PetscBool *useCone, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
99*6dd63270SBarry Smith {
100*6dd63270SBarry Smith   PetscInt *v = NULL;
101*6dd63270SBarry Smith   PetscInt  n;
102*6dd63270SBarry Smith 
103*6dd63270SBarry Smith   CHKFORTRANNULL(N);
104*6dd63270SBarry Smith   *ierr = DMPlexGetTransitiveClosure(*dm, *p, *useCone, &n, &v);
105*6dd63270SBarry Smith   if (*ierr) return;
106*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)v, MPIU_INT, 1, n * 2, ptr PETSC_F90_2PTR_PARAM(ptrd));
107*6dd63270SBarry Smith   if (N) *N = n;
108*6dd63270SBarry Smith }
109*6dd63270SBarry Smith 
110*6dd63270SBarry Smith PETSC_EXTERN void dmplexrestoretransitiveclosure_(DM *dm, PetscInt *p, PetscBool *useCone, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
111*6dd63270SBarry Smith {
112*6dd63270SBarry Smith   PetscInt *array;
113*6dd63270SBarry Smith 
114*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_INT, (void **)&array PETSC_F90_2PTR_PARAM(ptrd));
115*6dd63270SBarry Smith   if (*ierr) return;
116*6dd63270SBarry Smith   *ierr = DMPlexRestoreTransitiveClosure(*dm, *p, *useCone, NULL, &array);
117*6dd63270SBarry Smith   if (*ierr) return;
118*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
119*6dd63270SBarry Smith   if (*ierr) return;
120*6dd63270SBarry Smith }
121*6dd63270SBarry Smith 
122*6dd63270SBarry Smith PETSC_EXTERN void dmplexvecgetclosure_(DM *dm, PetscSection *section, Vec *x, PetscInt *point, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
123*6dd63270SBarry Smith {
124*6dd63270SBarry Smith   PetscScalar *v = NULL;
125*6dd63270SBarry Smith   PetscInt     n;
126*6dd63270SBarry Smith 
127*6dd63270SBarry Smith   CHKFORTRANNULL(N);
128*6dd63270SBarry Smith   *ierr = DMPlexVecGetClosure(*dm, *section, *x, *point, &n, &v);
129*6dd63270SBarry Smith   if (*ierr) return;
130*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)v, MPIU_SCALAR, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
131*6dd63270SBarry Smith   if (N) *N = n;
132*6dd63270SBarry Smith }
133*6dd63270SBarry Smith 
134*6dd63270SBarry Smith PETSC_EXTERN void dmplexvecrestoreclosure_(DM *dm, PetscSection *section, Vec *v, PetscInt *point, PetscInt *N, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
135*6dd63270SBarry Smith {
136*6dd63270SBarry Smith   PetscScalar *array;
137*6dd63270SBarry Smith 
138*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&array PETSC_F90_2PTR_PARAM(ptrd));
139*6dd63270SBarry Smith   if (*ierr) return;
140*6dd63270SBarry Smith   *ierr = DMPlexVecRestoreClosure(*dm, *section, *v, *point, NULL, &array);
141*6dd63270SBarry Smith   if (*ierr) return;
142*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
143*6dd63270SBarry Smith   if (*ierr) return;
144*6dd63270SBarry Smith }
145*6dd63270SBarry Smith 
146*6dd63270SBarry Smith PETSC_EXTERN void dmplexgetjoin_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
147*6dd63270SBarry Smith {
148*6dd63270SBarry Smith   const PetscInt *coveredPoints;
149*6dd63270SBarry Smith   PetscInt        n;
150*6dd63270SBarry Smith 
151*6dd63270SBarry Smith   CHKFORTRANNULL(N);
152*6dd63270SBarry Smith   *ierr = DMPlexGetJoin(*dm, *numPoints, points, &n, &coveredPoints);
153*6dd63270SBarry Smith   if (*ierr) return;
154*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
155*6dd63270SBarry Smith   if (N) *N = n;
156*6dd63270SBarry Smith }
157*6dd63270SBarry Smith 
158*6dd63270SBarry Smith PETSC_EXTERN void dmplexgetfulljoin_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
159*6dd63270SBarry Smith {
160*6dd63270SBarry Smith   const PetscInt *coveredPoints;
161*6dd63270SBarry Smith   PetscInt        n;
162*6dd63270SBarry Smith 
163*6dd63270SBarry Smith   CHKFORTRANNULL(N);
164*6dd63270SBarry Smith   *ierr = DMPlexGetFullJoin(*dm, *numPoints, points, &n, &coveredPoints);
165*6dd63270SBarry Smith   if (*ierr) return;
166*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
167*6dd63270SBarry Smith   if (N) *N = n;
168*6dd63270SBarry Smith }
169*6dd63270SBarry Smith 
170*6dd63270SBarry Smith PETSC_EXTERN void dmplexrestorejoin_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
171*6dd63270SBarry Smith {
172*6dd63270SBarry Smith   PetscInt *coveredPoints;
173*6dd63270SBarry Smith 
174*6dd63270SBarry Smith   *ierr = F90Array1dAccess(cptr, MPIU_INT, (void **)&coveredPoints PETSC_F90_2PTR_PARAM(cptrd));
175*6dd63270SBarry Smith   if (*ierr) return;
176*6dd63270SBarry Smith   *ierr = DMPlexRestoreJoin(*dm, 0, NULL, NULL, (const PetscInt **)&coveredPoints);
177*6dd63270SBarry Smith   if (*ierr) return;
178*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(cptr, MPIU_INT PETSC_F90_2PTR_PARAM(cptrd));
179*6dd63270SBarry Smith   if (*ierr) return;
180*6dd63270SBarry Smith }
181*6dd63270SBarry Smith 
182*6dd63270SBarry Smith PETSC_EXTERN void dmplexgetmeet_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
183*6dd63270SBarry Smith {
184*6dd63270SBarry Smith   const PetscInt *coveredPoints;
185*6dd63270SBarry Smith   PetscInt        n;
186*6dd63270SBarry Smith 
187*6dd63270SBarry Smith   CHKFORTRANNULL(N);
188*6dd63270SBarry Smith   *ierr = DMPlexGetMeet(*dm, *numPoints, points, &n, &coveredPoints);
189*6dd63270SBarry Smith   if (*ierr) return;
190*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
191*6dd63270SBarry Smith   if (N) *N = n;
192*6dd63270SBarry Smith }
193*6dd63270SBarry Smith 
194*6dd63270SBarry Smith PETSC_EXTERN void dmplexgetfullmeet_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
195*6dd63270SBarry Smith {
196*6dd63270SBarry Smith   const PetscInt *coveredPoints;
197*6dd63270SBarry Smith   PetscInt        n;
198*6dd63270SBarry Smith 
199*6dd63270SBarry Smith   CHKFORTRANNULL(N);
200*6dd63270SBarry Smith   if (*ierr) return;
201*6dd63270SBarry Smith   *ierr = DMPlexGetFullMeet(*dm, *numPoints, points, &n, &coveredPoints);
202*6dd63270SBarry Smith   if (*ierr) return;
203*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)coveredPoints, MPIU_INT, 1, n, cptr PETSC_F90_2PTR_PARAM(cptrd));
204*6dd63270SBarry Smith   if (N) *N = n;
205*6dd63270SBarry Smith }
206*6dd63270SBarry Smith 
207*6dd63270SBarry Smith PETSC_EXTERN void dmplexrestoremeet_(DM *dm, PetscInt *numPoints, PetscInt *points, PetscInt *N, F90Array1d *cptr, int *ierr PETSC_F90_2PTR_PROTO(cptrd))
208*6dd63270SBarry Smith {
209*6dd63270SBarry Smith   PetscInt *coveredPoints;
210*6dd63270SBarry Smith 
211*6dd63270SBarry Smith   *ierr = F90Array1dAccess(cptr, MPIU_INT, (void **)&coveredPoints PETSC_F90_2PTR_PARAM(cptrd));
212*6dd63270SBarry Smith   if (*ierr) return;
213*6dd63270SBarry Smith   *ierr = DMPlexRestoreMeet(*dm, 0, NULL, NULL, (const PetscInt **)&coveredPoints);
214*6dd63270SBarry Smith   if (*ierr) return;
215*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(cptr, MPIU_INT PETSC_F90_2PTR_PARAM(cptrd));
216*6dd63270SBarry Smith   if (*ierr) return;
217*6dd63270SBarry Smith }
218