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