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