xref: /petsc/src/dm/impls/plex/tests/ex25.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[]     = "Test DMCreateLocalVector_Plex, DMPlexGetCellFields and DMPlexRestoreCellFields work properly for 0 fields/cells/DS dimension\n\n";
2 static char FILENAME[] = "ex25.c";
3 
4 #include <petscdmplex.h>
5 #include <petscds.h>
6 #include <petscsnes.h>
7 
8 typedef struct {
9   PetscInt test;
10 } AppCtx;
11 
12 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
13   PetscFunctionBegin;
14   options->test = 0;
15   PetscOptionsBegin(comm, "", "Zero-sized DMPlexGetCellFields Test Options", "DMPLEX");
16   PetscCall(PetscOptionsBoundedInt("-test", "Test to run", FILENAME, options->test, &options->test, NULL, 0));
17   PetscOptionsEnd();
18   PetscFunctionReturn(0);
19 }
20 
21 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *options, DM *dm) {
22   PetscFunctionBegin;
23   PetscCall(DMCreate(comm, dm));
24   PetscCall(DMSetType(*dm, DMPLEX));
25   PetscCall(DMSetFromOptions(*dm));
26   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
27   PetscFunctionReturn(0);
28 }
29 
30 /* no discretization is given so DMGetNumFields yields 0 */
31 static PetscErrorCode test0(DM dm, AppCtx *options) {
32   Vec locX;
33 
34   PetscFunctionBegin;
35   PetscCall(DMGetLocalVector(dm, &locX));
36   PetscCall(DMRestoreLocalVector(dm, &locX));
37   PetscFunctionReturn(0);
38 }
39 
40 /* no discretization is given so DMGetNumFields and PetscDSGetTotalDimension yield 0 */
41 static PetscErrorCode test1(DM dm, AppCtx *options) {
42   IS           cells;
43   Vec          locX, locX_t, locA;
44   PetscScalar *u, *u_t, *a;
45 
46   PetscFunctionBegin;
47   PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &cells));
48   PetscCall(DMGetLocalVector(dm, &locX));
49   PetscCall(DMGetLocalVector(dm, &locX_t));
50   PetscCall(DMGetLocalVector(dm, &locA));
51   PetscCall(DMPlexGetCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
52   PetscCall(DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
53   PetscCall(DMRestoreLocalVector(dm, &locX));
54   PetscCall(DMRestoreLocalVector(dm, &locX_t));
55   PetscCall(DMRestoreLocalVector(dm, &locA));
56   PetscCall(ISDestroy(&cells));
57   PetscFunctionReturn(0);
58 }
59 
60 /* no discretization is given so DMGetNumFields and PetscDSGetTotalDimension yield 0 */
61 static PetscErrorCode test2(DM dm, AppCtx *options) {
62   IS           cells;
63   Vec          locX, locX_t, locA;
64   PetscScalar *u, *u_t, *a;
65   PetscMPIInt  rank;
66 
67   PetscFunctionBegin;
68   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
69   PetscCall(ISCreateStride(PETSC_COMM_SELF, rank ? 0 : 1, 0, 1, &cells));
70   PetscCall(DMGetLocalVector(dm, &locX));
71   PetscCall(DMGetLocalVector(dm, &locX_t));
72   PetscCall(DMGetLocalVector(dm, &locA));
73   PetscCall(DMPlexGetCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
74   PetscCall(DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
75   PetscCall(DMRestoreLocalVector(dm, &locX));
76   PetscCall(DMRestoreLocalVector(dm, &locX_t));
77   PetscCall(DMRestoreLocalVector(dm, &locA));
78   PetscCall(ISDestroy(&cells));
79   PetscFunctionReturn(0);
80 }
81 
82 static PetscErrorCode test3(DM dm, AppCtx *options) {
83   PetscDS   ds;
84   PetscFE   fe;
85   PetscInt  dim;
86   PetscBool simplex;
87 
88   PetscFunctionBegin;
89   PetscCall(DMGetDimension(dm, &dim));
90   PetscCall(DMPlexIsSimplex(dm, &simplex));
91   PetscCall(DMGetDS(dm, &ds));
92   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe));
93   PetscCall(PetscDSSetDiscretization(ds, 0, (PetscObject)fe));
94   PetscCall(PetscFEDestroy(&fe));
95   PetscCall(test1(dm, options));
96   PetscFunctionReturn(0);
97 }
98 
99 static PetscErrorCode test4(DM dm, AppCtx *options) {
100   PetscFE   fe;
101   PetscInt  dim;
102   PetscBool simplex;
103 
104   PetscFunctionBegin;
105   PetscCall(DMGetDimension(dm, &dim));
106   PetscCall(DMPlexIsSimplex(dm, &simplex));
107   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe));
108   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
109   PetscCall(PetscFEDestroy(&fe));
110   PetscCall(DMCreateDS(dm));
111   PetscCall(test2(dm, options));
112   PetscFunctionReturn(0);
113 }
114 
115 static PetscErrorCode test5(DM dm, AppCtx *options) {
116   IS           cells;
117   Vec          locX, locX_t, locA;
118   PetscScalar *u, *u_t, *a;
119 
120   PetscFunctionBegin;
121   locX_t = NULL;
122   locA   = NULL;
123   PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &cells));
124   PetscCall(DMGetLocalVector(dm, &locX));
125   PetscCall(DMPlexGetCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
126   PetscCall(DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
127   PetscCall(DMRestoreLocalVector(dm, &locX));
128   PetscCall(ISDestroy(&cells));
129   PetscFunctionReturn(0);
130 }
131 
132 static PetscErrorCode test6(DM dm, AppCtx *options) {
133   IS           cells;
134   Vec          locX, locX_t, locA;
135   PetscScalar *u, *u_t, *a;
136   PetscMPIInt  rank;
137 
138   PetscFunctionBegin;
139   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
140   locX_t = NULL;
141   locA   = NULL;
142   PetscCall(ISCreateStride(PETSC_COMM_SELF, rank ? 0 : 1, 0, 1, &cells));
143   PetscCall(DMGetLocalVector(dm, &locX));
144   PetscCall(DMPlexGetCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
145   PetscCall(DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a));
146   PetscCall(DMRestoreLocalVector(dm, &locX));
147   PetscCall(ISDestroy(&cells));
148   PetscFunctionReturn(0);
149 }
150 
151 static PetscErrorCode test7(DM dm, AppCtx *options) {
152   PetscFE   fe;
153   PetscInt  dim;
154   PetscBool simplex;
155 
156   PetscFunctionBegin;
157   PetscCall(DMGetDimension(dm, &dim));
158   PetscCall(DMPlexIsSimplex(dm, &simplex));
159   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe));
160   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
161   PetscCall(PetscFEDestroy(&fe));
162   PetscCall(DMCreateDS(dm));
163   PetscCall(test5(dm, options));
164   PetscFunctionReturn(0);
165 }
166 
167 static PetscErrorCode test8(DM dm, AppCtx *options) {
168   PetscFE   fe;
169   PetscInt  dim;
170   PetscBool simplex;
171 
172   PetscFunctionBegin;
173   PetscCall(DMGetDimension(dm, &dim));
174   PetscCall(DMPlexIsSimplex(dm, &simplex));
175   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe));
176   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
177   PetscCall(PetscFEDestroy(&fe));
178   PetscCall(DMCreateDS(dm));
179   PetscCall(test6(dm, options));
180   PetscFunctionReturn(0);
181 }
182 
183 int main(int argc, char **argv) {
184   MPI_Comm comm;
185   DM       dm;
186   AppCtx   options;
187 
188   PetscFunctionBeginUser;
189   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
190   comm = PETSC_COMM_WORLD;
191   PetscCall(ProcessOptions(comm, &options));
192   PetscCall(CreateMesh(comm, &options, &dm));
193 
194   switch (options.test) {
195   case 0: PetscCall(test0(dm, &options)); break;
196   case 1: PetscCall(test1(dm, &options)); break;
197   case 2: PetscCall(test2(dm, &options)); break;
198   case 3: PetscCall(test3(dm, &options)); break;
199   case 4: PetscCall(test4(dm, &options)); break;
200   case 5: PetscCall(test5(dm, &options)); break;
201   case 6: PetscCall(test6(dm, &options)); break;
202   case 7: PetscCall(test7(dm, &options)); break;
203   case 8: PetscCall(test8(dm, &options)); break;
204   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No such test: %" PetscInt_FMT, options.test);
205   }
206 
207   PetscCall(DMDestroy(&dm));
208   PetscCall(PetscFinalize());
209   return 0;
210 }
211 
212 /*TEST
213 
214   testset:
215     args: -dm_plex_dim 3 -dm_plex_interpolate 0
216 
217     test:
218       suffix: 0
219       requires: ctetgen
220       args: -test 0
221     test:
222       suffix: 1
223       requires: ctetgen
224       args: -test 1
225     test:
226       suffix: 2
227       requires: ctetgen
228       args: -test 2
229     test:
230       suffix: 3
231       requires: ctetgen
232       args: -test 3
233     test:
234       suffix: 4
235       requires: ctetgen
236       args: -test 4
237     test:
238       suffix: 5
239       requires: ctetgen
240       args: -test 5
241     test:
242       suffix: 6
243       requires: ctetgen
244       args: -test 6
245     test:
246       suffix: 7
247       requires: ctetgen
248       args: -test 7
249     test:
250       suffix: 8
251       requires: ctetgen
252       args: -test 8
253     test:
254       suffix: 9
255       requires: ctetgen
256       nsize: 2
257       args: -test 1
258     test:
259       suffix: 10
260       requires: ctetgen
261       nsize: 2
262       args: -test 2
263     test:
264       suffix: 11
265       requires: ctetgen
266       nsize: 2
267       args: -test 3
268     test:
269       suffix: 12
270       requires: ctetgen
271       nsize: 2
272       args: -test 4
273 
274 TEST*/
275