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