xref: /petsc/src/dm/impls/plex/tests/ex69.c (revision 98d129c30f3ee9fdddc40fdbc5a989b7be64f888)
1 static char help[] = "Tests for creation of cohesive meshes by transforms\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscsf.h>
5 
6 #include <petsc/private/dmpleximpl.h>
7 
8 PETSC_EXTERN char tri_2_cv[];
9 char              tri_2_cv[] = "\
10 2 4 6 3 1\n\
11 0 2 1\n\
12 1 2 3\n\
13 4 1 5\n\
14 4 0 1\n\
15 -1.0  0.0 0.0  1\n\
16  0.0  1.0 0.0 -1\n\
17  0.0 -1.0 0.0  1\n\
18  1.0  0.0 0.0 -1\n\
19 -2.0  1.0 0.0  1\n\
20 -1.0  2.0 0.0 -1";
21 
22 /* List of test meshes
23 
24 Test tri_0: triangle
25 
26  4-10--5      8-16--7-14--4
27  |\  1 |      |\     \  1 |
28  | \   |      | \     \   |
29  6  8  9  ->  9 12  2  11 13
30  |   \ |      |   \     \ |
31  | 0  \|      | 0  \     \|
32  2--7--3      3-10--6-15--5
33 
34 Test tri_1: triangle, not tensor
35 
36  4-10--5      8-10--7-16--4
37  |\  1 |      |\     \  1 |
38  | \   |      | \     \   |
39  6  8  9  -> 11 14  2  13 15
40  |   \ |      |   \     \ |
41  | 0  \|      | 0  \     \|
42  2--7--3      3-12--6--9--5
43 
44 Test tri_2: 4 triangles, non-oriented surface
45 
46            9
47           / \
48          /   \
49        17  2  16
50        /       \
51       /         \
52      8-----15----5
53       \         /|\
54        \       / | \
55        18  3  12 |  14
56          \   /   |   \
57           \ /    |    \
58            4  0 11  1  7
59             \    |    /
60              \   |   /
61              10  |  13
62                \ | /
63                 \|/
64                  6
65   becomes
66            8
67           / \
68          /   \
69         /     \
70       25   2   24
71       /         \
72      /           \
73    13-----18------9
74 28  |     5    26/ \
75    14----19----10   \
76      \         /|   |\
77       \       / |   | \
78       21  3  20 |   |  23
79         \   /   |   |   \
80          \ /    |   |    \
81           6  0 17 4 16 1  7
82            \    |   |    /
83             \   |   |   /
84             15  |   |  22
85               \ |   | /
86                \|   |/
87                12---11
88                  27
89 
90 Test tri_3: tri_2, in parallel
91 
92            6
93           / \
94          /   \
95         /     \
96       12   1   11
97       /         \
98      /           \
99     5-----10------2
100                    \
101     5-----9-----3   2
102      \         /|   |\
103       \       / |   | \
104       10  1  8  |   |  9
105         \   /   |   |   \
106          \ /    |   |    \
107           2  0  7   7  0  4
108            \    |   |    /
109             \   |   |   /
110              6  |   |  8
111               \ |   | /
112                \|   |/
113                 4   3
114   becomes
115                  11
116                 / \
117                /   \
118               /     \
119             19   1   18
120             /         \
121            /           \
122           8-----14------4
123         22 \     3       |
124             9------15    |\
125                     \    | \
126     9------14-----5  \  20 |
127   20\    3     18/ \  \/   |
128    10----15-----6   |  5   |
129      \         /|   |  |   |\
130       \       / |   |  |   | \
131       17  1 16  |   |  |   |  17
132         \   /   | 2 |  | 2 |   \
133          \ /    |   |  |   |    \
134           4  0  13 12  13  12 0 10
135            \    |   |  |   |    /
136             \   |   |  |   |   /
137             11  |   |  |   |  16
138               \ |   |  |   | /
139                \|   |  |   |/
140                 8---7  7---6
141                  19      21
142 
143 Test quad_0: quadrilateral
144 
145  5-10--6-11--7       5-12-10-20--9-14--6
146  |     |     |       |     |     |     |
147 12  0 13  1  14 --> 15  0 18  2 17  1  16
148  |     |     |       |     |     |     |
149  2--8--3--9--4       3-11--8-19--7-13--4
150 
151 Test quad_1: quadrilateral, not tensor
152 
153  5-10--6-11--7       5-14-10-12--9-16--6
154  |     |     |       |     |     |     |
155 12  0 13  1  14 --> 17  0 20  2 19  1  18
156  |     |     |       |     |     |     |
157  2--8--3--9--4       3-13--8-11--7-15--4
158 
159 Test quad_2: quadrilateral, 2 processes
160 
161  3--6--4  3--6--4       3--9--7-14--6   5-14--4--9--7
162  |     |  |     |       |     |     |   |     |     |
163  7  0  8  7  0  8  --> 10  0 12  1 11  12  1 11  0  10
164  |     |  |     |       |     |     |   |     |     |
165  1--5--2  1--5--2       2--8--5-13--4   3-13--2--8--6
166 
167 Test quad_3: quadrilateral, 4 processes, non-oriented surface
168 
169  3--6--4  3--6--4      3--9--7-14--6   5-14--4--9--7
170  |     |  |     |      |     |     |   |     |     |
171  7  0  8  7  0  8     10  0  12 1  11 12  1 11  0  10
172  |     |  |     |      |     |     |   |     |     |
173  1--5--2  1--5--2      2--8--5-13--4   3-13--2--8--6
174                    -->
175  3--6--4  3--6--4      3--9--7-14--6   5-14--4--9--7
176  |     |  |     |      |     |     |   |     |     |
177  7  0  8  7  0  8     10  0  12 1  11 12  1 11  0  10
178  |     |  |     |      |     |     |   |     |     |
179  1--5--2  1--5--2      2--8--5-13--4   3-13--2--8--6
180 
181 Test quad_4: embedded fault
182 
183 14-24-15-25-16-26--17
184  |     |     |     |
185 28  3 30  4 32  5  34
186  |     |     |     |
187 10-21-11-22-12-23--13
188  |     |     |     |
189 27  0 29  1 31  2  33
190  |     |     |     |
191  6-18--7-19--8-20--9
192 
193 becomes
194 
195  13-26-14-27-15-28--16
196   |     |     |     |
197  30  3 32  4 39  5  40
198   |     |     |     |
199  12-25-17-36-19-38--21
200         |     |     |
201        41  6 42  7  43
202         |     |     |
203  12-25-17-35-18-37--20
204   |     |     |     |
205  29  0 31  1 33  2  34
206   |     |     |     |
207   8-22--9-23-10-24--11
208 
209 Test quad_5: two faults
210 
211 14-24-15-25-16-26--17
212  |     |     |     |
213 28  3 30  4 32  5  34
214  |     |     |     |
215 10-21-11-22-12-23--13
216  |     |     |     |
217 27  0 29  1 31  2  33
218  |     |     |     |
219  6-18--7-19--8-20--9
220 
221 becomes
222 
223 12-26-13-27-14-28--15
224  |     |     |     |
225 37  4 31  3 33  5  40
226  |     |     |     |
227 17-36-18-25-19-39--21
228  |     |     |     |
229 43  6  44   41  7  42
230  |     |     |     |
231 16-35-18-25-19-38--20
232  |     |     |     |
233 29  0 30  1 32  2  34
234  |     |     |     |
235  8-22--9-23-10-24--11
236 
237 Test quad_6: T-junction
238 
239 14-24-15-25-16-26--17
240  |     |     |     |
241 28  3 30  4 32  5  34
242  |     |     |     |
243 10-21-11-22-12-23--13
244  |     |     |     |
245 27  0 29  1 31  2  33
246  |     |     |     |
247  6-18--7-19--8-20--9
248 
249 becomes
250 
251  13-26-14-27-15-28--16
252   |     |     |     |
253  30  3 32  4 39  5  40
254   |     |     |     |
255  12-25-17-36-19-38--21
256         |     |     |
257        41  6 42  7  43
258         |     |     |
259  12-25-17-35-18-37--20
260   |     |     |     |
261  29  0 31  1 33  2  34
262   |     |     |     |
263   8-22--9-23-10-24--11
264 
265 becomes
266 
267  14-28-15-41-21-44--20-29-16
268   |     |     |     |     |
269  31  3 33  5 43  8 42  4  40
270   |     |     |     |     |
271  13-27-17-37-23-46--23-39-19
272         |     |     |     |
273        47  6 48    48  7  49
274         |     |     |     |
275  13-27-17-36-22-45--22-38-18
276   |     |     |     |     |
277  30  0 32  1 34    34  2  35
278   |     |     |     |     |
279   9-24-10-25-11-----11-26-12
280 
281 List of future tests:
282 - Detect and error on intersecting faults
283 - 3D
284 */
285 
286 typedef struct {
287   PetscInt testNum; // The mesh to test
288 } AppCtx;
289 
290 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
291 {
292   PetscFunctionBegin;
293   options->testNum = 0;
294 
295   PetscOptionsBegin(comm, "", "Cohesive Meshing Options", "DMPLEX");
296   PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex5.c", options->testNum, &options->testNum, NULL, 0));
297   PetscOptionsEnd();
298   PetscFunctionReturn(PETSC_SUCCESS);
299 }
300 
301 static PetscErrorCode CreateQuadMesh1(MPI_Comm comm, AppCtx *user, DM *dm)
302 {
303   const PetscInt faces[2] = {1, 1};
304   PetscReal      lower[2], upper[2];
305   DMLabel        label;
306   PetscMPIInt    rank;
307   void          *get_tmp;
308   PetscInt64    *cidx;
309   PetscMPIInt    flg;
310 
311   PetscFunctionBeginUser;
312   PetscCallMPI(MPI_Comm_rank(comm, &rank));
313   // Create serial mesh
314   lower[0] = (PetscReal)(rank % 2);
315   lower[1] = (PetscReal)(rank / 2);
316   upper[0] = (PetscReal)(rank % 2) + 1.;
317   upper[1] = (PetscReal)(rank / 2) + 1.;
318   PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_SELF, 2, PETSC_FALSE, faces, lower, upper, NULL, PETSC_TRUE, dm));
319   PetscCall(PetscObjectSetName((PetscObject)*dm, "box"));
320   // Flip edges to make fault non-oriented
321   switch (rank) {
322   case 2:
323     PetscCall(DMPlexOrientPoint(*dm, 8, -1));
324     break;
325   case 3:
326     PetscCall(DMPlexOrientPoint(*dm, 7, -1));
327     break;
328   default:
329     break;
330   }
331   // Need this so that all procs create the cell types
332   PetscCall(DMPlexGetCellTypeLabel(*dm, &label));
333   // Replace comm in object (copied from PetscHeaderCreate/Destroy())
334   PetscCall(PetscCommDestroy(&(*dm)->hdr.comm));
335   PetscCall(PetscCommDuplicate(comm, &(*dm)->hdr.comm, &(*dm)->hdr.tag));
336   PetscCallMPI(MPI_Comm_get_attr((*dm)->hdr.comm, Petsc_CreationIdx_keyval, &get_tmp, &flg));
337   PetscCheck(flg, (*dm)->hdr.comm, PETSC_ERR_ARG_CORRUPT, "MPI_Comm does not have an object creation index");
338   cidx            = (PetscInt64 *)get_tmp;
339   (*dm)->hdr.cidx = (*cidx)++;
340   // Create new pointSF
341   {
342     PetscSF      sf;
343     PetscInt    *local  = NULL;
344     PetscSFNode *remote = NULL;
345     PetscInt     Nl;
346 
347     PetscCall(PetscSFCreate(comm, &sf));
348     switch (rank) {
349     case 0:
350       Nl = 5;
351       PetscCall(PetscMalloc1(Nl, &local));
352       PetscCall(PetscMalloc1(Nl, &remote));
353       local[0]        = 2;
354       remote[0].index = 1;
355       remote[0].rank  = 1;
356       local[1]        = 3;
357       remote[1].index = 1;
358       remote[1].rank  = 2;
359       local[2]        = 4;
360       remote[2].index = 1;
361       remote[2].rank  = 3;
362       local[3]        = 6;
363       remote[3].index = 5;
364       remote[3].rank  = 2;
365       local[4]        = 8;
366       remote[4].index = 7;
367       remote[4].rank  = 1;
368       break;
369     case 1:
370       Nl = 3;
371       PetscCall(PetscMalloc1(Nl, &local));
372       PetscCall(PetscMalloc1(Nl, &remote));
373       local[0]        = 3;
374       remote[0].index = 1;
375       remote[0].rank  = 3;
376       local[1]        = 4;
377       remote[1].index = 2;
378       remote[1].rank  = 3;
379       local[2]        = 6;
380       remote[2].index = 5;
381       remote[2].rank  = 3;
382       break;
383     case 2:
384       Nl = 3;
385       PetscCall(PetscMalloc1(Nl, &local));
386       PetscCall(PetscMalloc1(Nl, &remote));
387       local[0]        = 2;
388       remote[0].index = 1;
389       remote[0].rank  = 3;
390       local[1]        = 4;
391       remote[1].index = 3;
392       remote[1].rank  = 3;
393       local[2]        = 8;
394       remote[2].index = 7;
395       remote[2].rank  = 3;
396       break;
397     case 3:
398       Nl = 0;
399       break;
400     default:
401       SETERRQ(comm, PETSC_ERR_SUP, "This example only supports 4 ranks");
402     }
403     PetscCall(PetscSFSetGraph(sf, 9, Nl, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
404     PetscCall(DMSetPointSF(*dm, sf));
405     PetscCall(PetscSFDestroy(&sf));
406   }
407   // Create fault label
408   PetscCall(DMCreateLabel(*dm, "fault"));
409   PetscCall(DMGetLabel(*dm, "fault", &label));
410   switch (rank) {
411   case 0:
412   case 2:
413     PetscCall(DMLabelSetValue(label, 8, 1));
414     PetscCall(DMLabelSetValue(label, 2, 0));
415     PetscCall(DMLabelSetValue(label, 4, 0));
416     break;
417   case 1:
418   case 3:
419     PetscCall(DMLabelSetValue(label, 7, 1));
420     PetscCall(DMLabelSetValue(label, 1, 0));
421     PetscCall(DMLabelSetValue(label, 3, 0));
422     break;
423   default:
424     break;
425   }
426   PetscCall(DMPlexOrientLabel(*dm, label));
427   PetscCall(DMPlexLabelCohesiveComplete(*dm, label, NULL, 1, PETSC_FALSE, PETSC_FALSE, NULL));
428   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
429   PetscFunctionReturn(PETSC_SUCCESS);
430 }
431 
432 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
433 {
434   PetscFunctionBegin;
435   if (user->testNum) {
436     PetscCall(CreateQuadMesh1(comm, user, dm));
437   } else {
438     PetscCall(DMCreate(comm, dm));
439     PetscCall(DMSetType(*dm, DMPLEX));
440   }
441   PetscCall(DMSetFromOptions(*dm));
442   {
443     const char *prefix;
444 
445     // We cannot redistribute with cohesive cells in the SF
446     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
447     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)*dm, &prefix));
448     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "f0_"));
449     PetscCall(DMSetFromOptions(*dm));
450     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "f1_"));
451     PetscCall(DMSetFromOptions(*dm));
452     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, prefix));
453   }
454   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
455   PetscFunctionReturn(PETSC_SUCCESS);
456 }
457 
458 int main(int argc, char **argv)
459 {
460   DM     dm;
461   AppCtx user;
462 
463   PetscFunctionBeginUser;
464   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
465   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
466   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
467   PetscCall(DMDestroy(&dm));
468   PetscCall(PetscFinalize());
469   return 0;
470 }
471 
472 /*TEST
473 
474   testset:
475     requires: triangle
476     args: -dm_refine 1 -dm_plex_transform_type cohesive_extrude \
477             -dm_plex_transform_active fault \
478           -dm_view ::ascii_info_detail -coarse_dm_view ::ascii_info_detail
479 
480     test:
481       suffix: tri_0
482       args: -dm_plex_box_faces 1,1 -dm_plex_cohesive_label_fault 8
483     test:
484       suffix: tri_1
485       args: -dm_plex_box_faces 1,1 -dm_plex_cohesive_label_fault 8 \
486               -dm_plex_transform_extrude_use_tensor 0
487     test:
488       suffix: tri_2
489       args: -dm_plex_file_contents dat:tri_2_cv -dm_plex_cohesive_label_fault 11,15
490     test:
491       suffix: tri_3
492       nsize: 2
493       args: -dm_plex_file_contents dat:tri_2_cv -dm_plex_cohesive_label_fault 11,15 \
494               -petscpartitioner_type shell -petscpartitioner_shell_sizes 2,2 \
495               -petscpartitioner_shell_points 0,3,1,2
496 
497   testset:
498     args: -dm_plex_simplex 0 -dm_plex_box_faces 2,1 \
499           -dm_refine 1 -dm_plex_transform_type cohesive_extrude \
500             -dm_plex_transform_active fault -dm_plex_cohesive_label_fault 13 \
501           -dm_view ::ascii_info_detail -coarse_dm_view ::ascii_info_detail
502 
503     test:
504       suffix: quad_0
505     test:
506       suffix: quad_1
507       args: -dm_plex_transform_extrude_use_tensor 0
508     test:
509       suffix: quad_2
510       nsize: 2
511       args: -petscpartitioner_type simple
512 
513   test:
514     suffix: quad_3
515     nsize: 4
516     args: -test_num 1 \
517           -dm_refine 1 -dm_plex_transform_type cohesive_extrude \
518             -dm_plex_transform_active fault \
519           -dm_view ::ascii_info_detail -coarse_dm_view ::ascii_info_detail \
520           -orientation_view -orientation_view_synchronized
521 
522   test:
523     suffix: quad_4
524     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,2 \
525           -dm_refine 1 -dm_plex_transform_type cohesive_extrude \
526             -dm_plex_transform_active fault -dm_plex_cohesive_label_fault 22,23 \
527           -dm_view ::ascii_info_detail -coarse_dm_view ::ascii_info_detail
528 
529   test:
530     suffix: quad_5
531     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,2 \
532             -dm_plex_cohesive_label_fault0 21 \
533             -dm_plex_cohesive_label_fault1 23 \
534           -f0_dm_refine 1 -f0_dm_plex_transform_type cohesive_extrude \
535             -f0_dm_plex_transform_active fault0  -f0_coarse_dm_view ::ascii_info_detail \
536           -f1_dm_refine 1 -f1_dm_plex_transform_type cohesive_extrude \
537             -f1_dm_plex_transform_active fault1  -f1_coarse_dm_view ::ascii_info_detail \
538           -dm_view ::ascii_info_detail
539 
540   test:
541     suffix: quad_6
542     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,2 \
543             -dm_plex_cohesive_label_fault0 22,23 \
544             -dm_plex_cohesive_label_fault1 32 \
545           -f0_dm_refine 1 -f0_dm_plex_transform_type cohesive_extrude \
546             -f0_dm_plex_transform_active fault0  -f0_coarse_dm_view ::ascii_info_detail \
547           -f1_dm_refine 1 -f1_dm_plex_transform_type cohesive_extrude \
548             -f1_dm_plex_transform_active fault1  -f1_coarse_dm_view ::ascii_info_detail \
549           -dm_view ::ascii_info_detail
550 
551 TEST*/
552