xref: /petsc/src/dm/impls/plex/tests/ex5.c (revision 49664dce07a55e12bbee01dbb7f26988e081dbef)
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 
1553bdf08b3SMatthew G. Knepley Test 5:
1563bdf08b3SMatthew G. Knepley 
1573bdf08b3SMatthew G. Knepley   Fault only on points 2 and 5:
1583bdf08b3SMatthew G. Knepley 
1593bdf08b3SMatthew G. Knepley         6
1603bdf08b3SMatthew G. Knepley       / | \
1613bdf08b3SMatthew G. Knepley     13  |  17
1623bdf08b3SMatthew G. Knepley     /  15   \
1633bdf08b3SMatthew G. Knepley    7  0 | 1  9
1643bdf08b3SMatthew G. Knepley    |\   |   /|
1653bdf08b3SMatthew G. Knepley    | 14 | 16 |
1663bdf08b3SMatthew G. Knepley    |  \ | /  |
1673bdf08b3SMatthew G. Knepley  18| 2  8  3 |21
1683bdf08b3SMatthew G. Knepley    |  / | \  |
1693bdf08b3SMatthew G. Knepley    | 19 | 20 |
1703bdf08b3SMatthew G. Knepley    |/   |   \|
1713bdf08b3SMatthew G. Knepley   10  4 | 5  12
1723bdf08b3SMatthew G. Knepley     \  23   /
1733bdf08b3SMatthew G. Knepley     22  |  24
1743bdf08b3SMatthew G. Knepley       \ | /
1753bdf08b3SMatthew G. Knepley        11
1763bdf08b3SMatthew G. Knepley 
177c4762a1bSJed Brown Tetrahedron
178c4762a1bSJed Brown -----------
179c4762a1bSJed Brown Test 0:
180c4762a1bSJed Brown Two tets sharing a face
181c4762a1bSJed Brown 
182c4762a1bSJed Brown  cell   5 _______    cell
183c4762a1bSJed Brown  0    / | \      \       1
184c4762a1bSJed Brown     16  |  18     22
185c4762a1bSJed Brown     /8 19 10\      \
186c4762a1bSJed Brown    2-15-|----4--21--6
187c4762a1bSJed Brown     \  9| 7 /      /
188c4762a1bSJed Brown     14  |  17     20
189c4762a1bSJed Brown       \ | /      /
190c4762a1bSJed Brown         3-------
191c4762a1bSJed Brown 
192c4762a1bSJed Brown should become two tetrahedrons separated by a zero-volume cell with 3 faces/3 edges/6 vertices
193c4762a1bSJed Brown 
194c4762a1bSJed Brown  cell   6 ___36___10______    cell
195c4762a1bSJed Brown  0    / | \        |\      \     1
196c4762a1bSJed Brown     24  |  26      | 32     30
197c4762a1bSJed Brown     /12 27 14\    33  \      \
198c4762a1bSJed Brown    3-23-|----5--35-|---9--29--7
199c4762a1bSJed Brown     \ 13| 11/      |18 /      /
200c4762a1bSJed Brown     22  |  25      | 31     28
201c4762a1bSJed Brown       \ | /        |/      /
202c4762a1bSJed Brown         4----34----8------
203c4762a1bSJed Brown          cell 2
204c4762a1bSJed Brown 
205c4762a1bSJed Brown In parallel,
206c4762a1bSJed Brown 
207c4762a1bSJed Brown  cell   5 ___28____8      4______    cell
208c4762a1bSJed Brown  0    / | \        |\     |\      \     0
209c4762a1bSJed Brown     19  |   21     | 24   | 13  6  11
210c4762a1bSJed Brown     /10 22 12\    25  \   |8 \      \
211c4762a1bSJed Brown    2-18-|----4--27-|---7  14  3--10--1
212c4762a1bSJed Brown     \ 11| 9 /      |13 /  |  /      /
213c4762a1bSJed Brown     17  |  20      | 23   | 12  5  9
214c4762a1bSJed Brown       \ | /        |/     |/      /
215c4762a1bSJed Brown         3----26----6      2------
216c4762a1bSJed Brown          cell 1
217c4762a1bSJed Brown 
218c4762a1bSJed Brown Test 1:
219c4762a1bSJed Brown Four tets sharing two faces
220c4762a1bSJed Brown 
221c4762a1bSJed Brown Cells:    0-3,4-5
222c4762a1bSJed Brown Vertices: 6-15
223c4762a1bSJed Brown Faces:    16-29,30-34
224c4762a1bSJed Brown Edges:    35-52,53-56
225c4762a1bSJed Brown 
226c4762a1bSJed Brown Quadrilateral
227c4762a1bSJed Brown -------------
228c4762a1bSJed Brown Test 0:
229c4762a1bSJed Brown Two quads sharing a face
230c4762a1bSJed Brown 
231c4762a1bSJed Brown    5--10---4--14---7
232c4762a1bSJed Brown    |       |       |
233c4762a1bSJed Brown   11   0   9   1  13
234c4762a1bSJed Brown    |       |       |
235c4762a1bSJed Brown    2---8---3--12---6
236c4762a1bSJed Brown 
237c4762a1bSJed Brown should become two quads separated by a zero-volume cell with 4 vertices
238c4762a1bSJed Brown 
239c4762a1bSJed Brown    6--13---5-20-10--17---8    5--10---4-14--7  4---7---2
240c4762a1bSJed Brown    |       |     |       |    |       |     |  |       |
241c4762a1bSJed Brown   14   0  12  2 18   1  16   11   0   9  1 12  8   0   6
242c4762a1bSJed Brown    |       |     |       |    |       |     |  |       |
243c4762a1bSJed Brown    3--11---4-19--9--15---7    2---8---3-13--6  3---5---1
244c4762a1bSJed Brown 
245c4762a1bSJed Brown Test 1:
246c4762a1bSJed Brown 
247c4762a1bSJed Brown Original mesh with 9 cells,
248c4762a1bSJed Brown 
2495b9dfbb6SMatthew G. Knepley   9-----10-----11-----12
2505b9dfbb6SMatthew G. Knepley   |      |     ||      |
2515b9dfbb6SMatthew G. Knepley   |      |     ||      |
2525b9dfbb6SMatthew G. Knepley   |   0  |  1  ||  2   |
2535b9dfbb6SMatthew G. Knepley   |      |     ||      |
2545b9dfbb6SMatthew G. Knepley  13-----14-----15-----16
2555b9dfbb6SMatthew G. Knepley   |      |     ||      |
2565b9dfbb6SMatthew G. Knepley   |      |     ||      |
2575b9dfbb6SMatthew G. Knepley   |  3   |  4  ||  5   |
2585b9dfbb6SMatthew G. Knepley   |      |     ||      |
2595b9dfbb6SMatthew G. Knepley  17-----18-----19=====20
260c4762a1bSJed Brown   |      |      |      |
261c4762a1bSJed Brown   |      |      |      |
2625b9dfbb6SMatthew G. Knepley   |  6   |  7   |  8   |
263c4762a1bSJed Brown   |      |      |      |
2645b9dfbb6SMatthew G. Knepley  21-----22-----23-----24
265c4762a1bSJed Brown 
266c4762a1bSJed Brown After first fault,
267c4762a1bSJed Brown 
268c4762a1bSJed Brown  12 ----13 ----14-28 ----15
269c4762a1bSJed Brown   |      |      |  |      |
270c4762a1bSJed Brown   |  0   |  1   | 9|  2   |
271c4762a1bSJed Brown   |      |      |  |      |
272c4762a1bSJed Brown   |      |      |  |      |
273c4762a1bSJed Brown  16 ----17 ----18-29 ----19
274c4762a1bSJed Brown   |      |      |  |      |
275c4762a1bSJed Brown   |  3   |  4   |10|  5   |
276c4762a1bSJed Brown   |      |      |  |      |
277c4762a1bSJed Brown   |      |      |  |      |
278c4762a1bSJed Brown  20 ----21-----22-30 ----23
279c4762a1bSJed Brown   |      |      |  \--11- |
280c4762a1bSJed Brown   |  6   |  7   |     8   |
281c4762a1bSJed Brown   |      |      |         |
282c4762a1bSJed Brown   |      |      |         |
283c4762a1bSJed Brown  24 ----25 ----26--------27
284c4762a1bSJed Brown 
285c4762a1bSJed Brown After second fault,
286c4762a1bSJed Brown 
287c4762a1bSJed Brown  14 ----15 ----16-30 ----17
288c4762a1bSJed Brown   |      |      |  |      |
289c4762a1bSJed Brown   |  0   |  1   | 9|  2   |
290c4762a1bSJed Brown   |      |      |  |      |
291c4762a1bSJed Brown   |      |      |  |      |
292c4762a1bSJed Brown  18 ----19 ----20-31 ----21
293c4762a1bSJed Brown   |      |      |  |      |
294c4762a1bSJed Brown   |  3   |  4   |10|  5   |
295c4762a1bSJed Brown   |      |      |  |      |
296c4762a1bSJed Brown   |      |      |  |      |
297c4762a1bSJed Brown  33 ----34-----24-32 ----25
298c4762a1bSJed Brown   |  12  | 13 / |  \-11-- |
299c4762a1bSJed Brown  22 ----23---/  |         |
3005b9dfbb6SMatthew G. Knepley   |      |      |         |
3015b9dfbb6SMatthew G. Knepley   |  6   |   7  |     8   |
302c4762a1bSJed Brown   |      |      |         |
303c4762a1bSJed Brown   |      |      |         |
304c4762a1bSJed Brown  26 ----27 ----28--------29
305c4762a1bSJed Brown 
3065b9dfbb6SMatthew G. Knepley  Test 2:
3075b9dfbb6SMatthew G. Knepley  Two quads sharing a face in parallel
3085b9dfbb6SMatthew G. Knepley 
3095b9dfbb6SMatthew G. Knepley     4---7---3  2---8---4
3105b9dfbb6SMatthew G. Knepley     |       |  |       |
3115b9dfbb6SMatthew G. Knepley     8   0   6  5   0   7
3125b9dfbb6SMatthew G. Knepley     |       |  |       |
3135b9dfbb6SMatthew G. Knepley     1---5---2  1---6---3
3145b9dfbb6SMatthew G. Knepley 
3155b9dfbb6SMatthew G. Knepley  should become two quads separated by a zero-volume cell with 4 vertices
3165b9dfbb6SMatthew G. Knepley 
3175b9dfbb6SMatthew G. Knepley      4---7---3  3-14--7--11---5
3185b9dfbb6SMatthew G. Knepley      |       |  |     |       |
3195b9dfbb6SMatthew G. Knepley      8   0   6  8  1  12  0   10
3205b9dfbb6SMatthew G. Knepley      |       |  |     |       |
3215b9dfbb6SMatthew G. Knepley      1---5---2  2-13--6---9---4
3225b9dfbb6SMatthew G. Knepley 
3235b9dfbb6SMatthew G. Knepley  Test 3:
3245b9dfbb6SMatthew G. Knepley  Like Test 2, but with different partition
3255b9dfbb6SMatthew G. Knepley 
3265b9dfbb6SMatthew G. Knepley      5--10---4-14--7   2---8---4
3275b9dfbb6SMatthew G. Knepley      |       |     |   |       |
3285b9dfbb6SMatthew G. Knepley     11   0   9  1  12  5   0   7
3295b9dfbb6SMatthew G. Knepley      |       |     |   |       |
3305b9dfbb6SMatthew G. Knepley      2---8---3-13--6   1---6---3
3315b9dfbb6SMatthew G. Knepley 
332c4762a1bSJed Brown Hexahedron
333c4762a1bSJed Brown ----------
334c4762a1bSJed Brown Test 0:
335c4762a1bSJed Brown Two hexes sharing a face
336c4762a1bSJed Brown 
337c4762a1bSJed Brown cell   9-----31------8-----42------13 cell
338c4762a1bSJed Brown 0     /|            /|            /|     1
339c4762a1bSJed Brown     32 |   15      30|   21      41|
340c4762a1bSJed Brown     /  |          /  |          /  |
341c4762a1bSJed Brown    6-----29------7-----40------12  |
342c4762a1bSJed Brown    |   |     18  |   |     24  |   |
343c4762a1bSJed Brown    |  36         |  35         |   44
344c4762a1bSJed Brown    |19 |         |17 |         |23 |
345c4762a1bSJed Brown   33   |  16    34   |   22   43   |
346c4762a1bSJed Brown    |   5-----27--|---4-----39--|---11
347c4762a1bSJed Brown    |  /          |  /          |  /
348c4762a1bSJed Brown    | 28   14     | 26    20    | 38
349c4762a1bSJed Brown    |/            |/            |/
350c4762a1bSJed Brown    2-----25------3-----37------10
351c4762a1bSJed Brown 
352c4762a1bSJed Brown should become two hexes separated by a zero-volume cell with 8 vertices
353c4762a1bSJed Brown 
354c4762a1bSJed Brown                          cell 2
355c4762a1bSJed Brown cell  10-----41------9-----62------18----52------14 cell
356c4762a1bSJed Brown 0     /|            /|            /|            /|     1
357c4762a1bSJed Brown     42 |   20      40|  32       56|   26      51|
358c4762a1bSJed Brown     /  |          /  |          /  |          /  |
359c4762a1bSJed Brown    7-----39------8-----61------17--|-50------13  |
360c4762a1bSJed Brown    |   |     23  |   |         |   |     29  |   |
361c4762a1bSJed Brown    |  46         |  45         |   58        |   54
362c4762a1bSJed Brown    |24 |         |22 |         |30 |         |28 |
363c4762a1bSJed Brown   43   |  21    44   |        57   |   27   53   |
364c4762a1bSJed Brown    |   6-----37--|---5-----60--|---16----49--|---12
365c4762a1bSJed Brown    |  /          |  /          |  /          |  /
366c4762a1bSJed Brown    | 38   19     | 36   31     | 55    25    | 48
367c4762a1bSJed Brown    |/            |/            |/            |/
368c4762a1bSJed Brown    3-----35------4-----59------15----47------11
369c4762a1bSJed Brown 
370c4762a1bSJed Brown In parallel,
371c4762a1bSJed Brown 
372c4762a1bSJed Brown                          cell 2
373c4762a1bSJed Brown cell   9-----31------8-----44------13     8----20------4  cell
374c4762a1bSJed Brown 0     /|            /|            /|     /|           /|     1
375c4762a1bSJed Brown     32 |   15      30|  22       38|   24 |  10      19|
376c4762a1bSJed Brown     /  |          /  |          /  |   /  |         /  |
377c4762a1bSJed Brown    6-----29------7-----43------12  |  7----18------3   |
378c4762a1bSJed Brown    |   |     18  |   |         |   |  |   |    13  |   |
379c4762a1bSJed Brown    |  36         |  35         |   40 |  26        |   22
380c4762a1bSJed Brown    |19 |         |17 |         |20 |  |14 |        |12 |
381c4762a1bSJed Brown   33   |  16    34   |        39   |  25  |  11   21   |
382c4762a1bSJed Brown    |   5-----27--|---4-----42--|---11 |   6----17--|---2
383c4762a1bSJed Brown    |  /          |  /          |  /   |  /         |  /
384c4762a1bSJed Brown    | 28   14     | 26   21     | 37   |23     9    | 16
385c4762a1bSJed Brown    |/            |/            |/     |/           |/
386c4762a1bSJed Brown    2-----25------3-----41------10     5----15------1
387c4762a1bSJed Brown 
388c4762a1bSJed Brown Test 1:
389c4762a1bSJed Brown 
390c4762a1bSJed Brown */
391c4762a1bSJed Brown 
392c4762a1bSJed Brown typedef struct {
393c4762a1bSJed Brown   PetscInt  debug;          /* The debugging level */
394c4762a1bSJed Brown   PetscInt  dim;            /* The topological mesh dimension */
395c4762a1bSJed Brown   PetscBool cellSimplex;    /* Use simplices or hexes */
396c4762a1bSJed Brown   PetscBool testPartition;  /* Use a fixed partitioning for testing */
397c4762a1bSJed Brown   PetscInt  testNum;        /* The particular mesh to test */
398b7519becSMatthew G. Knepley   PetscInt  cohesiveFields; /* The number of cohesive fields */
399c4762a1bSJed Brown } AppCtx;
400c4762a1bSJed Brown 
401b253942bSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
402c4762a1bSJed Brown {
403c4762a1bSJed Brown   PetscFunctionBegin;
404c4762a1bSJed Brown   options->debug          = 0;
405c4762a1bSJed Brown   options->dim            = 2;
406c4762a1bSJed Brown   options->cellSimplex    = PETSC_TRUE;
407c4762a1bSJed Brown   options->testPartition  = PETSC_TRUE;
408c4762a1bSJed Brown   options->testNum        = 0;
409b7519becSMatthew G. Knepley   options->cohesiveFields = 1;
410c4762a1bSJed Brown 
411d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
4129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex5.c", options->debug, &options->debug, NULL,0));
4139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex5.c", options->dim, &options->dim, NULL,1,3));
4149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex5.c", options->cellSimplex, &options->cellSimplex, NULL));
4159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex5.c", options->testPartition, &options->testPartition, NULL));
4169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex5.c", options->testNum, &options->testNum, NULL,0));
4179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-cohesive_fields", "The number of cohesive fields", "ex5.c", options->cohesiveFields, &options->cohesiveFields, NULL, 0));
418d0609cedSBarry Smith   PetscOptionsEnd();
419c4762a1bSJed Brown   PetscFunctionReturn(0);
420c4762a1bSJed Brown }
421c4762a1bSJed Brown 
422b253942bSMatthew G. Knepley static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
423c4762a1bSJed Brown {
424c4762a1bSJed Brown   DM             idm;
425c4762a1bSJed Brown   PetscInt       p;
426c4762a1bSJed Brown   PetscMPIInt    rank;
427c4762a1bSJed Brown 
428c4762a1bSJed Brown   PetscFunctionBegin;
4299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
430dd400576SPatrick Sanan   if (rank == 0) {
431c4762a1bSJed Brown     switch (testNum) {
432c4762a1bSJed Brown     case 0:
433c4762a1bSJed Brown     {
434c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
435c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
436c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
437c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
438c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
439c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
440c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
441c4762a1bSJed Brown 
4429566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4439566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
4449566063dSJacob Faibussowitsch       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4459566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
4469566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
447c4762a1bSJed Brown     }
448c4762a1bSJed Brown     break;
449c4762a1bSJed Brown     case 1:
450c4762a1bSJed Brown     {
451c4762a1bSJed Brown       PetscInt    numPoints[2]         = {6, 4};
452c4762a1bSJed Brown       PetscInt    coneSize[10]         = {3, 3, 3, 3, 0, 0, 0, 0, 0, 0};
453c4762a1bSJed Brown       PetscInt    cones[12]            = {4, 6, 5,  5, 6, 7,  8, 5, 9,  8, 4, 5};
454c4762a1bSJed Brown       PetscInt    coneOrientations[12] = {0, 0, 0,  0, 0, 0,  0, 0, 0,  0, 0, 0};
455c4762a1bSJed 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};
456c4762a1bSJed Brown       PetscInt    markerPoints[6]      = {4, 1, 6, 1, 8, 1};
457c4762a1bSJed Brown       PetscInt    faultPoints[3]       = {5, 6, 8};
458c4762a1bSJed Brown 
4599566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4609566063dSJacob Faibussowitsch       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
4619566063dSJacob Faibussowitsch       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4629566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));PetscCall(DMSetLabelValue(*dm, "material", 3, 1));
4639566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));PetscCall(DMSetLabelValue(*dm, "material", 2, 2));
464c4762a1bSJed Brown     }
465c4762a1bSJed Brown     break;
4665b9dfbb6SMatthew G. Knepley     case 2:
4675b9dfbb6SMatthew G. Knepley     case 3:
4685b9dfbb6SMatthew G. Knepley     case 4:
4695b9dfbb6SMatthew G. Knepley     {
4705b9dfbb6SMatthew G. Knepley       PetscInt    numPoints[2]         = {8, 6};
4715b9dfbb6SMatthew G. Knepley       PetscInt    coneSize[14]         = {3, 3, 3, 3, 3, 3, 0, 0, 0, 0, 0, 0, 0, 0};
4725b9dfbb6SMatthew G. Knepley       PetscInt    cones[18]            = {6, 7,  9,   9, 12, 11,  7, 12,  9,
4735b9dfbb6SMatthew G. Knepley                                           7, 8, 10,  10, 13, 12,  7, 10, 12};
4745b9dfbb6SMatthew G. Knepley       PetscInt    coneOrientations[18] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
4755b9dfbb6SMatthew G. Knepley       PetscScalar vertexCoords[16]     = {-1., -1.,  0., -1.,  1., -1.,  -1., 0.,  1., 0.,
4765b9dfbb6SMatthew G. Knepley                                           -1.,  1.,  0.,  1.,  1.,  1.,};
4775b9dfbb6SMatthew G. Knepley       PetscInt    markerPoints[16]     = {6, 1, 7, 1, 8, 1, 9, 1, 10, 1, 11, 1, 12, 1, 13, 1};
4785b9dfbb6SMatthew G. Knepley       PetscInt    faultPoints[2]       = {7, 12};
4795b9dfbb6SMatthew G. Knepley 
4805b9dfbb6SMatthew G. Knepley       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4815b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
4825b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
4835b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
4845b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
4855b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 4, 2));
4865b9dfbb6SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
4875b9dfbb6SMatthew G. Knepley       for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
4885b9dfbb6SMatthew G. Knepley       if (testNum == 2) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4895b9dfbb6SMatthew G. Knepley       if (testNum == 3 || testNum == 4) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1));
4905b9dfbb6SMatthew G. Knepley     }
4915b9dfbb6SMatthew G. Knepley     break;
4923bdf08b3SMatthew G. Knepley     case 5:
4933bdf08b3SMatthew G. Knepley     {
4943bdf08b3SMatthew G. Knepley       PetscInt    numPoints[2]         = {7, 6};
4953bdf08b3SMatthew G. Knepley       PetscInt    coneSize[13]         = {3, 3, 3, 3, 3, 3, 0, 0, 0, 0, 0, 0, 0};
4963bdf08b3SMatthew G. Knepley       PetscInt    cones[18]            = {6, 7, 8,  8, 9, 6,  7, 10, 8,  9, 8, 12,  8, 10, 11,  11, 12, 8};
4973bdf08b3SMatthew G. Knepley       PetscInt    coneOrientations[18] = {0, 0, 0,  0, 0, 0,  0,  0, 0,  0, 0,  0,  0,  0,  0,   0,  0, 0};
4983bdf08b3SMatthew G. Knepley       PetscScalar vertexCoords[14]     = {0.0, 2.0,  -1.0, 1.0,  0.0, 0.0,  1.0, 1.0,  -1.0, -1.0,  0.0, -2.0,  1.0, -1.0};
4993bdf08b3SMatthew G. Knepley       PetscInt    faultPoints[2]       = {8, 11};
5003bdf08b3SMatthew G. Knepley       PetscInt    faultBdPoints[1]     = {8};
5013bdf08b3SMatthew G. Knepley 
5023bdf08b3SMatthew G. Knepley       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
5033bdf08b3SMatthew G. Knepley       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
5043bdf08b3SMatthew G. Knepley       for (p = 0; p < 1; ++p) PetscCall(DMSetLabelValue(*dm, "faultBd", faultBdPoints[p], 1));
5053bdf08b3SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 0, 0));PetscCall(DMSetLabelValue(*dm, "material", 2, 0));PetscCall(DMSetLabelValue(*dm, "material", 4, 0));
5063bdf08b3SMatthew G. Knepley       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));PetscCall(DMSetLabelValue(*dm, "material", 3, 2));PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
5073bdf08b3SMatthew G. Knepley     }
5083bdf08b3SMatthew G. Knepley     break;
509c4762a1bSJed Brown     default:
51063a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
511c4762a1bSJed Brown     }
512c4762a1bSJed Brown   } else {
513c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
514c4762a1bSJed Brown 
5159566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
5165b9dfbb6SMatthew G. Knepley     if (testNum == 3 || testNum == 4) PetscCall(DMCreateLabel(*dm, "pfault"));
5175b9dfbb6SMatthew G. Knepley     else                              PetscCall(DMCreateLabel(*dm, "fault"));
518c4762a1bSJed Brown   }
5199566063dSJacob Faibussowitsch   PetscCall(DMPlexInterpolate(*dm, &idm));
5209566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
5219566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
522c4762a1bSJed Brown   *dm  = idm;
523c4762a1bSJed Brown   PetscFunctionReturn(0);
524c4762a1bSJed Brown }
525c4762a1bSJed Brown 
526b253942bSMatthew G. Knepley static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, AppCtx *user, DM dm)
527c4762a1bSJed Brown {
528c4762a1bSJed Brown   PetscInt       depth = 3, testNum  = user->testNum, p;
529c4762a1bSJed Brown   PetscMPIInt    rank;
530c4762a1bSJed Brown 
531c4762a1bSJed Brown   PetscFunctionBegin;
5329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
533dd400576SPatrick Sanan   if (rank == 0) {
534c4762a1bSJed Brown     switch (testNum) {
535c4762a1bSJed Brown     case 0:
536c4762a1bSJed Brown     {
537c4762a1bSJed Brown       PetscInt    numPoints[4]         = {5, 7, 9, 2};
538c4762a1bSJed 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};
539c4762a1bSJed 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};
540b5a892a1SMatthew 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};
541c4762a1bSJed 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};
542c4762a1bSJed Brown       PetscInt    markerPoints[20]     = {2, 1, 3, 1, 4, 1, 5, 1, 14, 1, 15, 1, 16, 1, 17, 1, 18, 1, 19, 1};
543c4762a1bSJed Brown       PetscInt    faultPoints[3]      = {3, 4, 5};
544c4762a1bSJed Brown 
5459566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
546c4762a1bSJed Brown       for (p = 0; p < 10; ++p) {
5479566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
548c4762a1bSJed Brown       }
549c4762a1bSJed Brown       for (p = 0; p < 3; ++p) {
5509566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1));
551c4762a1bSJed Brown       }
5529566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 0, 1));
5539566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 1, 2));
554c4762a1bSJed Brown     }
555c4762a1bSJed Brown     break;
556c4762a1bSJed Brown     case 1:
557c4762a1bSJed Brown     {
558c4762a1bSJed Brown       PetscInt    numPoints[4]         = {6, 13, 12, 4};
559c4762a1bSJed 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};
560c4762a1bSJed 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};
561b5a892a1SMatthew 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};
562c4762a1bSJed 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};
563c4762a1bSJed Brown       PetscInt    markerPoints[14]     = {5, 1, 6, 1, 7, 1, 10, 1, 22, 1, 23, 1, 24, 1};
564c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {5, 6, 7, 8};
565c4762a1bSJed Brown 
5669566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
567c4762a1bSJed Brown       for (p = 0; p < 7; ++p) {
5689566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
569c4762a1bSJed Brown       }
570c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
5719566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1));
572c4762a1bSJed Brown       }
5739566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 0, 1));
5749566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 1, 2));
575c4762a1bSJed Brown     }
576c4762a1bSJed Brown     break;
577c4762a1bSJed Brown     default:
57863a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
579c4762a1bSJed Brown     }
580c4762a1bSJed Brown   } else {
581c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
582c4762a1bSJed Brown 
5839566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
5849566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "fault"));
585c4762a1bSJed Brown   }
586c4762a1bSJed Brown   PetscFunctionReturn(0);
587c4762a1bSJed Brown }
588c4762a1bSJed Brown 
589b253942bSMatthew G. Knepley static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
590c4762a1bSJed Brown {
591c4762a1bSJed Brown   DM             idm;
592c4762a1bSJed Brown   PetscInt       p;
593c4762a1bSJed Brown   PetscMPIInt    rank;
594c4762a1bSJed Brown 
595c4762a1bSJed Brown   PetscFunctionBegin;
5969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
597dd400576SPatrick Sanan   if (rank == 0) {
598c4762a1bSJed Brown     switch (testNum) {
599c4762a1bSJed Brown     case 0:
600c4762a1bSJed Brown     case 2:
6015b9dfbb6SMatthew G. Knepley     case 3:
602c4762a1bSJed Brown     {
603c4762a1bSJed Brown       PetscInt    numPoints[2]        = {6, 2};
604c4762a1bSJed Brown       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
605c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
606c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
607c4762a1bSJed 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};
608c4762a1bSJed Brown       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
609c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
610c4762a1bSJed Brown 
6119566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6129566063dSJacob Faibussowitsch       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
6139566063dSJacob Faibussowitsch       if (testNum == 0) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
6145b9dfbb6SMatthew G. Knepley       if (testNum == 2 || testNum == 3) for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1));
6159566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6169566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
617c4762a1bSJed Brown     }
618c4762a1bSJed Brown     break;
619c4762a1bSJed Brown     case 1:
620c4762a1bSJed Brown     {
621c4762a1bSJed Brown       PetscInt    numPoints[2]         = {16, 9};
622c4762a1bSJed 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};
623c4762a1bSJed Brown       PetscInt    cones[36]            = {9,  13, 14, 10,
624c4762a1bSJed Brown                                           10, 14, 15, 11,
625c4762a1bSJed Brown                                           11, 15, 16, 12,
626c4762a1bSJed Brown                                           13, 17, 18, 14,
627c4762a1bSJed Brown                                           14, 18, 19, 15,
628c4762a1bSJed Brown                                           15, 19, 20, 16,
629c4762a1bSJed Brown                                           17, 21, 22, 18,
630c4762a1bSJed Brown                                           18, 22, 23, 19,
631c4762a1bSJed Brown                                           19, 23, 24, 20};
632c4762a1bSJed 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};
633c4762a1bSJed 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,
634c4762a1bSJed 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};
6355b9dfbb6SMatthew G. Knepley       PetscInt    faultPoints[4]       = {11, 15, 19, 20};
636c4762a1bSJed Brown       PetscInt    fault2Points[2]      = {17, 18};
637c4762a1bSJed Brown 
6389566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6395b9dfbb6SMatthew G. Knepley       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "fault",  faultPoints[p], 1));
6405b9dfbb6SMatthew G. Knepley       for (p = 3; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "faultBd", faultPoints[p], 1));
6419566063dSJacob Faibussowitsch       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault2", fault2Points[p], 1));
6429566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6439566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
6449566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
6459566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 1));
6469566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 4, 1));
6479566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
6489566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 6, 2));
6499566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 7, 2));
6509566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 8, 2));
651c4762a1bSJed Brown     }
652c4762a1bSJed Brown     break;
653c4762a1bSJed Brown     default:
65463a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
655c4762a1bSJed Brown     }
656c4762a1bSJed Brown   } else {
657c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
658c4762a1bSJed Brown 
6599566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
6605b9dfbb6SMatthew G. Knepley     if (testNum == 2 || testNum == 3) PetscCall(DMCreateLabel(*dm, "pfault"));
6619566063dSJacob Faibussowitsch     else                              PetscCall(DMCreateLabel(*dm, "fault"));
662c4762a1bSJed Brown   }
6639566063dSJacob Faibussowitsch   PetscCall(DMPlexInterpolate(*dm, &idm));
6649566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
6659566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
666c4762a1bSJed Brown   *dm  = idm;
667c4762a1bSJed Brown   PetscFunctionReturn(0);
668c4762a1bSJed Brown }
669c4762a1bSJed Brown 
670b253942bSMatthew G. Knepley static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
671c4762a1bSJed Brown {
672c4762a1bSJed Brown   DM             idm;
673c4762a1bSJed Brown   PetscInt       depth = 3, p;
674c4762a1bSJed Brown   PetscMPIInt    rank;
675c4762a1bSJed Brown 
676c4762a1bSJed Brown   PetscFunctionBegin;
6779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
678dd400576SPatrick Sanan   if (rank == 0) {
679c4762a1bSJed Brown     switch (testNum) {
680c4762a1bSJed Brown     case 0:
681c4762a1bSJed Brown     {
682c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
683c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8,8, 0,0,0,0,0,0,0,0,0,0,0,0};
684c4762a1bSJed Brown       PetscInt    cones[16]            = {2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8};
685c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0,0,0,0,0,0,0,0,  0,0, 0,0 ,0, 0, 0,0};
686c4762a1bSJed 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,
687c4762a1bSJed 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,
688c4762a1bSJed 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};
689c4762a1bSJed Brown       PetscInt    markerPoints[52]     = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
690c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {3, 4, 7, 8};
691c4762a1bSJed Brown 
6929566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6939566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
6949566063dSJacob Faibussowitsch       for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(idm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
6959566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
6969566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6979566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
698c4762a1bSJed Brown     }
699c4762a1bSJed Brown     break;
700c4762a1bSJed Brown     case 1:
701c4762a1bSJed Brown     {
702c4762a1bSJed Brown       /* Cell Adjacency Graph:
703c4762a1bSJed Brown         0 -- { 8, 13, 21, 24} --> 1
704c4762a1bSJed Brown         0 -- {20, 21, 23, 24} --> 5 F
705c4762a1bSJed Brown         1 -- {10, 15, 21, 24} --> 2
706c4762a1bSJed Brown         1 -- {13, 14, 15, 24} --> 6
707c4762a1bSJed Brown         2 -- {21, 22, 24, 25} --> 4 F
708c4762a1bSJed Brown         3 -- {21, 24, 30, 35} --> 4
709c4762a1bSJed Brown         3 -- {21, 24, 28, 33} --> 5
710c4762a1bSJed Brown        */
711c4762a1bSJed Brown       PetscInt    numPoints[2]         = {30, 7};
712c4762a1bSJed 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};
713c4762a1bSJed Brown       PetscInt    cones[56]            = { 8, 21, 20,  7, 13, 12, 23, 24,
714c4762a1bSJed Brown                                           14, 15, 10,  9, 13,  8, 21, 24,
715c4762a1bSJed Brown                                           15, 16, 11, 10, 24, 21, 22, 25,
716c4762a1bSJed Brown                                           30, 29, 28, 21, 35, 24, 33, 34,
717c4762a1bSJed Brown                                           24, 21, 30, 35, 25, 36, 31, 22,
718c4762a1bSJed Brown                                           27, 20, 21, 28, 32, 33, 24, 23,
719c4762a1bSJed Brown                                           15, 24, 13, 14, 19, 18, 17, 26};
720c4762a1bSJed 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};
721c4762a1bSJed 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,
722c4762a1bSJed 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,
723c4762a1bSJed 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,
724c4762a1bSJed 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,
725c4762a1bSJed 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};
726c4762a1bSJed Brown       PetscInt    faultPoints[6]       = {20, 21, 22, 23, 24, 25};
727c4762a1bSJed Brown 
7289566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
7299566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
7309566063dSJacob Faibussowitsch       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
7319566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
7329566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
7339566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
7349566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
7359566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 4, 2));
7369566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
7379566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 6, 2));
738c4762a1bSJed Brown     }
739c4762a1bSJed Brown     break;
740c4762a1bSJed Brown     case 2:
741c4762a1bSJed Brown     {
742c4762a1bSJed Brown       /* Buried fault edge */
743c4762a1bSJed Brown       PetscInt    numPoints[2]         = {18, 4};
744c4762a1bSJed 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};
745c4762a1bSJed Brown       PetscInt    cones[32]            = { 4,  5,  8,  7, 13, 16, 17, 14,
746c4762a1bSJed Brown                                            5,  6,  9,  8, 14, 17, 18, 15,
747c4762a1bSJed Brown                                            7,  8, 11, 10, 16, 19, 20, 17,
748c4762a1bSJed Brown                                            8,  9, 12, 11, 17, 20, 21, 18};
749c4762a1bSJed 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};
750c4762a1bSJed 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,
751c4762a1bSJed 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,
752c4762a1bSJed 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};
753c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {7, 8, 16, 17};
754c4762a1bSJed Brown 
7559566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
7569566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
7579566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
7589566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
7599566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
7609566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 2));
7619566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
762c4762a1bSJed Brown     }
763c4762a1bSJed Brown     break;
76463a3b9bcSJacob Faibussowitsch     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
765c4762a1bSJed Brown     }
766c4762a1bSJed Brown   } else {
767c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
768c4762a1bSJed Brown 
7699566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
7709566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
7719566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(idm, "fault"));
772c4762a1bSJed Brown   }
7739566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
7749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
775c4762a1bSJed Brown   *dm  = idm;
776c4762a1bSJed Brown   PetscFunctionReturn(0);
777c4762a1bSJed Brown }
778c4762a1bSJed Brown 
779b253942bSMatthew G. Knepley static PetscErrorCode CreateFaultLabel(DM dm)
780b253942bSMatthew G. Knepley {
781b253942bSMatthew G. Knepley   DMLabel        label;
782b253942bSMatthew G. Knepley   PetscInt       dim, h, pStart, pEnd, pMax, p;
783b253942bSMatthew G. Knepley 
784b253942bSMatthew G. Knepley   PetscFunctionBegin;
7859566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7869566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "cohesive"));
7879566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &label));
788b253942bSMatthew G. Knepley   for (h = 0; h <= dim; ++h) {
7899566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(dm, h, NULL, &pMax));
7909566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &pStart, &pEnd));
7919566063dSJacob Faibussowitsch     for (p = pMax; p < pEnd; ++p) PetscCall(DMLabelSetValue(label, p, 1));
792b253942bSMatthew G. Knepley   }
793b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
794b253942bSMatthew G. Knepley }
795b253942bSMatthew G. Knepley 
796b253942bSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *user)
797b253942bSMatthew G. Knepley {
798b253942bSMatthew G. Knepley   PetscFE        fe;
799b253942bSMatthew G. Knepley   DMLabel        fault;
800b7519becSMatthew G. Knepley   PetscInt       dim, Ncf = user->cohesiveFields, f;
801b253942bSMatthew G. Knepley 
802b253942bSMatthew G. Knepley   PetscFunctionBegin;
8039566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
8049566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &fault));
8059566063dSJacob Faibussowitsch   PetscCall(DMLabelView(fault, PETSC_VIEWER_STDOUT_WORLD));
806b253942bSMatthew G. Knepley 
8079566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, user->cellSimplex, "displacement_", PETSC_DETERMINE, &fe));
8089566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(fe, "displacement"));
8099566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject) fe));
8109566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
811b253942bSMatthew G. Knepley 
812b7519becSMatthew G. Knepley   if (Ncf > 0) {
8139566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim-1, dim, user->cellSimplex, "faulttraction_", PETSC_DETERMINE, &fe));
8149566063dSJacob Faibussowitsch     PetscCall(PetscFESetName(fe, "fault traction"));
8159566063dSJacob Faibussowitsch     PetscCall(DMAddField(dm, fault, (PetscObject) fe));
8169566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
817b7519becSMatthew G. Knepley   }
818b7519becSMatthew G. Knepley   for (f = 1; f < Ncf; ++f) {
819b7519becSMatthew G. Knepley     char name[256], opt[256];
820b7519becSMatthew G. Knepley 
82163a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, 256, "fault field %" PetscInt_FMT, f));
82263a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(opt,  256, "faultfield_%" PetscInt_FMT "_", f));
8239566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim-1, dim, user->cellSimplex, opt, PETSC_DETERMINE, &fe));
8249566063dSJacob Faibussowitsch     PetscCall(PetscFESetName(fe, name));
8259566063dSJacob Faibussowitsch     PetscCall(DMAddField(dm, fault, (PetscObject) fe));
8269566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
827b7519becSMatthew G. Knepley   }
828b253942bSMatthew G. Knepley 
8299566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
830b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
831b253942bSMatthew G. Knepley }
832b253942bSMatthew G. Knepley 
833b253942bSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
834c4762a1bSJed Brown {
835c4762a1bSJed Brown   PetscInt       dim          = user->dim;
836c4762a1bSJed Brown   PetscBool      cellSimplex  = user->cellSimplex, hasFault, hasFault2, hasParallelFault;
837c4762a1bSJed Brown   PetscMPIInt    rank, size;
838dd96b1f9SToby Isaac   DMLabel        matLabel;
839c4762a1bSJed Brown 
840c4762a1bSJed Brown   PetscFunctionBegin;
8419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
8439566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
8449566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
8459566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
846c4762a1bSJed Brown   switch (dim) {
847c4762a1bSJed Brown   case 2:
848c4762a1bSJed Brown     if (cellSimplex) {
8499566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_2D(comm, user->testNum, dm));
850c4762a1bSJed Brown     } else {
8519566063dSJacob Faibussowitsch       PetscCall(CreateQuad_2D(comm, user->testNum, dm));
852c4762a1bSJed Brown     }
853c4762a1bSJed Brown     break;
854c4762a1bSJed Brown   case 3:
855c4762a1bSJed Brown     if (cellSimplex) {
8569566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_3D(comm, user, *dm));
857c4762a1bSJed Brown     } else {
8589566063dSJacob Faibussowitsch       PetscCall(CreateHex_3D(comm, user->testNum, dm));
859c4762a1bSJed Brown     }
860c4762a1bSJed Brown     break;
861c4762a1bSJed Brown   default:
86263a3b9bcSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %" PetscInt_FMT, dim);
863c4762a1bSJed Brown   }
8649566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "orig_"));
8659566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
8669566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
8679566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dm, "material", &matLabel));
868dd96b1f9SToby Isaac   if (matLabel) {
8699566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelComplete(*dm, matLabel));
870dd96b1f9SToby Isaac   }
8719566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
8729566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "fault", &hasFault));
873c4762a1bSJed Brown   if (hasFault) {
874b253942bSMatthew G. Knepley     DM      dmHybrid = NULL, dmInterface = NULL;
875b253942bSMatthew G. Knepley     DMLabel faultLabel, faultBdLabel, hybridLabel, splitLabel;
876c4762a1bSJed Brown 
8779566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
8789566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "faultBd", &faultBdLabel));
879caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, &splitLabel, &dmInterface, &dmHybrid));
8809566063dSJacob Faibussowitsch     PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD));
8819566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
8829566063dSJacob Faibussowitsch     PetscCall(DMLabelView(splitLabel, PETSC_VIEWER_STDOUT_WORLD));
8839566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&splitLabel));
8849566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_interface_view"));
8859566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmInterface));
8869566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
887c4762a1bSJed Brown     *dm  = dmHybrid;
888c4762a1bSJed Brown   }
8899566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "fault2", &hasFault2));
890c4762a1bSJed Brown   if (hasFault2) {
891c4762a1bSJed Brown     DM      dmHybrid = NULL;
892c4762a1bSJed Brown     DMLabel faultLabel, faultBdLabel, hybridLabel;
893c4762a1bSJed Brown 
8949566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "faulted_"));
8959566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
8969566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
8979566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
8989566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
8999566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault2", &faultLabel));
9009566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault2Bd", &faultBdLabel));
901caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, NULL, &dmHybrid));
9029566063dSJacob Faibussowitsch     PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD));
9039566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
9049566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
905c4762a1bSJed Brown     *dm  = dmHybrid;
906c4762a1bSJed Brown   }
907c4762a1bSJed Brown   if (user->testPartition && size > 1) {
908c4762a1bSJed Brown     PetscPartitioner part;
909c4762a1bSJed Brown     PetscInt *sizes  = NULL;
910c4762a1bSJed Brown     PetscInt *points = NULL;
911c4762a1bSJed Brown 
912dd400576SPatrick Sanan     if (rank == 0) {
913c4762a1bSJed Brown       if (dim == 2 && cellSimplex && size == 2) {
914c4762a1bSJed Brown         switch (user->testNum) {
915c4762a1bSJed Brown         case 0: {
916c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {1, 2};
917c4762a1bSJed Brown           PetscInt triPoints_p2[3] = {0, 1, 2};
918c4762a1bSJed Brown 
9199566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9209566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  triSizes_p2, 2));
9219566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, triPoints_p2, 3));break;}
9225b9dfbb6SMatthew G. Knepley         case 3: {
9235b9dfbb6SMatthew G. Knepley           PetscInt triSizes_p2[2]  = {3, 3};
9245b9dfbb6SMatthew G. Knepley           PetscInt triPoints_p2[6] = {0, 1, 2,  3, 4, 5};
9255b9dfbb6SMatthew G. Knepley 
9265b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(2, &sizes, 6, &points));
9275b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes,  triSizes_p2, 2));
9285b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(points, triPoints_p2, 6));break;}
929c4762a1bSJed Brown         default:
9305b9dfbb6SMatthew 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);
9315b9dfbb6SMatthew G. Knepley         }
9325b9dfbb6SMatthew G. Knepley       } else if (dim == 2 && cellSimplex && size == 6) {
9335b9dfbb6SMatthew G. Knepley         switch (user->testNum) {
9345b9dfbb6SMatthew G. Knepley         case 4: {
9355b9dfbb6SMatthew G. Knepley           PetscInt triSizes_p6[6]  = {1, 1, 1, 1, 1, 1};
9365b9dfbb6SMatthew G. Knepley           PetscInt triPoints_p6[6] = {0, 1, 2, 3, 4, 5};
9375b9dfbb6SMatthew G. Knepley 
9385b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(6, &sizes, 6, &points));
9395b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes,  triSizes_p6, 6));
9405b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(points, triPoints_p6, 6));break;}
9415b9dfbb6SMatthew G. Knepley         default:
9425b9dfbb6SMatthew 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);
943c4762a1bSJed Brown         }
944c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && size == 2) {
945c4762a1bSJed Brown         switch (user->testNum) {
946c4762a1bSJed Brown         case 0: {
947c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 2};
948c4762a1bSJed Brown           PetscInt quadPoints_p2[3] = {0, 1, 2};
949c4762a1bSJed Brown 
9509566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9519566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  quadSizes_p2, 2));
9529566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, quadPoints_p2, 3));break;}
953c4762a1bSJed Brown         case 2: {
954c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 1};
955c4762a1bSJed Brown           PetscInt quadPoints_p2[2] = {0, 1};
956c4762a1bSJed Brown 
9579566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
9589566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  quadSizes_p2, 2));
9599566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));break;}
9605b9dfbb6SMatthew G. Knepley         case 3: {
9615b9dfbb6SMatthew G. Knepley           PetscInt quadSizes_p2[2]  = {1, 1};
9625b9dfbb6SMatthew G. Knepley           PetscInt quadPoints_p2[2] = {1, 0};
9635b9dfbb6SMatthew G. Knepley 
9645b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
9655b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes,  quadSizes_p2, 2));
9665b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));break;}
967c4762a1bSJed Brown         default:
9685b9dfbb6SMatthew 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);
969c4762a1bSJed Brown         }
970c4762a1bSJed Brown       } else if (dim == 3 && cellSimplex && size == 2) {
971c4762a1bSJed Brown         switch (user->testNum) {
972c4762a1bSJed Brown         case 0: {
973c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 2};
974c4762a1bSJed Brown           PetscInt tetPoints_p2[3] = {0, 1, 2};
975c4762a1bSJed Brown 
9769566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9779566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  tetSizes_p2, 2));
9789566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, tetPoints_p2, 3));break;}
979c4762a1bSJed Brown         default:
9805b9dfbb6SMatthew 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);
981c4762a1bSJed Brown         }
982c4762a1bSJed Brown       } else if (dim == 3 && !cellSimplex && size == 2) {
983c4762a1bSJed Brown         switch (user->testNum) {
984c4762a1bSJed Brown         case 0: {
985c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {1, 2};
986c4762a1bSJed Brown           PetscInt hexPoints_p2[3] = {0, 1, 2};
987c4762a1bSJed Brown 
9889566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9899566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  hexSizes_p2, 2));
9909566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, hexPoints_p2, 3));break;}
991c4762a1bSJed Brown         default:
9925b9dfbb6SMatthew 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);
993c4762a1bSJed Brown         }
994c4762a1bSJed Brown       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
995c4762a1bSJed Brown     }
9969566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
9979566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
9989566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
9999566063dSJacob Faibussowitsch     PetscCall(PetscFree2(sizes, points));
1000c4762a1bSJed Brown   }
1001c4762a1bSJed Brown   {
1002c4762a1bSJed Brown     DM pdm = NULL;
1003c4762a1bSJed Brown 
1004c4762a1bSJed Brown     /* Distribute mesh over processes */
10059566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
1006c4762a1bSJed Brown     if (pdm) {
10079566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
10089566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
1009c4762a1bSJed Brown       *dm  = pdm;
1010c4762a1bSJed Brown     }
1011c4762a1bSJed Brown   }
10129566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "pfault", &hasParallelFault));
1013c4762a1bSJed Brown   if (hasParallelFault) {
10145b9dfbb6SMatthew G. Knepley     DM      dmHybrid = NULL, dmInterface;
1015c4762a1bSJed Brown     DMLabel faultLabel, faultBdLabel, hybridLabel;
1016c4762a1bSJed Brown 
10179566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "pfault", &faultLabel));
10189566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "pfaultBd", &faultBdLabel));
1019caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, &dmInterface, &dmHybrid));
10205b9dfbb6SMatthew G. Knepley     PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_fault_view"));
10215b9dfbb6SMatthew G. Knepley     {
10225b9dfbb6SMatthew G. Knepley       PetscViewer viewer;
10235b9dfbb6SMatthew G. Knepley       PetscMPIInt rank;
10245b9dfbb6SMatthew G. Knepley 
10255b9dfbb6SMatthew G. Knepley       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank));
10265b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
10275b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "Rank %d\n", rank));
10285b9dfbb6SMatthew G. Knepley       PetscCall(DMLabelView(hybridLabel, viewer));
10295b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
10305b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
10315b9dfbb6SMatthew G. Knepley     }
10329566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
10335b9dfbb6SMatthew G. Knepley     PetscCall(DMDestroy(&dmInterface));
10349566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
1035c4762a1bSJed Brown     *dm  = dmHybrid;
1036c4762a1bSJed Brown   }
10379566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, "Hybrid Mesh"));
10389566063dSJacob Faibussowitsch   PetscCall(CreateFaultLabel(*dm));
10399566063dSJacob Faibussowitsch   PetscCall(CreateDiscretization(*dm, user));
10409566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
10419566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
10429566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
10439566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1044c4762a1bSJed Brown   PetscFunctionReturn(0);
1045c4762a1bSJed Brown }
1046c4762a1bSJed Brown 
1047b253942bSMatthew G. Knepley static PetscErrorCode TestMesh(DM dm, AppCtx *user)
1048b253942bSMatthew G. Knepley {
1049b253942bSMatthew G. Knepley   PetscFunctionBegin;
10509566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckSymmetry(dm));
10519566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckSkeleton(dm, 0));
10529566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckFaces(dm, 0));
1053b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
1054b253942bSMatthew G. Knepley }
1055b253942bSMatthew G. Knepley 
1056b253942bSMatthew G. Knepley static PetscErrorCode TestDiscretization(DM dm, AppCtx *user)
1057b253942bSMatthew G. Knepley {
1058b253942bSMatthew G. Knepley   PetscSection   s;
1059b253942bSMatthew G. Knepley 
1060b253942bSMatthew G. Knepley   PetscFunctionBegin;
10619566063dSJacob Faibussowitsch   PetscCall(DMGetSection(dm, &s));
10629566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject) s, NULL, "-local_section_view"));
1063b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
1064b253942bSMatthew G. Knepley }
1065b253942bSMatthew G. Knepley 
1066b253942bSMatthew G. Knepley static PetscErrorCode r(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1067b253942bSMatthew G. Knepley {
1068b253942bSMatthew G. Knepley   PetscInt d;
1069b253942bSMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = x[d];
1070b253942bSMatthew G. Knepley   return 0;
1071b253942bSMatthew G. Knepley }
1072b253942bSMatthew G. Knepley 
1073b253942bSMatthew G. Knepley static PetscErrorCode rp1(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1074b253942bSMatthew G. Knepley {
1075b253942bSMatthew G. Knepley   PetscInt d;
1076b253942bSMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = x[d] + (d > 0 ? 1.0 : 0.0);
1077b253942bSMatthew G. Knepley   return 0;
1078b253942bSMatthew G. Knepley }
1079b253942bSMatthew G. Knepley 
1080b253942bSMatthew G. Knepley static PetscErrorCode phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1081b253942bSMatthew G. Knepley {
1082b253942bSMatthew G. Knepley   PetscInt d;
1083b253942bSMatthew G. Knepley   u[0] = -x[1];
1084b253942bSMatthew G. Knepley   u[1] =  x[0];
1085b253942bSMatthew G. Knepley   for (d = 2; d < dim; ++d) u[d] = x[d];
1086b253942bSMatthew G. Knepley   return 0;
1087b253942bSMatthew G. Knepley }
1088b253942bSMatthew G. Knepley 
1089b253942bSMatthew G. Knepley /* \lambda \cdot (\psi_u^- - \psi_u^+) */
1090b253942bSMatthew G. Knepley static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1091b253942bSMatthew G. Knepley                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1092b253942bSMatthew G. Knepley                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1093b253942bSMatthew G. Knepley                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1094b253942bSMatthew G. Knepley {
1095b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[1]-uOff[0];
1096b253942bSMatthew G. Knepley   PetscInt       c;
1097b253942bSMatthew G. Knepley   for (c = 0;  c < Nc;   ++c) {
1098b253942bSMatthew G. Knepley     f0[c]    =   u[Nc*2+c] + x[Nc-c-1];
1099b253942bSMatthew G. Knepley     f0[Nc+c] = -(u[Nc*2+c] + x[Nc-c-1]);
1100b253942bSMatthew G. Knepley   }
1101b253942bSMatthew G. Knepley }
1102b253942bSMatthew G. Knepley 
1103b253942bSMatthew G. Knepley /* (d - u^+ + u^-) \cdot \psi_\lambda */
1104b253942bSMatthew G. Knepley static void f0_bd_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1105b253942bSMatthew G. Knepley                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1106b253942bSMatthew G. Knepley                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1107b253942bSMatthew G. Knepley                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1108b253942bSMatthew G. Knepley {
1109b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[2]-uOff[1];
1110b253942bSMatthew G. Knepley   PetscInt       c;
1111b253942bSMatthew G. Knepley 
1112b253942bSMatthew G. Knepley   for (c = 0; c < Nc; ++c) f0[c] = (c > 0 ? 1.0 : 0.0) + u[c] - u[Nc+c];
1113b253942bSMatthew G. Knepley }
1114b253942bSMatthew G. Knepley 
1115b253942bSMatthew G. Knepley /* \psi_lambda \cdot (\psi_u^- - \psi_u^+) */
1116b253942bSMatthew G. Knepley static void g0_bd_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1117b253942bSMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1118b253942bSMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1119b253942bSMatthew G. Knepley                      PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1120b253942bSMatthew G. Knepley {
1121b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[1]-uOff[0];
1122b253942bSMatthew G. Knepley   PetscInt       c;
1123b253942bSMatthew G. Knepley 
1124b253942bSMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
1125b253942bSMatthew G. Knepley     g0[(0 +c)*Nc+c] =  1.0;
1126b253942bSMatthew G. Knepley     g0[(Nc+c)*Nc+c] = -1.0;
1127b253942bSMatthew G. Knepley   }
1128b253942bSMatthew G. Knepley }
1129b253942bSMatthew G. Knepley 
1130b253942bSMatthew G. Knepley /* (-\psi_u^+ + \psi_u^-) \cdot \psi_\lambda */
1131b253942bSMatthew G. Knepley static void g0_bd_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1132b253942bSMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1133b253942bSMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1134b253942bSMatthew G. Knepley                      PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1135b253942bSMatthew G. Knepley {
1136b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[2]-uOff[1];
1137b253942bSMatthew G. Knepley   PetscInt       c;
1138b253942bSMatthew G. Knepley 
1139b253942bSMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
1140b253942bSMatthew G. Knepley     g0[c*Nc*2+c]    =  1.0;
1141b253942bSMatthew G. Knepley     g0[c*Nc*2+Nc+c] = -1.0;
1142b253942bSMatthew G. Knepley   }
1143b253942bSMatthew G. Knepley }
1144b253942bSMatthew G. Knepley 
1145b253942bSMatthew G. Knepley static PetscErrorCode TestAssembly(DM dm, AppCtx *user)
1146b253942bSMatthew G. Knepley {
1147b253942bSMatthew G. Knepley   Mat              J;
1148b253942bSMatthew G. Knepley   Vec              locX, locF;
1149b253942bSMatthew G. Knepley   PetscDS          probh;
1150b253942bSMatthew G. Knepley   DMLabel          fault, material;
1151b253942bSMatthew G. Knepley   IS               cohesiveCells;
1152b7519becSMatthew G. Knepley   PetscWeakForm    wf;
115306ad1575SMatthew G. Knepley   PetscFormKey     keys[3];
1154b253942bSMatthew G. Knepley   PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
1155b253942bSMatthew G. Knepley   PetscInt         dim, Nf, cMax, cEnd, id;
1156b7519becSMatthew G. Knepley   PetscMPIInt      rank;
1157b253942bSMatthew G. Knepley 
1158b253942bSMatthew G. Knepley   PetscFunctionBegin;
11599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank));
11609566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
11619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, NULL, &cMax));
11629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
11639566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, cEnd - cMax, cMax, 1, &cohesiveCells));
11649566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &fault));
11659566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
11669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) locX, "Local Solution"));
11679566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
11689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) locF, "Local Residual"));
11699566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
11709566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) J, "Jacobian"));
1171b253942bSMatthew G. Knepley 
1172b253942bSMatthew G. Knepley   /* The initial guess has displacement shifted by one unit in each fault parallel direction across the fault */
11739566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "material", &material));
1174b253942bSMatthew G. Knepley   id   = 1;
1175b253942bSMatthew G. Knepley   initialGuess[0] = r;
1176b253942bSMatthew G. Knepley   initialGuess[1] = NULL;
11779566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
1178b253942bSMatthew G. Knepley   id   = 2;
1179b253942bSMatthew G. Knepley   initialGuess[0] = rp1;
1180b253942bSMatthew G. Knepley   initialGuess[1] = NULL;
11819566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
1182b253942bSMatthew G. Knepley   id   = 1;
1183b253942bSMatthew G. Knepley   initialGuess[0] = NULL;
1184b253942bSMatthew G. Knepley   initialGuess[1] = phi;
11859566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, fault, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
11869566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(locX, NULL, "-local_solution_view"));
1187b253942bSMatthew G. Knepley 
11889566063dSJacob Faibussowitsch   PetscCall(DMGetCellDS(dm, cMax, &probh));
11899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(probh, &wf));
11909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(probh, &Nf));
11919566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 1, 0, 0, 0, f0_bd_u, 0, NULL));
11929566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 2, 0, 0, 0, f0_bd_u, 0, NULL));
11939566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 1, 0, 1, 0, 0, g0_bd_ul, 0, NULL, 0, NULL, 0, NULL));
11949566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 2, 0, 1, 0, 0, g0_bd_ul, 0, NULL, 0, NULL, 0, NULL));
1195b7519becSMatthew G. Knepley   if (Nf > 1) {
11969566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdResidual(wf, fault, 1, 1, 0, 0, f0_bd_l, 0, NULL));
11979566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdJacobian(wf, fault, 1, 1, 0, 0, 0, g0_bd_lu, 0, NULL, 0, NULL, 0, NULL));
1198b7519becSMatthew G. Knepley   }
1199c5853193SPierre Jolivet   if (rank == 0) PetscCall(PetscDSView(probh, NULL));
1200b253942bSMatthew G. Knepley 
1201dd96b1f9SToby Isaac   keys[0].label = NULL;
1202dd96b1f9SToby Isaac   keys[0].value = 0;
12036528b96dSMatthew G. Knepley   keys[0].field = 0;
1204dd96b1f9SToby Isaac   keys[0].part  = 0;
1205b7519becSMatthew G. Knepley   keys[1].label = material;
1206b7519becSMatthew G. Knepley   keys[1].value = 2;
12076528b96dSMatthew G. Knepley   keys[1].field = 0;
1208dd96b1f9SToby Isaac   keys[1].part  = 0;
1209b7519becSMatthew G. Knepley   keys[2].label = fault;
1210b7519becSMatthew G. Knepley   keys[2].value = 1;
1211b7519becSMatthew G. Knepley   keys[2].field = 1;
1212dd96b1f9SToby Isaac   keys[2].part  = 0;
12139566063dSJacob Faibussowitsch   PetscCall(VecSet(locF, 0.));
12149566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeResidual_Hybrid_Internal(dm, keys, cohesiveCells, 0.0, locX, NULL, 0.0, locF, user));
12159566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(locF, NULL, "-local_residual_view"));
12169566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(J));
12179566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeJacobian_Hybrid_Internal(dm, keys, cohesiveCells, 0.0, 0.0, locX, NULL, J, J, user));
1218*49664dceSMatthew G. Knepley   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1219*49664dceSMatthew G. Knepley   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
12209566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-local_jacobian_view"));
1221b253942bSMatthew G. Knepley 
12229566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
12239566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
12249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
12259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cohesiveCells));
1226b253942bSMatthew G. Knepley   PetscFunctionReturn(0);
1227b253942bSMatthew G. Knepley }
1228b253942bSMatthew G. Knepley 
1229c4762a1bSJed Brown int main(int argc, char **argv)
1230c4762a1bSJed Brown {
1231c4762a1bSJed Brown   DM             dm;
1232c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
1233c4762a1bSJed Brown 
12349566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
12359566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
12369566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
12379566063dSJacob Faibussowitsch   PetscCall(TestMesh(dm, &user));
12389566063dSJacob Faibussowitsch   PetscCall(TestDiscretization(dm, &user));
12399566063dSJacob Faibussowitsch   PetscCall(TestAssembly(dm, &user));
12409566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
12419566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1242b122ec5aSJacob Faibussowitsch   return 0;
1243c4762a1bSJed Brown }
1244c4762a1bSJed Brown 
1245c4762a1bSJed Brown /*TEST
1246c4762a1bSJed Brown   testset:
1247b253942bSMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1248b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1 -local_section_view \
1249b253942bSMatthew G. Knepley           -local_solution_view -local_residual_view -local_jacobian_view
1250c4762a1bSJed Brown     test:
1251c4762a1bSJed Brown       suffix: tri_0
1252c4762a1bSJed Brown       args: -dim 2
12531c6715b8SMatthew 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"
1254c4762a1bSJed Brown     test:
1255c4762a1bSJed Brown       suffix: tri_t1_0
1256c4762a1bSJed Brown       args: -dim 2 -test_num 1
12571c6715b8SMatthew 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"
1258c4762a1bSJed Brown     test:
12595b9dfbb6SMatthew G. Knepley       suffix: tri_t2_0
12605b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 2
12615b9dfbb6SMatthew 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"
12625b9dfbb6SMatthew G. Knepley     test:
12633bdf08b3SMatthew G. Knepley       suffix: tri_t5_0
12643bdf08b3SMatthew G. Knepley       args: -dim 2 -test_num 5
12653bdf08b3SMatthew 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"
12663bdf08b3SMatthew G. Knepley     test:
1267c4762a1bSJed Brown       suffix: tet_0
1268c4762a1bSJed Brown       args: -dim 3
12691c6715b8SMatthew 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"
1270c4762a1bSJed Brown     test:
1271b253942bSMatthew G. Knepley       suffix: tet_t1_0
1272b253942bSMatthew G. Knepley       args: -dim 3 -test_num 1
12731c6715b8SMatthew 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"
1274b253942bSMatthew G. Knepley 
1275b253942bSMatthew G. Knepley   testset:
1276b253942bSMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1277b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1
1278b253942bSMatthew G. Knepley     test:
1279c4762a1bSJed Brown       suffix: tet_1
1280c4762a1bSJed Brown       nsize: 2
1281c4762a1bSJed Brown       args: -dim 3
12821c6715b8SMatthew 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"
1283c4762a1bSJed Brown     test:
1284b253942bSMatthew G. Knepley       suffix: tri_1
1285b253942bSMatthew G. Knepley       nsize: 2
1286b253942bSMatthew G. Knepley       args: -dim 2
12871c6715b8SMatthew 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"
12885b9dfbb6SMatthew G. Knepley     test:
12895b9dfbb6SMatthew G. Knepley       suffix: tri_t3_0
12905b9dfbb6SMatthew G. Knepley       nsize: 2
12915b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 3
12925b9dfbb6SMatthew 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"
12935b9dfbb6SMatthew G. Knepley     test:
12945b9dfbb6SMatthew G. Knepley       suffix: tri_t4_0
12955b9dfbb6SMatthew G. Knepley       nsize: 6
12965b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 4
12975b9dfbb6SMatthew 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"
1298c4762a1bSJed Brown 
1299c4762a1bSJed Brown   testset:
1300ecfb78b5SMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1301b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1
1302c4762a1bSJed Brown     # 2D Quads
1303c4762a1bSJed Brown     test:
1304c4762a1bSJed Brown       suffix: quad_0
1305c4762a1bSJed Brown       args: -dim 2 -cell_simplex 0
13061c6715b8SMatthew 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"
1307c4762a1bSJed Brown     test:
1308c4762a1bSJed Brown       suffix: quad_1
1309c4762a1bSJed Brown       nsize: 2
1310c4762a1bSJed Brown       args: -dim 2 -cell_simplex 0
13111c6715b8SMatthew 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"
1312c4762a1bSJed Brown     test:
1313c4762a1bSJed Brown       suffix: quad_t1_0
131454fcfd0cSMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 1 -faulted_dm_plex_check_all
13151c6715b8SMatthew 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"
13165b9dfbb6SMatthew G. Knepley     test:
13175b9dfbb6SMatthew G. Knepley       suffix: quad_t2_0
13185b9dfbb6SMatthew G. Knepley       nsize: 2
13195b9dfbb6SMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 2
13205b9dfbb6SMatthew 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"
13215b9dfbb6SMatthew G. Knepley     test:
13225b9dfbb6SMatthew G. Knepley       # TODO: The PetscSF is wrong here (connects to wrong side of split)
13235b9dfbb6SMatthew G. Knepley       suffix: quad_t3_0
13245b9dfbb6SMatthew G. Knepley       nsize: 2
13255b9dfbb6SMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 3
13265b9dfbb6SMatthew 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"
1327c4762a1bSJed Brown     # 3D Hex
1328c4762a1bSJed Brown     test:
1329c4762a1bSJed Brown       suffix: hex_0
1330c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0
13311c6715b8SMatthew 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"
1332c4762a1bSJed Brown     test:
1333c4762a1bSJed Brown       suffix: hex_1
1334c4762a1bSJed Brown       nsize: 2
1335c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0
13361c6715b8SMatthew 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"
1337c4762a1bSJed Brown     test:
1338c4762a1bSJed Brown       suffix: hex_t1_0
1339c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0 -test_num 1
13401c6715b8SMatthew 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"
1341c4762a1bSJed Brown     test:
1342c4762a1bSJed Brown       suffix: hex_t2_0
1343c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0 -test_num 2
13441c6715b8SMatthew 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"
1345c4762a1bSJed Brown 
1346c4762a1bSJed Brown TEST*/
1347