xref: /petsc/src/dm/impls/plex/tests/ex5.c (revision caf9e14d8ca99dc20b2119107fcf65b01f7b29cf)
1c4762a1bSJed Brown static char help[] = "Tests for creation of hybrid meshes\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /* TODO
4c4762a1bSJed Brown  - Propagate hybridSize with distribution
5c4762a1bSJed Brown  - Test with multiple fault segments
6c4762a1bSJed Brown  - Test with embedded fault
7c4762a1bSJed Brown  - Test with multiple faults
8c4762a1bSJed Brown  - Move over all PyLith tests?
9c4762a1bSJed Brown */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown #include <petscdmplex.h>
12b253942bSMatthew G. Knepley #include <petscds.h>
13b253942bSMatthew G. Knepley #include <petsc/private/dmpleximpl.h>
14c4762a1bSJed Brown /* List of test meshes
15c4762a1bSJed Brown 
16c4762a1bSJed Brown Triangle
17c4762a1bSJed Brown --------
18c4762a1bSJed Brown Test 0:
19c4762a1bSJed Brown Two triangles sharing a face
20c4762a1bSJed Brown 
21c4762a1bSJed Brown         4
22c4762a1bSJed Brown       / | \
23c4762a1bSJed Brown      8  |  9
24c4762a1bSJed Brown     /   |   \
25c4762a1bSJed Brown    2  0 7 1  5
26c4762a1bSJed Brown     \   |   /
27c4762a1bSJed Brown      6  |  10
28c4762a1bSJed Brown       \ | /
29c4762a1bSJed Brown         3
30c4762a1bSJed Brown 
31c4762a1bSJed Brown should become two triangles separated by a zero-volume cell with 4 vertices
32c4762a1bSJed Brown 
33c4762a1bSJed Brown         5--16--8              4--12--6 3
34c4762a1bSJed Brown       / |      | \          / |      | | \
35c4762a1bSJed Brown     11  |      |  12       9  |      | |  4
36c4762a1bSJed Brown     /   |      |   \      /   |      | |   \
37c4762a1bSJed Brown    3  0 10  2 14 1  6    2  0 8  1  10 6 0  1
38c4762a1bSJed Brown     \   |      |   /      \   |      | |   /
39c4762a1bSJed Brown      9  |      |  13       7  |      | |  5
40c4762a1bSJed Brown       \ |      | /          \ |      | | /
41c4762a1bSJed Brown         4--15--7              3--11--5 2
42c4762a1bSJed Brown 
43c4762a1bSJed Brown Test 1:
44c4762a1bSJed Brown Four triangles sharing two faces which are oriented against each other
45c4762a1bSJed Brown 
46c4762a1bSJed Brown           9
47c4762a1bSJed Brown          / \
48c4762a1bSJed Brown         /   \
49c4762a1bSJed Brown       17  2  16
50c4762a1bSJed Brown       /       \
51c4762a1bSJed Brown      /         \
52c4762a1bSJed Brown     8-----15----5
53c4762a1bSJed Brown      \         /|\
54c4762a1bSJed Brown       \       / | \
55c4762a1bSJed Brown       18  3  12 |  14
56c4762a1bSJed Brown         \   /   |   \
57c4762a1bSJed Brown          \ /    |    \
58c4762a1bSJed Brown           4  0 11  1  7
59c4762a1bSJed Brown            \    |    /
60c4762a1bSJed Brown             \   |   /
61c4762a1bSJed Brown             10  |  13
62c4762a1bSJed Brown               \ | /
63c4762a1bSJed Brown                \|/
64c4762a1bSJed Brown                 6
65c4762a1bSJed Brown 
66c4762a1bSJed Brown Fault mesh
67c4762a1bSJed Brown 
68c4762a1bSJed Brown 0 --> 0
69c4762a1bSJed Brown 1 --> 1
70c4762a1bSJed Brown 2 --> 2
71c4762a1bSJed Brown 3 --> 3
72c4762a1bSJed Brown 4 --> 5
73c4762a1bSJed Brown 5 --> 6
74c4762a1bSJed Brown 6 --> 8
75c4762a1bSJed Brown 7 --> 11
76c4762a1bSJed Brown 8 --> 15
77c4762a1bSJed Brown 
78c4762a1bSJed Brown        2
79c4762a1bSJed Brown        |
80c4762a1bSJed Brown   6----8----4
81c4762a1bSJed Brown        |    |
82c4762a1bSJed Brown        3    |
83c4762a1bSJed Brown           0-7-1
84c4762a1bSJed Brown             |
85c4762a1bSJed Brown             |
86c4762a1bSJed Brown             5
87c4762a1bSJed Brown 
88c4762a1bSJed Brown should become four triangles separated by two zero-volume cells with 4 vertices
89c4762a1bSJed Brown 
90c4762a1bSJed Brown           11
91c4762a1bSJed Brown           / \
92c4762a1bSJed Brown          /   \
93c4762a1bSJed Brown         /     \
94c4762a1bSJed Brown       22   2   21
95c4762a1bSJed Brown       /         \
96c4762a1bSJed Brown      /           \
97c4762a1bSJed Brown    10-----20------7
98c4762a1bSJed Brown 28  |     5    26/ \
99c4762a1bSJed Brown    14----25----12   \
100c4762a1bSJed Brown      \         /|   |\
101c4762a1bSJed Brown       \       / |   | \
102c4762a1bSJed Brown       23  3  17 |   |  19
103c4762a1bSJed Brown         \   /   |   |   \
104c4762a1bSJed Brown          \ /    |   |    \
105c4762a1bSJed Brown           6  0 24 4 16 1  9
106c4762a1bSJed Brown            \    |   |    /
107c4762a1bSJed Brown             \   |   |   /
108c4762a1bSJed Brown             15  |   |  18
109c4762a1bSJed Brown               \ |   | /
110c4762a1bSJed Brown                \|   |/
111c4762a1bSJed Brown                13---8
112c4762a1bSJed Brown                  27
113c4762a1bSJed Brown 
1145b9dfbb6SMatthew G. Knepley Test 2:
1155b9dfbb6SMatthew G. Knepley Six triangles sharing one face
1165b9dfbb6SMatthew G. Knepley 
1175b9dfbb6SMatthew G. Knepley 11-----12------13
1185b9dfbb6SMatthew G. Knepley  |     /|\     |
1195b9dfbb6SMatthew G. Knepley  | 1  / | \ 4  |
1205b9dfbb6SMatthew G. Knepley  |   /  |  \   |
1215b9dfbb6SMatthew G. Knepley  |  /   |   \  |
1225b9dfbb6SMatthew G. Knepley  | /    |    \ |
1235b9dfbb6SMatthew G. Knepley  |/     |     \|
1245b9dfbb6SMatthew G. Knepley  9  2   |   5  10
1255b9dfbb6SMatthew G. Knepley  |\     |     /|
1265b9dfbb6SMatthew G. Knepley  | \    |    / |
1275b9dfbb6SMatthew G. Knepley  |  \   |   /  |
1285b9dfbb6SMatthew G. Knepley  |   \  |  /   |
1295b9dfbb6SMatthew G. Knepley  | 0  \ | / 3  |
1305b9dfbb6SMatthew G. Knepley  |     \|/     |
1315b9dfbb6SMatthew G. Knepley  6------7------8
1325b9dfbb6SMatthew G. Knepley 
1335b9dfbb6SMatthew G. Knepley Test 3:
1345b9dfbb6SMatthew G. Knepley This is Test 2 on two processes. After the fault, we have
1355b9dfbb6SMatthew G. Knepley 
1365b9dfbb6SMatthew G. Knepley  6--12--7    7--20-10--16-8
1375b9dfbb6SMatthew G. Knepley  |     /|    |     |\     |
1385b9dfbb6SMatthew G. Knepley  | 1  / |    |     | \  1 |
1395b9dfbb6SMatthew G. Knepley 13  11  |    |     |  17  15
1405b9dfbb6SMatthew G. Knepley  |  /   |    |     |   \  |
1415b9dfbb6SMatthew G. Knepley  | /    |    |     |    \ |
1425b9dfbb6SMatthew G. Knepley  |/     |    |     |     \|
1435b9dfbb6SMatthew G. Knepley  5   2  14  11  3 18  2   6
1445b9dfbb6SMatthew G. Knepley  |\     |    |     |     /|
1455b9dfbb6SMatthew G. Knepley  | \    |    |     |    / |
1465b9dfbb6SMatthew G. Knepley  |  \   |    |     |   /  |
1475b9dfbb6SMatthew G. Knepley 10   9  |    |     |  14  13
1485b9dfbb6SMatthew G. Knepley  | 0  \ |    |     | /  0 |
1495b9dfbb6SMatthew G. Knepley  |     \|    |     |/     |
1505b9dfbb6SMatthew G. Knepley  3---8--4    4--19-9--12--5
1515b9dfbb6SMatthew G. Knepley 
1525b9dfbb6SMatthew G. Knepley Test 4:
1535b9dfbb6SMatthew G. Knepley This is Test 2 on six processes. After the fault, we have
1545b9dfbb6SMatthew G. Knepley 
155c4762a1bSJed Brown Tetrahedron
156c4762a1bSJed Brown -----------
157c4762a1bSJed Brown Test 0:
158c4762a1bSJed Brown Two tets sharing a face
159c4762a1bSJed Brown 
160c4762a1bSJed Brown  cell   5 _______    cell
161c4762a1bSJed Brown  0    / | \      \       1
162c4762a1bSJed Brown     16  |  18     22
163c4762a1bSJed Brown     /8 19 10\      \
164c4762a1bSJed Brown    2-15-|----4--21--6
165c4762a1bSJed Brown     \  9| 7 /      /
166c4762a1bSJed Brown     14  |  17     20
167c4762a1bSJed Brown       \ | /      /
168c4762a1bSJed Brown         3-------
169c4762a1bSJed Brown 
170c4762a1bSJed Brown should become two tetrahedrons separated by a zero-volume cell with 3 faces/3 edges/6 vertices
171c4762a1bSJed Brown 
172c4762a1bSJed Brown  cell   6 ___36___10______    cell
173c4762a1bSJed Brown  0    / | \        |\      \     1
174c4762a1bSJed Brown     24  |  26      | 32     30
175c4762a1bSJed Brown     /12 27 14\    33  \      \
176c4762a1bSJed Brown    3-23-|----5--35-|---9--29--7
177c4762a1bSJed Brown     \ 13| 11/      |18 /      /
178c4762a1bSJed Brown     22  |  25      | 31     28
179c4762a1bSJed Brown       \ | /        |/      /
180c4762a1bSJed Brown         4----34----8------
181c4762a1bSJed Brown          cell 2
182c4762a1bSJed Brown 
183c4762a1bSJed Brown In parallel,
184c4762a1bSJed Brown 
185c4762a1bSJed Brown  cell   5 ___28____8      4______    cell
186c4762a1bSJed Brown  0    / | \        |\     |\      \     0
187c4762a1bSJed Brown     19  |   21     | 24   | 13  6  11
188c4762a1bSJed Brown     /10 22 12\    25  \   |8 \      \
189c4762a1bSJed Brown    2-18-|----4--27-|---7  14  3--10--1
190c4762a1bSJed Brown     \ 11| 9 /      |13 /  |  /      /
191c4762a1bSJed Brown     17  |  20      | 23   | 12  5  9
192c4762a1bSJed Brown       \ | /        |/     |/      /
193c4762a1bSJed Brown         3----26----6      2------
194c4762a1bSJed Brown          cell 1
195c4762a1bSJed Brown 
196c4762a1bSJed Brown Test 1:
197c4762a1bSJed Brown Four tets sharing two faces
198c4762a1bSJed Brown 
199c4762a1bSJed Brown Cells:    0-3,4-5
200c4762a1bSJed Brown Vertices: 6-15
201c4762a1bSJed Brown Faces:    16-29,30-34
202c4762a1bSJed Brown Edges:    35-52,53-56
203c4762a1bSJed Brown 
204c4762a1bSJed Brown Quadrilateral
205c4762a1bSJed Brown -------------
206c4762a1bSJed Brown Test 0:
207c4762a1bSJed Brown Two quads sharing a face
208c4762a1bSJed Brown 
209c4762a1bSJed Brown    5--10---4--14---7
210c4762a1bSJed Brown    |       |       |
211c4762a1bSJed Brown   11   0   9   1  13
212c4762a1bSJed Brown    |       |       |
213c4762a1bSJed Brown    2---8---3--12---6
214c4762a1bSJed Brown 
215c4762a1bSJed Brown should become two quads separated by a zero-volume cell with 4 vertices
216c4762a1bSJed Brown 
217c4762a1bSJed Brown    6--13---5-20-10--17---8    5--10---4-14--7  4---7---2
218c4762a1bSJed Brown    |       |     |       |    |       |     |  |       |
219c4762a1bSJed Brown   14   0  12  2 18   1  16   11   0   9  1 12  8   0   6
220c4762a1bSJed Brown    |       |     |       |    |       |     |  |       |
221c4762a1bSJed Brown    3--11---4-19--9--15---7    2---8---3-13--6  3---5---1
222c4762a1bSJed Brown 
223c4762a1bSJed Brown Test 1:
224c4762a1bSJed Brown 
225c4762a1bSJed Brown Original mesh with 9 cells,
226c4762a1bSJed Brown 
2275b9dfbb6SMatthew G. Knepley   9-----10-----11-----12
2285b9dfbb6SMatthew G. Knepley   |      |     ||      |
2295b9dfbb6SMatthew G. Knepley   |      |     ||      |
2305b9dfbb6SMatthew G. Knepley   |   0  |  1  ||  2   |
2315b9dfbb6SMatthew G. Knepley   |      |     ||      |
2325b9dfbb6SMatthew G. Knepley  13-----14-----15-----16
2335b9dfbb6SMatthew G. Knepley   |      |     ||      |
2345b9dfbb6SMatthew G. Knepley   |      |     ||      |
2355b9dfbb6SMatthew G. Knepley   |  3   |  4  ||  5   |
2365b9dfbb6SMatthew G. Knepley   |      |     ||      |
2375b9dfbb6SMatthew G. Knepley  17-----18-----19=====20
238c4762a1bSJed Brown   |      |      |      |
239c4762a1bSJed Brown   |      |      |      |
2405b9dfbb6SMatthew G. Knepley   |  6   |  7   |  8   |
241c4762a1bSJed Brown   |      |      |      |
2425b9dfbb6SMatthew G. Knepley  21-----22-----23-----24
243c4762a1bSJed Brown 
244c4762a1bSJed Brown After first fault,
245c4762a1bSJed Brown 
246c4762a1bSJed Brown  12 ----13 ----14-28 ----15
247c4762a1bSJed Brown   |      |      |  |      |
248c4762a1bSJed Brown   |  0   |  1   | 9|  2   |
249c4762a1bSJed Brown   |      |      |  |      |
250c4762a1bSJed Brown   |      |      |  |      |
251c4762a1bSJed Brown  16 ----17 ----18-29 ----19
252c4762a1bSJed Brown   |      |      |  |      |
253c4762a1bSJed Brown   |  3   |  4   |10|  5   |
254c4762a1bSJed Brown   |      |      |  |      |
255c4762a1bSJed Brown   |      |      |  |      |
256c4762a1bSJed Brown  20 ----21-----22-30 ----23
257c4762a1bSJed Brown   |      |      |  \--11- |
258c4762a1bSJed Brown   |  6   |  7   |     8   |
259c4762a1bSJed Brown   |      |      |         |
260c4762a1bSJed Brown   |      |      |         |
261c4762a1bSJed Brown  24 ----25 ----26--------27
262c4762a1bSJed Brown 
263c4762a1bSJed Brown After second fault,
264c4762a1bSJed Brown 
265c4762a1bSJed Brown  14 ----15 ----16-30 ----17
266c4762a1bSJed Brown   |      |      |  |      |
267c4762a1bSJed Brown   |  0   |  1   | 9|  2   |
268c4762a1bSJed Brown   |      |      |  |      |
269c4762a1bSJed Brown   |      |      |  |      |
270c4762a1bSJed Brown  18 ----19 ----20-31 ----21
271c4762a1bSJed Brown   |      |      |  |      |
272c4762a1bSJed Brown   |  3   |  4   |10|  5   |
273c4762a1bSJed Brown   |      |      |  |      |
274c4762a1bSJed Brown   |      |      |  |      |
275c4762a1bSJed Brown  33 ----34-----24-32 ----25
276c4762a1bSJed Brown   |  12  | 13 / |  \-11-- |
277c4762a1bSJed Brown  22 ----23---/  |         |
2785b9dfbb6SMatthew G. Knepley   |      |      |         |
2795b9dfbb6SMatthew G. Knepley   |  6   |   7  |     8   |
280c4762a1bSJed Brown   |      |      |         |
281c4762a1bSJed Brown   |      |      |         |
282c4762a1bSJed Brown  26 ----27 ----28--------29
283c4762a1bSJed Brown 
2845b9dfbb6SMatthew G. Knepley  Test 2:
2855b9dfbb6SMatthew G. Knepley  Two quads sharing a face in parallel
2865b9dfbb6SMatthew G. Knepley 
2875b9dfbb6SMatthew G. Knepley     4---7---3  2---8---4
2885b9dfbb6SMatthew G. Knepley     |       |  |       |
2895b9dfbb6SMatthew G. Knepley     8   0   6  5   0   7
2905b9dfbb6SMatthew G. Knepley     |       |  |       |
2915b9dfbb6SMatthew G. Knepley     1---5---2  1---6---3
2925b9dfbb6SMatthew G. Knepley 
2935b9dfbb6SMatthew G. Knepley  should become two quads separated by a zero-volume cell with 4 vertices
2945b9dfbb6SMatthew G. Knepley 
2955b9dfbb6SMatthew G. Knepley      4---7---3  3-14--7--11---5
2965b9dfbb6SMatthew G. Knepley      |       |  |     |       |
2975b9dfbb6SMatthew G. Knepley      8   0   6  8  1  12  0   10
2985b9dfbb6SMatthew G. Knepley      |       |  |     |       |
2995b9dfbb6SMatthew G. Knepley      1---5---2  2-13--6---9---4
3005b9dfbb6SMatthew G. Knepley 
3015b9dfbb6SMatthew G. Knepley  Test 3:
3025b9dfbb6SMatthew G. Knepley  Like Test 2, but with different partition
3035b9dfbb6SMatthew G. Knepley 
3045b9dfbb6SMatthew G. Knepley      5--10---4-14--7   2---8---4
3055b9dfbb6SMatthew G. Knepley      |       |     |   |       |
3065b9dfbb6SMatthew G. Knepley     11   0   9  1  12  5   0   7
3075b9dfbb6SMatthew G. Knepley      |       |     |   |       |
3085b9dfbb6SMatthew G. Knepley      2---8---3-13--6   1---6---3
3095b9dfbb6SMatthew G. Knepley 
310c4762a1bSJed Brown Hexahedron
311c4762a1bSJed Brown ----------
312c4762a1bSJed Brown Test 0:
313c4762a1bSJed Brown Two hexes sharing a face
314c4762a1bSJed Brown 
315c4762a1bSJed Brown cell   9-----31------8-----42------13 cell
316c4762a1bSJed Brown 0     /|            /|            /|     1
317c4762a1bSJed Brown     32 |   15      30|   21      41|
318c4762a1bSJed Brown     /  |          /  |          /  |
319c4762a1bSJed Brown    6-----29------7-----40------12  |
320c4762a1bSJed Brown    |   |     18  |   |     24  |   |
321c4762a1bSJed Brown    |  36         |  35         |   44
322c4762a1bSJed Brown    |19 |         |17 |         |23 |
323c4762a1bSJed Brown   33   |  16    34   |   22   43   |
324c4762a1bSJed Brown    |   5-----27--|---4-----39--|---11
325c4762a1bSJed Brown    |  /          |  /          |  /
326c4762a1bSJed Brown    | 28   14     | 26    20    | 38
327c4762a1bSJed Brown    |/            |/            |/
328c4762a1bSJed Brown    2-----25------3-----37------10
329c4762a1bSJed Brown 
330c4762a1bSJed Brown should become two hexes separated by a zero-volume cell with 8 vertices
331c4762a1bSJed Brown 
332c4762a1bSJed Brown                          cell 2
333c4762a1bSJed Brown cell  10-----41------9-----62------18----52------14 cell
334c4762a1bSJed Brown 0     /|            /|            /|            /|     1
335c4762a1bSJed Brown     42 |   20      40|  32       56|   26      51|
336c4762a1bSJed Brown     /  |          /  |          /  |          /  |
337c4762a1bSJed Brown    7-----39------8-----61------17--|-50------13  |
338c4762a1bSJed Brown    |   |     23  |   |         |   |     29  |   |
339c4762a1bSJed Brown    |  46         |  45         |   58        |   54
340c4762a1bSJed Brown    |24 |         |22 |         |30 |         |28 |
341c4762a1bSJed Brown   43   |  21    44   |        57   |   27   53   |
342c4762a1bSJed Brown    |   6-----37--|---5-----60--|---16----49--|---12
343c4762a1bSJed Brown    |  /          |  /          |  /          |  /
344c4762a1bSJed Brown    | 38   19     | 36   31     | 55    25    | 48
345c4762a1bSJed Brown    |/            |/            |/            |/
346c4762a1bSJed Brown    3-----35------4-----59------15----47------11
347c4762a1bSJed Brown 
348c4762a1bSJed Brown In parallel,
349c4762a1bSJed Brown 
350c4762a1bSJed Brown                          cell 2
351c4762a1bSJed Brown cell   9-----31------8-----44------13     8----20------4  cell
352c4762a1bSJed Brown 0     /|            /|            /|     /|           /|     1
353c4762a1bSJed Brown     32 |   15      30|  22       38|   24 |  10      19|
354c4762a1bSJed Brown     /  |          /  |          /  |   /  |         /  |
355c4762a1bSJed Brown    6-----29------7-----43------12  |  7----18------3   |
356c4762a1bSJed Brown    |   |     18  |   |         |   |  |   |    13  |   |
357c4762a1bSJed Brown    |  36         |  35         |   40 |  26        |   22
358c4762a1bSJed Brown    |19 |         |17 |         |20 |  |14 |        |12 |
359c4762a1bSJed Brown   33   |  16    34   |        39   |  25  |  11   21   |
360c4762a1bSJed Brown    |   5-----27--|---4-----42--|---11 |   6----17--|---2
361c4762a1bSJed Brown    |  /          |  /          |  /   |  /         |  /
362c4762a1bSJed Brown    | 28   14     | 26   21     | 37   |23     9    | 16
363c4762a1bSJed Brown    |/            |/            |/     |/           |/
364c4762a1bSJed Brown    2-----25------3-----41------10     5----15------1
365c4762a1bSJed Brown 
366c4762a1bSJed Brown Test 1:
367c4762a1bSJed Brown 
368c4762a1bSJed Brown */
369c4762a1bSJed Brown 
370c4762a1bSJed Brown typedef struct {
371c4762a1bSJed Brown   PetscInt  debug;          /* The debugging level */
372c4762a1bSJed Brown   PetscInt  dim;            /* The topological mesh dimension */
373c4762a1bSJed Brown   PetscBool cellSimplex;    /* Use simplices or hexes */
374c4762a1bSJed Brown   PetscBool testPartition;  /* Use a fixed partitioning for testing */
375c4762a1bSJed Brown   PetscInt  testNum;        /* The particular mesh to test */
376b7519becSMatthew G. Knepley   PetscInt  cohesiveFields; /* The number of cohesive fields */
377c4762a1bSJed Brown } AppCtx;
378c4762a1bSJed Brown 
379b253942bSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
380c4762a1bSJed Brown {
381c4762a1bSJed Brown   PetscFunctionBegin;
382c4762a1bSJed Brown   options->debug          = 0;
383c4762a1bSJed Brown   options->dim            = 2;
384c4762a1bSJed Brown   options->cellSimplex    = PETSC_TRUE;
385c4762a1bSJed Brown   options->testPartition  = PETSC_TRUE;
386c4762a1bSJed Brown   options->testNum        = 0;
387b7519becSMatthew G. Knepley   options->cohesiveFields = 1;
388c4762a1bSJed Brown 
389d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
3909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex5.c", options->debug, &options->debug, NULL,0));
3919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex5.c", options->dim, &options->dim, NULL,1,3));
3929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex5.c", options->cellSimplex, &options->cellSimplex, NULL));
3939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex5.c", options->testPartition, &options->testPartition, NULL));
3949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex5.c", options->testNum, &options->testNum, NULL,0));
3959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-cohesive_fields", "The number of cohesive fields", "ex5.c", options->cohesiveFields, &options->cohesiveFields, NULL, 0));
396d0609cedSBarry Smith   PetscOptionsEnd();
397c4762a1bSJed Brown   PetscFunctionReturn(0);
398c4762a1bSJed Brown }
399c4762a1bSJed Brown 
400b253942bSMatthew G. Knepley static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
401c4762a1bSJed Brown {
402c4762a1bSJed Brown   DM             idm;
403c4762a1bSJed Brown   PetscInt       p;
404c4762a1bSJed Brown   PetscMPIInt    rank;
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   PetscFunctionBegin;
4079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
408dd400576SPatrick Sanan   if (rank == 0) {
409c4762a1bSJed Brown     switch (testNum) {
410c4762a1bSJed Brown     case 0:
411c4762a1bSJed Brown     {
412c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
413c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
414c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
415c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
416c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
417c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
418c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
419c4762a1bSJed Brown 
4209566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4219566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
4229566063dSJacob Faibussowitsch       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4239566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
4249566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
425c4762a1bSJed Brown     }
426c4762a1bSJed Brown     break;
427c4762a1bSJed Brown     case 1:
428c4762a1bSJed Brown     {
429c4762a1bSJed Brown       PetscInt    numPoints[2]         = {6, 4};
430c4762a1bSJed Brown       PetscInt    coneSize[10]         = {3, 3, 3, 3, 0, 0, 0, 0, 0, 0};
431c4762a1bSJed Brown       PetscInt    cones[12]            = {4, 6, 5,  5, 6, 7,  8, 5, 9,  8, 4, 5};
432c4762a1bSJed Brown       PetscInt    coneOrientations[12] = {0, 0, 0,  0, 0, 0,  0, 0, 0,  0, 0, 0};
433c4762a1bSJed Brown       PetscScalar vertexCoords[12]     = {-1.0, 0.0, 0.0, 1.0, 0.0, -1.0, 1.0, 0.0, -2.0, 1.0, -1.0, 2.0};
434c4762a1bSJed Brown       PetscInt    markerPoints[6]      = {4, 1, 6, 1, 8, 1};
435c4762a1bSJed Brown       PetscInt    faultPoints[3]       = {5, 6, 8};
436c4762a1bSJed Brown 
4379566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4389566063dSJacob Faibussowitsch       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
4399566063dSJacob Faibussowitsch       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4409566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));PetscCall(DMSetLabelValue(*dm, "material", 3, 1));
4419566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));PetscCall(DMSetLabelValue(*dm, "material", 2, 2));
442c4762a1bSJed Brown     }
443c4762a1bSJed Brown     break;
4445b9dfbb6SMatthew G. Knepley     case 2:
4455b9dfbb6SMatthew G. Knepley     case 3:
4465b9dfbb6SMatthew G. Knepley     case 4:
4475b9dfbb6SMatthew G. Knepley     {
4485b9dfbb6SMatthew G. Knepley       PetscInt    numPoints[2]         = {8, 6};
4495b9dfbb6SMatthew G. Knepley       PetscInt    coneSize[14]         = {3, 3, 3, 3, 3, 3, 0, 0, 0, 0, 0, 0, 0, 0};
4505b9dfbb6SMatthew G. Knepley       PetscInt    cones[18]            = {6, 7,  9,   9, 12, 11,  7, 12,  9,
4515b9dfbb6SMatthew G. Knepley                                           7, 8, 10,  10, 13, 12,  7, 10, 12};
4525b9dfbb6SMatthew G. Knepley       PetscInt    coneOrientations[18] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
4535b9dfbb6SMatthew G. Knepley       PetscScalar vertexCoords[16]     = {-1., -1.,  0., -1.,  1., -1.,  -1., 0.,  1., 0.,
4545b9dfbb6SMatthew G. Knepley                                           -1.,  1.,  0.,  1.,  1.,  1.,};
4555b9dfbb6SMatthew G. Knepley       PetscInt    markerPoints[16]     = {6, 1, 7, 1, 8, 1, 9, 1, 10, 1, 11, 1, 12, 1, 13, 1};
4565b9dfbb6SMatthew G. Knepley       PetscInt    faultPoints[2]       = {7, 12};
4575b9dfbb6SMatthew G. Knepley 
4585b9dfbb6SMatthew G. Knepley       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4595b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
4605b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
4615b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
4625b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
4635b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 4, 2));
4645b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
4655b9dfbb6SMatthew G. Knepley       for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
4665b9dfbb6SMatthew G. Knepley       if (testNum == 2) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4675b9dfbb6SMatthew G. Knepley       if (testNum == 3 || testNum == 4) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1));
4685b9dfbb6SMatthew G. Knepley     }
4695b9dfbb6SMatthew G. Knepley     break;
470c4762a1bSJed Brown     default:
47163a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
472c4762a1bSJed Brown     }
473c4762a1bSJed Brown   } else {
474c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
475c4762a1bSJed Brown 
4769566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
4775b9dfbb6SMatthew G. Knepley     if (testNum == 3 || testNum == 4) PetscCall(DMCreateLabel(*dm, "pfault"));
4785b9dfbb6SMatthew G. Knepley     else                              PetscCall(DMCreateLabel(*dm, "fault"));
479c4762a1bSJed Brown   }
4809566063dSJacob Faibussowitsch   PetscCall(DMPlexInterpolate(*dm, &idm));
4819566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
4829566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
483c4762a1bSJed Brown   *dm  = idm;
484c4762a1bSJed Brown   PetscFunctionReturn(0);
485c4762a1bSJed Brown }
486c4762a1bSJed Brown 
487b253942bSMatthew G. Knepley static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, AppCtx *user, DM dm)
488c4762a1bSJed Brown {
489c4762a1bSJed Brown   PetscInt       depth = 3, testNum  = user->testNum, p;
490c4762a1bSJed Brown   PetscMPIInt    rank;
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   PetscFunctionBegin;
4939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
494dd400576SPatrick Sanan   if (rank == 0) {
495c4762a1bSJed Brown     switch (testNum) {
496c4762a1bSJed Brown     case 0:
497c4762a1bSJed Brown     {
498c4762a1bSJed Brown       PetscInt    numPoints[4]         = {5, 7, 9, 2};
499c4762a1bSJed Brown       PetscInt    coneSize[23]         = {4, 4, 0, 0, 0, 0, 0, 3, 3, 3, 3, 3, 3, 3, 2, 2, 2, 2, 2, 2, 2, 2, 2};
500c4762a1bSJed Brown       PetscInt    cones[47]            = {7, 8, 9, 10,  11, 10, 13, 12,  15, 17, 14,  16, 18, 15,  14, 19, 16,  17, 18, 19,  17, 21, 20,  18, 22, 21,  22, 19, 20,   2, 3,  2, 4,  2, 5,  3, 4,  4, 5,  5, 3,  3, 6,  4, 6,  5, 6};
501b5a892a1SMatthew G. Knepley       PetscInt    coneOrientations[47] = {0, 0, 0,  0,   0, -2,  2,  2,   0, -1, -1,   0, -1, -1,   0, -1, -1,   0,  0,  0,   0,  0, -1,   0,  0, -1,  -1,  0,  0,   0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0};
502c4762a1bSJed Brown       PetscScalar vertexCoords[15]     = {0.0, 0.0, -0.5,  0.0, -0.5, 0.0,  1.0, 0.0, 0.0,  0.0, 0.5, 0.0,  0.0, 0.0, 0.5};
503c4762a1bSJed Brown       PetscInt    markerPoints[20]     = {2, 1, 3, 1, 4, 1, 5, 1, 14, 1, 15, 1, 16, 1, 17, 1, 18, 1, 19, 1};
504c4762a1bSJed Brown       PetscInt    faultPoints[3]      = {3, 4, 5};
505c4762a1bSJed Brown 
5069566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
507c4762a1bSJed Brown       for (p = 0; p < 10; ++p) {
5089566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
509c4762a1bSJed Brown       }
510c4762a1bSJed Brown       for (p = 0; p < 3; ++p) {
5119566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1));
512c4762a1bSJed Brown       }
5139566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 0, 1));
5149566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 1, 2));
515c4762a1bSJed Brown     }
516c4762a1bSJed Brown     break;
517c4762a1bSJed Brown     case 1:
518c4762a1bSJed Brown     {
519c4762a1bSJed Brown       PetscInt    numPoints[4]         = {6, 13, 12, 4};
520c4762a1bSJed Brown       PetscInt    coneSize[35]         = {4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2};
521c4762a1bSJed Brown       PetscInt    cones[78]            = {10, 11, 12, 13,  10, 15, 16, 14,  17, 18, 14, 19,  20, 13, 19, 21,  22, 23, 24,  25, 26, 22,  24, 27, 25,  23, 26, 27,  28, 29, 23,  24, 30, 28,  22, 29, 30,   31, 32, 28,  29, 33, 31,  32, 33, 23,  26, 34, 33,  34, 27, 32,  6, 5,  5, 7,  7, 6,  6, 4,  4, 5,  7, 4,  7, 9,  9, 5,  6, 9,  9, 8,  8, 7,  5, 8,  4, 8};
522b5a892a1SMatthew G. Knepley       PetscInt    coneOrientations[78] = { 0,  0,  0,  0,  -2,  1,  0,  2,   0,  0, -3,  0,   0, -3, -1,  0,   0,  0,  0,   0,  0, -1,  -1,  0, -1,  -1, -1, -1,   0,  0,  0,   0,  0, -1,   0, -1, -1,    0,  0,  0,   0,  0, -1,  -1, -1,  0,  -1,  0, -1,  -1, -1, -1,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0};
523c4762a1bSJed Brown       PetscScalar vertexCoords[18]     = {-1.0, 0.0, 0.0,  0.0, -1.0, 0.0,  0.0, 0.0, 1.0,  0.0, 1.0, 0.0,  0.0, 0.0, -1.0,  1.0, 0.0, 0.0};
524c4762a1bSJed Brown       PetscInt    markerPoints[14]     = {5, 1, 6, 1, 7, 1, 10, 1, 22, 1, 23, 1, 24, 1};
525c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {5, 6, 7, 8};
526c4762a1bSJed Brown 
5279566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
528c4762a1bSJed Brown       for (p = 0; p < 7; ++p) {
5299566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
530c4762a1bSJed Brown       }
531c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
5329566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1));
533c4762a1bSJed Brown       }
5349566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 0, 1));
5359566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 1, 2));
536c4762a1bSJed Brown     }
537c4762a1bSJed Brown     break;
538c4762a1bSJed Brown     default:
53963a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
540c4762a1bSJed Brown     }
541c4762a1bSJed Brown   } else {
542c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
543c4762a1bSJed Brown 
5449566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
5459566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "fault"));
546c4762a1bSJed Brown   }
547c4762a1bSJed Brown   PetscFunctionReturn(0);
548c4762a1bSJed Brown }
549c4762a1bSJed Brown 
550b253942bSMatthew G. Knepley static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
551c4762a1bSJed Brown {
552c4762a1bSJed Brown   DM             idm;
553c4762a1bSJed Brown   PetscInt       p;
554c4762a1bSJed Brown   PetscMPIInt    rank;
555c4762a1bSJed Brown 
556c4762a1bSJed Brown   PetscFunctionBegin;
5579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
558dd400576SPatrick Sanan   if (rank == 0) {
559c4762a1bSJed Brown     switch (testNum) {
560c4762a1bSJed Brown     case 0:
561c4762a1bSJed Brown     case 2:
5625b9dfbb6SMatthew G. Knepley     case 3:
563c4762a1bSJed Brown     {
564c4762a1bSJed Brown       PetscInt    numPoints[2]        = {6, 2};
565c4762a1bSJed Brown       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
566c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
567c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
568c4762a1bSJed Brown       PetscScalar vertexCoords[12]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0, 0.5, 1.0};
569c4762a1bSJed Brown       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
570c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
571c4762a1bSJed Brown 
5729566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
5739566063dSJacob Faibussowitsch       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
5749566063dSJacob Faibussowitsch       if (testNum == 0) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
5755b9dfbb6SMatthew G. Knepley       if (testNum == 2 || testNum == 3) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1));
5769566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
5779566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
578c4762a1bSJed Brown     }
579c4762a1bSJed Brown     break;
580c4762a1bSJed Brown     case 1:
581c4762a1bSJed Brown     {
582c4762a1bSJed Brown       PetscInt    numPoints[2]         = {16, 9};
583c4762a1bSJed Brown       PetscInt    coneSize[25]         = {4, 4, 4, 4, 4, 4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
584c4762a1bSJed Brown       PetscInt    cones[36]            = {9,  13, 14, 10,
585c4762a1bSJed Brown                                           10, 14, 15, 11,
586c4762a1bSJed Brown                                           11, 15, 16, 12,
587c4762a1bSJed Brown                                           13, 17, 18, 14,
588c4762a1bSJed Brown                                           14, 18, 19, 15,
589c4762a1bSJed Brown                                           15, 19, 20, 16,
590c4762a1bSJed Brown                                           17, 21, 22, 18,
591c4762a1bSJed Brown                                           18, 22, 23, 19,
592c4762a1bSJed Brown                                           19, 23, 24, 20};
593c4762a1bSJed Brown       PetscInt    coneOrientations[36] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
594c4762a1bSJed Brown       PetscScalar vertexCoords[32]     = {-3.0,  3.0,  -1.0,  3.0,  1.0,  3.0,  3.0,  3.0,  -3.0,  1.0,  -1.0,  1.0,  1.0,  1.0,  3.0,  1.0,
595c4762a1bSJed Brown                                           -3.0, -1.0,  -1.0, -1.0,  1.0, -1.0,  3.0, -1.0,  -3.0, -3.0,  -1.0, -3.0,  1.0, -3.0,  3.0, -3.0};
5965b9dfbb6SMatthew G. Knepley       PetscInt    faultPoints[4]       = {11, 15, 19, 20};
597c4762a1bSJed Brown       PetscInt    fault2Points[2]      = {17, 18};
598c4762a1bSJed Brown 
5999566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6005b9dfbb6SMatthew G. Knepley       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "fault",  faultPoints[p], 1));
6015b9dfbb6SMatthew G. Knepley       for (p = 3; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "faultBd", faultPoints[p], 1));
6029566063dSJacob Faibussowitsch       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault2", fault2Points[p], 1));
6039566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6049566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
6059566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
6069566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 1));
6079566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 4, 1));
6089566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
6099566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 6, 2));
6109566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 7, 2));
6119566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 8, 2));
612c4762a1bSJed Brown     }
613c4762a1bSJed Brown     break;
614c4762a1bSJed Brown     default:
61563a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
616c4762a1bSJed Brown     }
617c4762a1bSJed Brown   } else {
618c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
619c4762a1bSJed Brown 
6209566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
6215b9dfbb6SMatthew G. Knepley     if (testNum == 2 || testNum == 3) PetscCall(DMCreateLabel(*dm, "pfault"));
6229566063dSJacob Faibussowitsch     else                              PetscCall(DMCreateLabel(*dm, "fault"));
623c4762a1bSJed Brown   }
6249566063dSJacob Faibussowitsch   PetscCall(DMPlexInterpolate(*dm, &idm));
6259566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
6269566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
627c4762a1bSJed Brown   *dm  = idm;
628c4762a1bSJed Brown   PetscFunctionReturn(0);
629c4762a1bSJed Brown }
630c4762a1bSJed Brown 
631b253942bSMatthew G. Knepley static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
632c4762a1bSJed Brown {
633c4762a1bSJed Brown   DM             idm;
634c4762a1bSJed Brown   PetscInt       depth = 3, p;
635c4762a1bSJed Brown   PetscMPIInt    rank;
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   PetscFunctionBegin;
6389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
639dd400576SPatrick Sanan   if (rank == 0) {
640c4762a1bSJed Brown     switch (testNum) {
641c4762a1bSJed Brown     case 0:
642c4762a1bSJed Brown     {
643c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
644c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8,8, 0,0,0,0,0,0,0,0,0,0,0,0};
645c4762a1bSJed Brown       PetscInt    cones[16]            = {2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8};
646c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0,0,0,0,0,0,0,0,  0,0, 0,0 ,0, 0, 0,0};
647c4762a1bSJed Brown       PetscScalar vertexCoords[36]     = {-0.5,0.0,0.0, 0.0,0.0,0.0, 0.0,1.0,0.0, -0.5,1.0,0.0,
648c4762a1bSJed Brown                                           -0.5,0.0,1.0, 0.0,0.0,1.0, 0.0,1.0,1.0, -0.5,1.0,1.0,
649c4762a1bSJed Brown                                            0.5,0.0,0.0, 0.5,1.0,0.0, 0.5,0.0,1.0,  0.5,1.0,1.0};
650c4762a1bSJed Brown       PetscInt    markerPoints[52]     = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
651c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {3, 4, 7, 8};
652c4762a1bSJed Brown 
6539566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6549566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
6559566063dSJacob Faibussowitsch       for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(idm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
6569566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
6579566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6589566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
659c4762a1bSJed Brown     }
660c4762a1bSJed Brown     break;
661c4762a1bSJed Brown     case 1:
662c4762a1bSJed Brown     {
663c4762a1bSJed Brown       /* Cell Adjacency Graph:
664c4762a1bSJed Brown         0 -- { 8, 13, 21, 24} --> 1
665c4762a1bSJed Brown         0 -- {20, 21, 23, 24} --> 5 F
666c4762a1bSJed Brown         1 -- {10, 15, 21, 24} --> 2
667c4762a1bSJed Brown         1 -- {13, 14, 15, 24} --> 6
668c4762a1bSJed Brown         2 -- {21, 22, 24, 25} --> 4 F
669c4762a1bSJed Brown         3 -- {21, 24, 30, 35} --> 4
670c4762a1bSJed Brown         3 -- {21, 24, 28, 33} --> 5
671c4762a1bSJed Brown        */
672c4762a1bSJed Brown       PetscInt    numPoints[2]         = {30, 7};
673c4762a1bSJed Brown       PetscInt    coneSize[37]         = {8,8,8,8,8,8,8, 0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
674c4762a1bSJed Brown       PetscInt    cones[56]            = { 8, 21, 20,  7, 13, 12, 23, 24,
675c4762a1bSJed Brown                                           14, 15, 10,  9, 13,  8, 21, 24,
676c4762a1bSJed Brown                                           15, 16, 11, 10, 24, 21, 22, 25,
677c4762a1bSJed Brown                                           30, 29, 28, 21, 35, 24, 33, 34,
678c4762a1bSJed Brown                                           24, 21, 30, 35, 25, 36, 31, 22,
679c4762a1bSJed Brown                                           27, 20, 21, 28, 32, 33, 24, 23,
680c4762a1bSJed Brown                                           15, 24, 13, 14, 19, 18, 17, 26};
681c4762a1bSJed Brown       PetscInt    coneOrientations[56] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
682c4762a1bSJed Brown       PetscScalar vertexCoords[90]     = {-2.0, -2.0, -2.0,  -2.0, -1.0, -2.0,  -3.0,  0.0, -2.0,  -2.0,  1.0, -2.0,  -2.0,  2.0, -2.0,  -2.0, -2.0,  0.0,
683c4762a1bSJed Brown                                           -2.0, -1.0,  0.0,  -3.0,  0.0,  0.0,  -2.0,  1.0,  0.0,  -2.0,  2.0,  0.0,  -2.0, -1.0,  2.0,  -3.0,  0.0,  2.0,
684c4762a1bSJed Brown                                           -2.0,  1.0,  2.0,   0.0, -2.0, -2.0,   0.0,  0.0, -2.0,   0.0,  2.0, -2.0,   0.0, -2.0,  0.0,   0.0,  0.0,  0.0,
685c4762a1bSJed Brown                                            0.0,  2.0,  0.0,   0.0,  0.0,  2.0,   2.0, -2.0, -2.0,   2.0, -1.0, -2.0,   3.0,  0.0, -2.0,   2.0,  1.0, -2.0,
686c4762a1bSJed Brown                                            2.0,  2.0, -2.0,   2.0, -2.0,  0.0,   2.0, -1.0,  0.0,   3.0,  0.0,  0.0,   2.0,  1.0,  0.0,   2.0,  2.0,  0.0};
687c4762a1bSJed Brown       PetscInt    faultPoints[6]       = {20, 21, 22, 23, 24, 25};
688c4762a1bSJed Brown 
6899566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6909566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
6919566063dSJacob Faibussowitsch       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
6929566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6939566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
6949566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
6959566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
6969566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 4, 2));
6979566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
6989566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 6, 2));
699c4762a1bSJed Brown     }
700c4762a1bSJed Brown     break;
701c4762a1bSJed Brown     case 2:
702c4762a1bSJed Brown     {
703c4762a1bSJed Brown       /* Buried fault edge */
704c4762a1bSJed Brown       PetscInt    numPoints[2]         = {18, 4};
705c4762a1bSJed Brown       PetscInt    coneSize[22]         = {8,8,8,8, 0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
706c4762a1bSJed Brown       PetscInt    cones[32]            = { 4,  5,  8,  7, 13, 16, 17, 14,
707c4762a1bSJed Brown                                            5,  6,  9,  8, 14, 17, 18, 15,
708c4762a1bSJed Brown                                            7,  8, 11, 10, 16, 19, 20, 17,
709c4762a1bSJed Brown                                            8,  9, 12, 11, 17, 20, 21, 18};
710c4762a1bSJed Brown       PetscInt    coneOrientations[32] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
711c4762a1bSJed Brown       PetscScalar vertexCoords[54]     = {-2.0, -2.0,  0.0,  -2.0,  0.0,  0.0,  -2.0,  2.0,  0.0,   0.0, -2.0,  0.0,   0.0,  0.0,  0.0,   0.0,  2.0,  0.0,
712c4762a1bSJed Brown                                            2.0, -2.0,  0.0,   2.0,  0.0,  0.0,   2.0,  2.0,  0.0,  -2.0, -2.0,  2.0,  -2.0,  0.0,  2.0,  -2.0,  2.0,  2.0,
713c4762a1bSJed Brown                                            0.0, -2.0,  2.0,   0.0,  0.0,  2.0,   0.0,  2.0,  2.0,   2.0, -2.0,  2.0,   2.0,  0.0,  2.0,   2.0,  2.0,  2.0};
714c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {7, 8, 16, 17};
715c4762a1bSJed Brown 
7169566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
7179566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
7189566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
7199566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
7209566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
7219566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 2));
7229566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
723c4762a1bSJed Brown     }
724c4762a1bSJed Brown     break;
72563a3b9bcSJacob Faibussowitsch     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
726c4762a1bSJed Brown     }
727c4762a1bSJed Brown   } else {
728c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
729c4762a1bSJed Brown 
7309566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
7319566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
7329566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(idm, "fault"));
733c4762a1bSJed Brown   }
7349566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
7359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
736c4762a1bSJed Brown   *dm  = idm;
737c4762a1bSJed Brown   PetscFunctionReturn(0);
738c4762a1bSJed Brown }
739c4762a1bSJed Brown 
740b253942bSMatthew G. Knepley static PetscErrorCode CreateFaultLabel(DM dm)
741b253942bSMatthew G. Knepley {
742b253942bSMatthew G. Knepley   DMLabel        label;
743b253942bSMatthew G. Knepley   PetscInt       dim, h, pStart, pEnd, pMax, p;
744b253942bSMatthew G. Knepley 
745b253942bSMatthew G. Knepley   PetscFunctionBegin;
7469566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7479566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "cohesive"));
7489566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &label));
749b253942bSMatthew G. Knepley   for (h = 0; h <= dim; ++h) {
7509566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(dm, h, NULL, &pMax));
7519566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &pStart, &pEnd));
7529566063dSJacob Faibussowitsch     for (p = pMax; p < pEnd; ++p) PetscCall(DMLabelSetValue(label, p, 1));
753b253942bSMatthew G. Knepley   }
754b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
755b253942bSMatthew G. Knepley }
756b253942bSMatthew G. Knepley 
757b253942bSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *user)
758b253942bSMatthew G. Knepley {
759b253942bSMatthew G. Knepley   PetscFE        fe;
760b253942bSMatthew G. Knepley   DMLabel        fault;
761b7519becSMatthew G. Knepley   PetscInt       dim, Ncf = user->cohesiveFields, f;
762b253942bSMatthew G. Knepley 
763b253942bSMatthew G. Knepley   PetscFunctionBegin;
7649566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7659566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &fault));
7669566063dSJacob Faibussowitsch   PetscCall(DMLabelView(fault, PETSC_VIEWER_STDOUT_WORLD));
767b253942bSMatthew G. Knepley 
7689566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, user->cellSimplex, "displacement_", PETSC_DETERMINE, &fe));
7699566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(fe, "displacement"));
7709566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject) fe));
7719566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
772b253942bSMatthew G. Knepley 
773b7519becSMatthew G. Knepley   if (Ncf > 0) {
7749566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim-1, dim, user->cellSimplex, "faulttraction_", PETSC_DETERMINE, &fe));
7759566063dSJacob Faibussowitsch     PetscCall(PetscFESetName(fe, "fault traction"));
7769566063dSJacob Faibussowitsch     PetscCall(DMAddField(dm, fault, (PetscObject) fe));
7779566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
778b7519becSMatthew G. Knepley   }
779b7519becSMatthew G. Knepley   for (f = 1; f < Ncf; ++f) {
780b7519becSMatthew G. Knepley     char name[256], opt[256];
781b7519becSMatthew G. Knepley 
78263a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, 256, "fault field %" PetscInt_FMT, f));
78363a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(opt,  256, "faultfield_%" PetscInt_FMT "_", f));
7849566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim-1, dim, user->cellSimplex, opt, PETSC_DETERMINE, &fe));
7859566063dSJacob Faibussowitsch     PetscCall(PetscFESetName(fe, name));
7869566063dSJacob Faibussowitsch     PetscCall(DMAddField(dm, fault, (PetscObject) fe));
7879566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
788b7519becSMatthew G. Knepley   }
789b253942bSMatthew G. Knepley 
7909566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
791b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
792b253942bSMatthew G. Knepley }
793b253942bSMatthew G. Knepley 
794b253942bSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
795c4762a1bSJed Brown {
796c4762a1bSJed Brown   PetscInt       dim          = user->dim;
797c4762a1bSJed Brown   PetscBool      cellSimplex  = user->cellSimplex, hasFault, hasFault2, hasParallelFault;
798c4762a1bSJed Brown   PetscMPIInt    rank, size;
799dd96b1f9SToby Isaac   DMLabel        matLabel;
800c4762a1bSJed Brown 
801c4762a1bSJed Brown   PetscFunctionBegin;
8029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
8049566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
8059566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
8069566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
807c4762a1bSJed Brown   switch (dim) {
808c4762a1bSJed Brown   case 2:
809c4762a1bSJed Brown     if (cellSimplex) {
8109566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_2D(comm, user->testNum, dm));
811c4762a1bSJed Brown     } else {
8129566063dSJacob Faibussowitsch       PetscCall(CreateQuad_2D(comm, user->testNum, dm));
813c4762a1bSJed Brown     }
814c4762a1bSJed Brown     break;
815c4762a1bSJed Brown   case 3:
816c4762a1bSJed Brown     if (cellSimplex) {
8179566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_3D(comm, user, *dm));
818c4762a1bSJed Brown     } else {
8199566063dSJacob Faibussowitsch       PetscCall(CreateHex_3D(comm, user->testNum, dm));
820c4762a1bSJed Brown     }
821c4762a1bSJed Brown     break;
822c4762a1bSJed Brown   default:
82363a3b9bcSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %" PetscInt_FMT, dim);
824c4762a1bSJed Brown   }
8259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "orig_"));
8269566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
8279566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
8289566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dm, "material", &matLabel));
829dd96b1f9SToby Isaac   if (matLabel) {
8309566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelComplete(*dm, matLabel));
831dd96b1f9SToby Isaac   }
8329566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
8339566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "fault", &hasFault));
834c4762a1bSJed Brown   if (hasFault) {
835b253942bSMatthew G. Knepley     DM      dmHybrid = NULL, dmInterface = NULL;
836b253942bSMatthew G. Knepley     DMLabel faultLabel, faultBdLabel, hybridLabel, splitLabel;
837c4762a1bSJed Brown 
8389566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
8399566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "faultBd", &faultBdLabel));
840*caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, &splitLabel, &dmInterface, &dmHybrid));
8419566063dSJacob Faibussowitsch     PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD));
8429566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
8439566063dSJacob Faibussowitsch     PetscCall(DMLabelView(splitLabel, PETSC_VIEWER_STDOUT_WORLD));
8449566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&splitLabel));
8459566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_interface_view"));
8469566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmInterface));
8479566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
848c4762a1bSJed Brown     *dm  = dmHybrid;
849c4762a1bSJed Brown   }
8509566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "fault2", &hasFault2));
851c4762a1bSJed Brown   if (hasFault2) {
852c4762a1bSJed Brown     DM      dmHybrid = NULL;
853c4762a1bSJed Brown     DMLabel faultLabel, faultBdLabel, hybridLabel;
854c4762a1bSJed Brown 
8559566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "faulted_"));
8569566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
8579566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
8589566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
8599566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
8609566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault2", &faultLabel));
8619566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault2Bd", &faultBdLabel));
862*caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, NULL, &dmHybrid));
8639566063dSJacob Faibussowitsch     PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD));
8649566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
8659566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
866c4762a1bSJed Brown     *dm  = dmHybrid;
867c4762a1bSJed Brown   }
868c4762a1bSJed Brown   if (user->testPartition && size > 1) {
869c4762a1bSJed Brown     PetscPartitioner part;
870c4762a1bSJed Brown     PetscInt *sizes  = NULL;
871c4762a1bSJed Brown     PetscInt *points = NULL;
872c4762a1bSJed Brown 
873dd400576SPatrick Sanan     if (rank == 0) {
874c4762a1bSJed Brown       if (dim == 2 && cellSimplex && size == 2) {
875c4762a1bSJed Brown         switch (user->testNum) {
876c4762a1bSJed Brown         case 0: {
877c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {1, 2};
878c4762a1bSJed Brown           PetscInt triPoints_p2[3] = {0, 1, 2};
879c4762a1bSJed Brown 
8809566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
8819566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  triSizes_p2, 2));
8829566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, triPoints_p2, 3));break;}
8835b9dfbb6SMatthew G. Knepley         case 3: {
8845b9dfbb6SMatthew G. Knepley           PetscInt triSizes_p2[2]  = {3, 3};
8855b9dfbb6SMatthew G. Knepley           PetscInt triPoints_p2[6] = {0, 1, 2,  3, 4, 5};
8865b9dfbb6SMatthew G. Knepley 
8875b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(2, &sizes, 6, &points));
8885b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes,  triSizes_p2, 2));
8895b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(points, triPoints_p2, 6));break;}
890c4762a1bSJed Brown         default:
8915b9dfbb6SMatthew G. Knepley           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 2 procs", user->testNum);
8925b9dfbb6SMatthew G. Knepley         }
8935b9dfbb6SMatthew G. Knepley       } else if (dim == 2 && cellSimplex && size == 6) {
8945b9dfbb6SMatthew G. Knepley         switch (user->testNum) {
8955b9dfbb6SMatthew G. Knepley         case 4: {
8965b9dfbb6SMatthew G. Knepley           PetscInt triSizes_p6[6]  = {1, 1, 1, 1, 1, 1};
8975b9dfbb6SMatthew G. Knepley           PetscInt triPoints_p6[6] = {0, 1, 2, 3, 4, 5};
8985b9dfbb6SMatthew G. Knepley 
8995b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(6, &sizes, 6, &points));
9005b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes,  triSizes_p6, 6));
9015b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(points, triPoints_p6, 6));break;}
9025b9dfbb6SMatthew G. Knepley         default:
9035b9dfbb6SMatthew G. Knepley           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 6 procs", user->testNum);
904c4762a1bSJed Brown         }
905c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && size == 2) {
906c4762a1bSJed Brown         switch (user->testNum) {
907c4762a1bSJed Brown         case 0: {
908c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 2};
909c4762a1bSJed Brown           PetscInt quadPoints_p2[3] = {0, 1, 2};
910c4762a1bSJed Brown 
9119566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9129566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  quadSizes_p2, 2));
9139566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, quadPoints_p2, 3));break;}
914c4762a1bSJed Brown         case 2: {
915c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 1};
916c4762a1bSJed Brown           PetscInt quadPoints_p2[2] = {0, 1};
917c4762a1bSJed Brown 
9189566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
9199566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  quadSizes_p2, 2));
9209566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));break;}
9215b9dfbb6SMatthew G. Knepley         case 3: {
9225b9dfbb6SMatthew G. Knepley           PetscInt quadSizes_p2[2]  = {1, 1};
9235b9dfbb6SMatthew G. Knepley           PetscInt quadPoints_p2[2] = {1, 0};
9245b9dfbb6SMatthew G. Knepley 
9255b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
9265b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes,  quadSizes_p2, 2));
9275b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));break;}
928c4762a1bSJed Brown         default:
9295b9dfbb6SMatthew G. Knepley           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral mesh on 2 procs", user->testNum);
930c4762a1bSJed Brown         }
931c4762a1bSJed Brown       } else if (dim == 3 && cellSimplex && size == 2) {
932c4762a1bSJed Brown         switch (user->testNum) {
933c4762a1bSJed Brown         case 0: {
934c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 2};
935c4762a1bSJed Brown           PetscInt tetPoints_p2[3] = {0, 1, 2};
936c4762a1bSJed Brown 
9379566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9389566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  tetSizes_p2, 2));
9399566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, tetPoints_p2, 3));break;}
940c4762a1bSJed Brown         default:
9415b9dfbb6SMatthew G. Knepley           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for teterehedral mesh on 2 procs", user->testNum);
942c4762a1bSJed Brown         }
943c4762a1bSJed Brown       } else if (dim == 3 && !cellSimplex && size == 2) {
944c4762a1bSJed Brown         switch (user->testNum) {
945c4762a1bSJed Brown         case 0: {
946c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {1, 2};
947c4762a1bSJed Brown           PetscInt hexPoints_p2[3] = {0, 1, 2};
948c4762a1bSJed Brown 
9499566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9509566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  hexSizes_p2, 2));
9519566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, hexPoints_p2, 3));break;}
952c4762a1bSJed Brown         default:
9535b9dfbb6SMatthew G. Knepley           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral mesh on 2 procs", user->testNum);
954c4762a1bSJed Brown         }
955c4762a1bSJed Brown       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
956c4762a1bSJed Brown     }
9579566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
9589566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
9599566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
9609566063dSJacob Faibussowitsch     PetscCall(PetscFree2(sizes, points));
961c4762a1bSJed Brown   }
962c4762a1bSJed Brown   {
963c4762a1bSJed Brown     DM pdm = NULL;
964c4762a1bSJed Brown 
965c4762a1bSJed Brown     /* Distribute mesh over processes */
9669566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
967c4762a1bSJed Brown     if (pdm) {
9689566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
9699566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
970c4762a1bSJed Brown       *dm  = pdm;
971c4762a1bSJed Brown     }
972c4762a1bSJed Brown   }
9739566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "pfault", &hasParallelFault));
974c4762a1bSJed Brown   if (hasParallelFault) {
9755b9dfbb6SMatthew G. Knepley     DM      dmHybrid = NULL, dmInterface;
976c4762a1bSJed Brown     DMLabel faultLabel, faultBdLabel, hybridLabel;
977c4762a1bSJed Brown 
9789566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "pfault", &faultLabel));
9799566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "pfaultBd", &faultBdLabel));
980*caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, &dmInterface, &dmHybrid));
9815b9dfbb6SMatthew G. Knepley     PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_fault_view"));
9825b9dfbb6SMatthew G. Knepley     {
9835b9dfbb6SMatthew G. Knepley       PetscViewer viewer;
9845b9dfbb6SMatthew G. Knepley       PetscMPIInt rank;
9855b9dfbb6SMatthew G. Knepley 
9865b9dfbb6SMatthew G. Knepley       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank));
9875b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
9885b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "Rank %d\n", rank));
9895b9dfbb6SMatthew G. Knepley       PetscCall(DMLabelView(hybridLabel, viewer));
9905b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
9915b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
9925b9dfbb6SMatthew G. Knepley     }
9939566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
9945b9dfbb6SMatthew G. Knepley     PetscCall(DMDestroy(&dmInterface));
9959566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
996c4762a1bSJed Brown     *dm  = dmHybrid;
997c4762a1bSJed Brown   }
9989566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, "Hybrid Mesh"));
9999566063dSJacob Faibussowitsch   PetscCall(CreateFaultLabel(*dm));
10009566063dSJacob Faibussowitsch   PetscCall(CreateDiscretization(*dm, user));
10019566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
10029566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
10039566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
10049566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1005c4762a1bSJed Brown   PetscFunctionReturn(0);
1006c4762a1bSJed Brown }
1007c4762a1bSJed Brown 
1008b253942bSMatthew G. Knepley static PetscErrorCode TestMesh(DM dm, AppCtx *user)
1009b253942bSMatthew G. Knepley {
1010b253942bSMatthew G. Knepley   PetscFunctionBegin;
10119566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckSymmetry(dm));
10129566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckSkeleton(dm, 0));
10139566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckFaces(dm, 0));
1014b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
1015b253942bSMatthew G. Knepley }
1016b253942bSMatthew G. Knepley 
1017b253942bSMatthew G. Knepley static PetscErrorCode TestDiscretization(DM dm, AppCtx *user)
1018b253942bSMatthew G. Knepley {
1019b253942bSMatthew G. Knepley   PetscSection   s;
1020b253942bSMatthew G. Knepley 
1021b253942bSMatthew G. Knepley   PetscFunctionBegin;
10229566063dSJacob Faibussowitsch   PetscCall(DMGetSection(dm, &s));
10239566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject) s, NULL, "-local_section_view"));
1024b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
1025b253942bSMatthew G. Knepley }
1026b253942bSMatthew G. Knepley 
1027b253942bSMatthew G. Knepley static PetscErrorCode r(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1028b253942bSMatthew G. Knepley {
1029b253942bSMatthew G. Knepley   PetscInt d;
1030b253942bSMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = x[d];
1031b253942bSMatthew G. Knepley   return 0;
1032b253942bSMatthew G. Knepley }
1033b253942bSMatthew G. Knepley 
1034b253942bSMatthew G. Knepley static PetscErrorCode rp1(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1035b253942bSMatthew G. Knepley {
1036b253942bSMatthew G. Knepley   PetscInt d;
1037b253942bSMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = x[d] + (d > 0 ? 1.0 : 0.0);
1038b253942bSMatthew G. Knepley   return 0;
1039b253942bSMatthew G. Knepley }
1040b253942bSMatthew G. Knepley 
1041b253942bSMatthew G. Knepley static PetscErrorCode phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1042b253942bSMatthew G. Knepley {
1043b253942bSMatthew G. Knepley   PetscInt d;
1044b253942bSMatthew G. Knepley   u[0] = -x[1];
1045b253942bSMatthew G. Knepley   u[1] =  x[0];
1046b253942bSMatthew G. Knepley   for (d = 2; d < dim; ++d) u[d] = x[d];
1047b253942bSMatthew G. Knepley   return 0;
1048b253942bSMatthew G. Knepley }
1049b253942bSMatthew G. Knepley 
1050b253942bSMatthew G. Knepley /* \lambda \cdot (\psi_u^- - \psi_u^+) */
1051b253942bSMatthew G. Knepley static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1052b253942bSMatthew G. Knepley                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1053b253942bSMatthew G. Knepley                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1054b253942bSMatthew G. Knepley                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1055b253942bSMatthew G. Knepley {
1056b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[1]-uOff[0];
1057b253942bSMatthew G. Knepley   PetscInt       c;
1058b253942bSMatthew G. Knepley   for (c = 0;  c < Nc;   ++c) {
1059b253942bSMatthew G. Knepley     f0[c]    =   u[Nc*2+c] + x[Nc-c-1];
1060b253942bSMatthew G. Knepley     f0[Nc+c] = -(u[Nc*2+c] + x[Nc-c-1]);
1061b253942bSMatthew G. Knepley   }
1062b253942bSMatthew G. Knepley }
1063b253942bSMatthew G. Knepley 
1064b253942bSMatthew G. Knepley /* (d - u^+ + u^-) \cdot \psi_\lambda */
1065b253942bSMatthew G. Knepley static void f0_bd_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1066b253942bSMatthew G. Knepley                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1067b253942bSMatthew G. Knepley                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1068b253942bSMatthew G. Knepley                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1069b253942bSMatthew G. Knepley {
1070b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[2]-uOff[1];
1071b253942bSMatthew G. Knepley   PetscInt       c;
1072b253942bSMatthew G. Knepley 
1073b253942bSMatthew G. Knepley   for (c = 0; c < Nc; ++c) f0[c] = (c > 0 ? 1.0 : 0.0) + u[c] - u[Nc+c];
1074b253942bSMatthew G. Knepley }
1075b253942bSMatthew G. Knepley 
1076b253942bSMatthew G. Knepley /* \psi_lambda \cdot (\psi_u^- - \psi_u^+) */
1077b253942bSMatthew G. Knepley static void g0_bd_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1078b253942bSMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1079b253942bSMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1080b253942bSMatthew G. Knepley                      PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1081b253942bSMatthew G. Knepley {
1082b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[1]-uOff[0];
1083b253942bSMatthew G. Knepley   PetscInt       c;
1084b253942bSMatthew G. Knepley 
1085b253942bSMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
1086b253942bSMatthew G. Knepley     g0[(0 +c)*Nc+c] =  1.0;
1087b253942bSMatthew G. Knepley     g0[(Nc+c)*Nc+c] = -1.0;
1088b253942bSMatthew G. Knepley   }
1089b253942bSMatthew G. Knepley }
1090b253942bSMatthew G. Knepley 
1091b253942bSMatthew G. Knepley /* (-\psi_u^+ + \psi_u^-) \cdot \psi_\lambda */
1092b253942bSMatthew G. Knepley static void g0_bd_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1093b253942bSMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1094b253942bSMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1095b253942bSMatthew G. Knepley                      PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1096b253942bSMatthew G. Knepley {
1097b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[2]-uOff[1];
1098b253942bSMatthew G. Knepley   PetscInt       c;
1099b253942bSMatthew G. Knepley 
1100b253942bSMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
1101b253942bSMatthew G. Knepley     g0[c*Nc*2+c]    =  1.0;
1102b253942bSMatthew G. Knepley     g0[c*Nc*2+Nc+c] = -1.0;
1103b253942bSMatthew G. Knepley   }
1104b253942bSMatthew G. Knepley }
1105b253942bSMatthew G. Knepley 
1106b253942bSMatthew G. Knepley static PetscErrorCode TestAssembly(DM dm, AppCtx *user)
1107b253942bSMatthew G. Knepley {
1108b253942bSMatthew G. Knepley   Mat              J;
1109b253942bSMatthew G. Knepley   Vec              locX, locF;
1110b253942bSMatthew G. Knepley   PetscDS          probh;
1111b253942bSMatthew G. Knepley   DMLabel          fault, material;
1112b253942bSMatthew G. Knepley   IS               cohesiveCells;
1113b7519becSMatthew G. Knepley   PetscWeakForm    wf;
111406ad1575SMatthew G. Knepley   PetscFormKey     keys[3];
1115b253942bSMatthew G. Knepley   PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
1116b253942bSMatthew G. Knepley   PetscInt         dim, Nf, cMax, cEnd, id;
1117b7519becSMatthew G. Knepley   PetscMPIInt      rank;
1118b253942bSMatthew G. Knepley 
1119b253942bSMatthew G. Knepley   PetscFunctionBegin;
11209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank));
11219566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
11229566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, NULL, &cMax));
11239566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
11249566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, cEnd - cMax, cMax, 1, &cohesiveCells));
11259566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &fault));
11269566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
11279566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) locX, "Local Solution"));
11289566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
11299566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) locF, "Local Residual"));
11309566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
11319566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) J, "Jacobian"));
1132b253942bSMatthew G. Knepley 
1133b253942bSMatthew G. Knepley   /* The initial guess has displacement shifted by one unit in each fault parallel direction across the fault */
11349566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "material", &material));
1135b253942bSMatthew G. Knepley   id   = 1;
1136b253942bSMatthew G. Knepley   initialGuess[0] = r;
1137b253942bSMatthew G. Knepley   initialGuess[1] = NULL;
11389566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
1139b253942bSMatthew G. Knepley   id   = 2;
1140b253942bSMatthew G. Knepley   initialGuess[0] = rp1;
1141b253942bSMatthew G. Knepley   initialGuess[1] = NULL;
11429566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
1143b253942bSMatthew G. Knepley   id   = 1;
1144b253942bSMatthew G. Knepley   initialGuess[0] = NULL;
1145b253942bSMatthew G. Knepley   initialGuess[1] = phi;
11469566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, fault, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
11479566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(locX, NULL, "-local_solution_view"));
1148b253942bSMatthew G. Knepley 
11499566063dSJacob Faibussowitsch   PetscCall(DMGetCellDS(dm, cMax, &probh));
11509566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(probh, &wf));
11519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(probh, &Nf));
11529566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 1, 0, 0, 0, f0_bd_u, 0, NULL));
11539566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 2, 0, 0, 0, f0_bd_u, 0, NULL));
11549566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 1, 0, 1, 0, 0, g0_bd_ul, 0, NULL, 0, NULL, 0, NULL));
11559566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 2, 0, 1, 0, 0, g0_bd_ul, 0, NULL, 0, NULL, 0, NULL));
1156b7519becSMatthew G. Knepley   if (Nf > 1) {
11579566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdResidual(wf, fault, 1, 1, 0, 0, f0_bd_l, 0, NULL));
11589566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdJacobian(wf, fault, 1, 1, 0, 0, 0, g0_bd_lu, 0, NULL, 0, NULL, 0, NULL));
1159b7519becSMatthew G. Knepley   }
11609566063dSJacob Faibussowitsch   if (!rank) PetscCall(PetscDSView(probh, NULL));
1161b253942bSMatthew G. Knepley 
1162dd96b1f9SToby Isaac   keys[0].label = NULL;
1163dd96b1f9SToby Isaac   keys[0].value = 0;
11646528b96dSMatthew G. Knepley   keys[0].field = 0;
1165dd96b1f9SToby Isaac   keys[0].part  = 0;
1166b7519becSMatthew G. Knepley   keys[1].label = material;
1167b7519becSMatthew G. Knepley   keys[1].value = 2;
11686528b96dSMatthew G. Knepley   keys[1].field = 0;
1169dd96b1f9SToby Isaac   keys[1].part  = 0;
1170b7519becSMatthew G. Knepley   keys[2].label = fault;
1171b7519becSMatthew G. Knepley   keys[2].value = 1;
1172b7519becSMatthew G. Knepley   keys[2].field = 1;
1173dd96b1f9SToby Isaac   keys[2].part  = 0;
11749566063dSJacob Faibussowitsch   PetscCall(VecSet(locF, 0.));
11759566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeResidual_Hybrid_Internal(dm, keys, cohesiveCells, 0.0, locX, NULL, 0.0, locF, user));
11769566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(locF, NULL, "-local_residual_view"));
11779566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(J));
11789566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeJacobian_Hybrid_Internal(dm, keys, cohesiveCells, 0.0, 0.0, locX, NULL, J, J, user));
11799566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-local_jacobian_view"));
1180b253942bSMatthew G. Knepley 
11819566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
11829566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
11839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
11849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cohesiveCells));
1185b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
1186b253942bSMatthew G. Knepley }
1187b253942bSMatthew G. Knepley 
1188c4762a1bSJed Brown int main(int argc, char **argv)
1189c4762a1bSJed Brown {
1190c4762a1bSJed Brown   DM             dm;
1191c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
1192c4762a1bSJed Brown 
11939566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
11949566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
11959566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
11969566063dSJacob Faibussowitsch   PetscCall(TestMesh(dm, &user));
11979566063dSJacob Faibussowitsch   PetscCall(TestDiscretization(dm, &user));
11989566063dSJacob Faibussowitsch   PetscCall(TestAssembly(dm, &user));
11999566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
12009566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1201b122ec5aSJacob Faibussowitsch   return 0;
1202c4762a1bSJed Brown }
1203c4762a1bSJed Brown 
1204c4762a1bSJed Brown /*TEST
1205c4762a1bSJed Brown   testset:
1206b253942bSMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1207b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1 -local_section_view \
1208b253942bSMatthew G. Knepley           -local_solution_view -local_residual_view -local_jacobian_view
1209c4762a1bSJed Brown     test:
1210c4762a1bSJed Brown       suffix: tri_0
1211c4762a1bSJed Brown       args: -dim 2
12121c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1213c4762a1bSJed Brown     test:
1214c4762a1bSJed Brown       suffix: tri_t1_0
1215c4762a1bSJed Brown       args: -dim 2 -test_num 1
12161c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1217c4762a1bSJed Brown     test:
12185b9dfbb6SMatthew G. Knepley       suffix: tri_t2_0
12195b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 2
12205b9dfbb6SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
12215b9dfbb6SMatthew G. Knepley     test:
1222c4762a1bSJed Brown       suffix: tet_0
1223c4762a1bSJed Brown       args: -dim 3
12241c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1225c4762a1bSJed Brown     test:
1226b253942bSMatthew G. Knepley       suffix: tet_t1_0
1227b253942bSMatthew G. Knepley       args: -dim 3 -test_num 1
12281c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1229b253942bSMatthew G. Knepley 
1230b253942bSMatthew G. Knepley   testset:
1231b253942bSMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1232b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1
1233b253942bSMatthew G. Knepley     test:
1234c4762a1bSJed Brown       suffix: tet_1
1235c4762a1bSJed Brown       nsize: 2
1236c4762a1bSJed Brown       args: -dim 3
12371c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1238c4762a1bSJed Brown     test:
1239b253942bSMatthew G. Knepley       suffix: tri_1
1240b253942bSMatthew G. Knepley       nsize: 2
1241b253942bSMatthew G. Knepley       args: -dim 2
12421c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
12435b9dfbb6SMatthew G. Knepley     test:
12445b9dfbb6SMatthew G. Knepley       suffix: tri_t3_0
12455b9dfbb6SMatthew G. Knepley       nsize: 2
12465b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 3
12475b9dfbb6SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
12485b9dfbb6SMatthew G. Knepley     test:
12495b9dfbb6SMatthew G. Knepley       suffix: tri_t4_0
12505b9dfbb6SMatthew G. Knepley       nsize: 6
12515b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 4
12525b9dfbb6SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1253c4762a1bSJed Brown 
1254c4762a1bSJed Brown   testset:
1255ecfb78b5SMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1256b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1
1257c4762a1bSJed Brown     # 2D Quads
1258c4762a1bSJed Brown     test:
1259c4762a1bSJed Brown       suffix: quad_0
1260c4762a1bSJed Brown       args: -dim 2 -cell_simplex 0
12611c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1262c4762a1bSJed Brown     test:
1263c4762a1bSJed Brown       suffix: quad_1
1264c4762a1bSJed Brown       nsize: 2
1265c4762a1bSJed Brown       args: -dim 2 -cell_simplex 0
12661c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1267c4762a1bSJed Brown     test:
1268c4762a1bSJed Brown       suffix: quad_t1_0
126954fcfd0cSMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 1 -faulted_dm_plex_check_all
12701c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
12715b9dfbb6SMatthew G. Knepley     test:
12725b9dfbb6SMatthew G. Knepley       suffix: quad_t2_0
12735b9dfbb6SMatthew G. Knepley       nsize: 2
12745b9dfbb6SMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 2
12755b9dfbb6SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
12765b9dfbb6SMatthew G. Knepley     test:
12775b9dfbb6SMatthew G. Knepley       # TODO: The PetscSF is wrong here (connects to wrong side of split)
12785b9dfbb6SMatthew G. Knepley       suffix: quad_t3_0
12795b9dfbb6SMatthew G. Knepley       nsize: 2
12805b9dfbb6SMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 3
12815b9dfbb6SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1282c4762a1bSJed Brown     # 3D Hex
1283c4762a1bSJed Brown     test:
1284c4762a1bSJed Brown       suffix: hex_0
1285c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0
12861c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1287c4762a1bSJed Brown     test:
1288c4762a1bSJed Brown       suffix: hex_1
1289c4762a1bSJed Brown       nsize: 2
1290c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0
12911c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1292c4762a1bSJed Brown     test:
1293c4762a1bSJed Brown       suffix: hex_t1_0
1294c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0 -test_num 1
12951c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1296c4762a1bSJed Brown     test:
1297c4762a1bSJed Brown       suffix: hex_t2_0
1298c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0 -test_num 2
12991c6715b8SMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul//g" -e "s/g0_bd_lu//g"
1300c4762a1bSJed Brown 
1301c4762a1bSJed Brown TEST*/
1302