xref: /petsc/src/dm/impls/plex/tests/ex5.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 */
39767de2c8aSMatthew G. Knepley   PetscBool testAssembly;   // Flag for assembly test
398c4762a1bSJed Brown   PetscInt  testNum;        /* The particular mesh to test */
399b7519becSMatthew G. Knepley   PetscInt  cohesiveFields; /* The number of cohesive fields */
400c4762a1bSJed Brown } AppCtx;
401c4762a1bSJed Brown 
ProcessOptions(MPI_Comm comm,AppCtx * options)402d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
403d71ae5a4SJacob Faibussowitsch {
404c4762a1bSJed Brown   PetscFunctionBegin;
405c4762a1bSJed Brown   options->debug          = 0;
406c4762a1bSJed Brown   options->dim            = 2;
407c4762a1bSJed Brown   options->cellSimplex    = PETSC_TRUE;
408c4762a1bSJed Brown   options->testPartition  = PETSC_TRUE;
40967de2c8aSMatthew G. Knepley   options->testAssembly   = PETSC_TRUE;
410c4762a1bSJed Brown   options->testNum        = 0;
411b7519becSMatthew G. Knepley   options->cohesiveFields = 1;
412c4762a1bSJed Brown 
413d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
4149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex5.c", options->debug, &options->debug, NULL, 0));
4159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex5.c", options->dim, &options->dim, NULL, 1, 3));
4169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex5.c", options->cellSimplex, &options->cellSimplex, NULL));
4179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex5.c", options->testPartition, &options->testPartition, NULL));
41867de2c8aSMatthew G. Knepley   PetscCall(PetscOptionsBool("-test_assembly", "Run the assembly test", "ex5.c", options->testAssembly, &options->testAssembly, NULL));
4199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex5.c", options->testNum, &options->testNum, NULL, 0));
4209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-cohesive_fields", "The number of cohesive fields", "ex5.c", options->cohesiveFields, &options->cohesiveFields, NULL, 0));
421d0609cedSBarry Smith   PetscOptionsEnd();
4223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
423c4762a1bSJed Brown }
424c4762a1bSJed Brown 
CreateSimplex_2D(MPI_Comm comm,PetscInt testNum,DM * dm)425d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
426d71ae5a4SJacob Faibussowitsch {
427c4762a1bSJed Brown   DM          idm;
428c4762a1bSJed Brown   PetscInt    p;
429c4762a1bSJed Brown   PetscMPIInt rank;
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
433dd400576SPatrick Sanan   if (rank == 0) {
434c4762a1bSJed Brown     switch (testNum) {
4359371c9d4SSatish Balay     case 0: {
436c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
437c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
438c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4, 5, 4, 3};
439c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
440c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
441c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
442c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
443c4762a1bSJed Brown 
4449566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4459566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
4469566063dSJacob Faibussowitsch       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4479566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
4489566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
4499371c9d4SSatish Balay     } break;
4509371c9d4SSatish Balay     case 1: {
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));
4629371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
4639371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 3, 1));
4649371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
4659371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 2, 2));
4669371c9d4SSatish Balay     } break;
4675b9dfbb6SMatthew G. Knepley     case 2:
4685b9dfbb6SMatthew G. Knepley     case 3:
4699371c9d4SSatish Balay     case 4: {
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};
4729371c9d4SSatish Balay       PetscInt    cones[18]            = {6, 7, 9, 9, 12, 11, 7, 12, 9, 7, 8, 10, 10, 13, 12, 7, 10, 12};
4735b9dfbb6SMatthew G. Knepley       PetscInt    coneOrientations[18] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
4749371c9d4SSatish Balay       PetscScalar vertexCoords[16]     = {
4759371c9d4SSatish Balay         -1., -1., 0., -1., 1., -1., -1., 0., 1., 0., -1., 1., 0., 1., 1., 1.,
4769371c9d4SSatish Balay       };
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]));
4889371c9d4SSatish Balay       if (testNum == 2)
4899371c9d4SSatish Balay         for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4909371c9d4SSatish Balay       if (testNum == 3 || testNum == 4)
4919371c9d4SSatish Balay         for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1));
4929371c9d4SSatish Balay     } break;
4939371c9d4SSatish Balay     case 5: {
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));
5059371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 0, 0));
5069371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 2, 0));
5079371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 4, 0));
5089371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
5099371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
5109371c9d4SSatish Balay       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
5119371c9d4SSatish Balay     } break;
512d71ae5a4SJacob Faibussowitsch     default:
513d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
514c4762a1bSJed Brown     }
515c4762a1bSJed Brown   } else {
516c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
517c4762a1bSJed Brown 
5189566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
5195b9dfbb6SMatthew G. Knepley     if (testNum == 3 || testNum == 4) PetscCall(DMCreateLabel(*dm, "pfault"));
5205b9dfbb6SMatthew G. Knepley     else PetscCall(DMCreateLabel(*dm, "fault"));
521c4762a1bSJed Brown   }
5229566063dSJacob Faibussowitsch   PetscCall(DMPlexInterpolate(*dm, &idm));
5239566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
5249566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
525c4762a1bSJed Brown   *dm = idm;
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
527c4762a1bSJed Brown }
528c4762a1bSJed Brown 
CreateSimplex_3D(MPI_Comm comm,AppCtx * user,DM dm)529d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, AppCtx *user, DM dm)
530d71ae5a4SJacob Faibussowitsch {
531c4762a1bSJed Brown   PetscInt    depth = 3, testNum = user->testNum, p;
532c4762a1bSJed Brown   PetscMPIInt rank;
533c4762a1bSJed Brown 
534c4762a1bSJed Brown   PetscFunctionBegin;
5359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
536dd400576SPatrick Sanan   if (rank == 0) {
537c4762a1bSJed Brown     switch (testNum) {
5389371c9d4SSatish Balay     case 0: {
539c4762a1bSJed Brown       PetscInt    numPoints[4]         = {5, 7, 9, 2};
540c4762a1bSJed 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};
541c4762a1bSJed 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};
542b5a892a1SMatthew 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};
543c4762a1bSJed 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};
544c4762a1bSJed Brown       PetscInt    markerPoints[20]     = {2, 1, 3, 1, 4, 1, 5, 1, 14, 1, 15, 1, 16, 1, 17, 1, 18, 1, 19, 1};
545c4762a1bSJed Brown       PetscInt    faultPoints[3]       = {3, 4, 5};
546c4762a1bSJed Brown 
5479566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
54848a46eb9SPierre Jolivet       for (p = 0; p < 10; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
54948a46eb9SPierre Jolivet       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1));
5509566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 0, 1));
5519566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 1, 2));
5529371c9d4SSatish Balay     } break;
5539371c9d4SSatish Balay     case 1: {
554c4762a1bSJed Brown       PetscInt    numPoints[4]         = {6, 13, 12, 4};
555c4762a1bSJed 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};
5569371c9d4SSatish Balay       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,
5579371c9d4SSatish Balay                                           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};
5589371c9d4SSatish Balay       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,
5599371c9d4SSatish Balay                                           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};
560c4762a1bSJed 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};
561c4762a1bSJed Brown       PetscInt    markerPoints[14]     = {5, 1, 6, 1, 7, 1, 10, 1, 22, 1, 23, 1, 24, 1};
562c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {5, 6, 7, 8};
563c4762a1bSJed Brown 
5649566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
56548a46eb9SPierre Jolivet       for (p = 0; p < 7; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
56648a46eb9SPierre Jolivet       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1));
5679566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 0, 1));
5689566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(dm, "material", 1, 2));
5699371c9d4SSatish Balay     } break;
570d71ae5a4SJacob Faibussowitsch     default:
571d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
572c4762a1bSJed Brown     }
573c4762a1bSJed Brown   } else {
574c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
575c4762a1bSJed Brown 
5769566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
5779566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "fault"));
578c4762a1bSJed Brown   }
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
580c4762a1bSJed Brown }
581c4762a1bSJed Brown 
CreateQuad_2D(MPI_Comm comm,PetscInt testNum,DM * dm)582d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
583d71ae5a4SJacob Faibussowitsch {
584c4762a1bSJed Brown   DM          idm;
585c4762a1bSJed Brown   PetscInt    p;
586c4762a1bSJed Brown   PetscMPIInt rank;
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   PetscFunctionBegin;
5899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
590dd400576SPatrick Sanan   if (rank == 0) {
591c4762a1bSJed Brown     switch (testNum) {
592c4762a1bSJed Brown     case 0:
593c4762a1bSJed Brown     case 2:
5949371c9d4SSatish Balay     case 3: {
595c4762a1bSJed Brown       PetscInt    numPoints[2]        = {6, 2};
596c4762a1bSJed Brown       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
597c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 3, 4, 5, 3, 6, 7, 4};
598c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
599c4762a1bSJed 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};
600c4762a1bSJed Brown       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
601c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
602c4762a1bSJed Brown 
6039566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6049566063dSJacob Faibussowitsch       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
6059371c9d4SSatish Balay       if (testNum == 0)
6069371c9d4SSatish Balay         for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
6079371c9d4SSatish Balay       if (testNum == 2 || testNum == 3)
6089371c9d4SSatish Balay         for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1));
6099566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6109566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
6119371c9d4SSatish Balay     } break;
6129371c9d4SSatish Balay     case 1: {
613c4762a1bSJed Brown       PetscInt    numPoints[2]         = {16, 9};
614c4762a1bSJed 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};
6159371c9d4SSatish Balay       PetscInt    cones[36]            = {9, 13, 14, 10, 10, 14, 15, 11, 11, 15, 16, 12, 13, 17, 18, 14, 14, 18, 19, 15, 15, 19, 20, 16, 17, 21, 22, 18, 18, 22, 23, 19, 19, 23, 24, 20};
616c4762a1bSJed 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};
6179371c9d4SSatish Balay       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, -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};
6185b9dfbb6SMatthew G. Knepley       PetscInt    faultPoints[4]       = {11, 15, 19, 20};
619c4762a1bSJed Brown       PetscInt    fault2Points[2]      = {17, 18};
620c4762a1bSJed Brown 
6219566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6225b9dfbb6SMatthew G. Knepley       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
6235b9dfbb6SMatthew G. Knepley       for (p = 3; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "faultBd", faultPoints[p], 1));
6249566063dSJacob Faibussowitsch       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault2", fault2Points[p], 1));
6259566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6269566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
6279566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
6289566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 1));
6299566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 4, 1));
6309566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
6319566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 6, 2));
6329566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 7, 2));
6339566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 8, 2));
6349371c9d4SSatish Balay     } break;
635d71ae5a4SJacob Faibussowitsch     default:
636d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
637c4762a1bSJed Brown     }
638c4762a1bSJed Brown   } else {
639c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
640c4762a1bSJed Brown 
6419566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
6425b9dfbb6SMatthew G. Knepley     if (testNum == 2 || testNum == 3) PetscCall(DMCreateLabel(*dm, "pfault"));
6439566063dSJacob Faibussowitsch     else PetscCall(DMCreateLabel(*dm, "fault"));
644c4762a1bSJed Brown   }
6459566063dSJacob Faibussowitsch   PetscCall(DMPlexInterpolate(*dm, &idm));
6469566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
6479566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
648c4762a1bSJed Brown   *dm = idm;
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
650c4762a1bSJed Brown }
651c4762a1bSJed Brown 
CreateHex_3D(MPI_Comm comm,PetscInt testNum,DM * dm)652d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
653d71ae5a4SJacob Faibussowitsch {
654c4762a1bSJed Brown   DM          idm;
655c4762a1bSJed Brown   PetscInt    depth = 3, p;
656c4762a1bSJed Brown   PetscMPIInt rank;
657c4762a1bSJed Brown 
658c4762a1bSJed Brown   PetscFunctionBegin;
6599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
660dd400576SPatrick Sanan   if (rank == 0) {
661c4762a1bSJed Brown     switch (testNum) {
6629371c9d4SSatish Balay     case 0: {
663c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
664c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
665c4762a1bSJed Brown       PetscInt    cones[16]            = {2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8};
666c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
6679371c9d4SSatish Balay       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, -0.5, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, -0.5, 1.0, 1.0, 0.5, 0.0, 0.0, 0.5, 1.0, 0.0, 0.5, 0.0, 1.0, 0.5, 1.0, 1.0};
668c4762a1bSJed Brown       PetscInt    markerPoints[52]     = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1};
669c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {3, 4, 7, 8};
670c4762a1bSJed Brown 
6719566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6729566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
6739566063dSJacob Faibussowitsch       for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(idm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
6749566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
6759566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
6769566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 2));
6779371c9d4SSatish Balay     } break;
6789371c9d4SSatish Balay     case 1: {
679c4762a1bSJed Brown       /* Cell Adjacency Graph:
680c4762a1bSJed Brown         0 -- { 8, 13, 21, 24} --> 1
681c4762a1bSJed Brown         0 -- {20, 21, 23, 24} --> 5 F
682c4762a1bSJed Brown         1 -- {10, 15, 21, 24} --> 2
683c4762a1bSJed Brown         1 -- {13, 14, 15, 24} --> 6
684c4762a1bSJed Brown         2 -- {21, 22, 24, 25} --> 4 F
685c4762a1bSJed Brown         3 -- {21, 24, 30, 35} --> 4
686c4762a1bSJed Brown         3 -- {21, 24, 28, 33} --> 5
687c4762a1bSJed Brown        */
688c4762a1bSJed Brown       PetscInt numPoints[2] = {30, 7};
689c4762a1bSJed 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};
6909371c9d4SSatish Balay       PetscInt cones[56] = {8, 21, 20, 7, 13, 12, 23, 24, 14, 15, 10, 9, 13, 8, 21, 24, 15, 16, 11, 10, 24, 21, 22, 25, 30, 29, 28, 21, 35, 24, 33, 34, 24, 21, 30, 35, 25, 36, 31, 22, 27, 20, 21, 28, 32, 33, 24, 23, 15, 24, 13, 14, 19, 18, 17, 26};
691c4762a1bSJed 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};
6929371c9d4SSatish Balay       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,  -2.0, -1.0, 0.0, -3.0, 0.0, 0.0, -2.0, 1.0, 0.0, -2.0, 2.0, 0.0,
6939371c9d4SSatish Balay                                        -2.0, -1.0, 2.0,  -3.0, 0.0,  2.0,  -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, 0.0,  2.0, 0.0, 0.0,  0.0, 2.0,
6949371c9d4SSatish Balay                                        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,  2.0,  -1.0, 0.0, 3.0,  0.0, 0.0, 2.0,  1.0, 0.0, 2.0,  2.0, 0.0};
695c4762a1bSJed Brown       PetscInt    faultPoints[6]    = {20, 21, 22, 23, 24, 25};
696c4762a1bSJed Brown 
6979566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
6989566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
6999566063dSJacob Faibussowitsch       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
7009566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
7019566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
7029566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 1));
7039566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
7049566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 4, 2));
7059566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 5, 2));
7069566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 6, 2));
7079371c9d4SSatish Balay     } break;
7089371c9d4SSatish Balay     case 2: {
709c4762a1bSJed Brown       /* Buried fault edge */
710c4762a1bSJed Brown       PetscInt    numPoints[2]         = {18, 4};
711c4762a1bSJed 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};
7129371c9d4SSatish Balay       PetscInt    cones[32]            = {4, 5, 8, 7, 13, 16, 17, 14, 5, 6, 9, 8, 14, 17, 18, 15, 7, 8, 11, 10, 16, 19, 20, 17, 8, 9, 12, 11, 17, 20, 21, 18};
713c4762a1bSJed 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};
7149371c9d4SSatish Balay       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, 2.0, -2.0, 0.0, 2.0, 0.0, 0.0, 2.0, 2.0, 0.0,
7159371c9d4SSatish Balay                                           -2.0, -2.0, 2.0, -2.0, 0.0, 2.0, -2.0, 2.0, 2.0, 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};
716c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {7, 8, 16, 17};
717c4762a1bSJed Brown 
7189566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
7199566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
7209566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1));
7219566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 0, 1));
7229566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 1, 1));
7239566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 2, 2));
7249566063dSJacob Faibussowitsch       PetscCall(DMSetLabelValue(*dm, "material", 3, 2));
7259371c9d4SSatish Balay     } break;
726d71ae5a4SJacob Faibussowitsch     default:
727d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
728c4762a1bSJed Brown     }
729c4762a1bSJed Brown   } else {
730c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
731c4762a1bSJed Brown 
7329566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
7339566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
7349566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(idm, "fault"));
735c4762a1bSJed Brown   }
7369566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
7379566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
738c4762a1bSJed Brown   *dm = idm;
7393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
740c4762a1bSJed Brown }
741c4762a1bSJed Brown 
CreateFaultLabel(DM dm)742d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFaultLabel(DM dm)
743d71ae5a4SJacob Faibussowitsch {
744b253942bSMatthew G. Knepley   DMLabel  label;
745b253942bSMatthew G. Knepley   PetscInt dim, h, pStart, pEnd, pMax, p;
746b253942bSMatthew G. Knepley 
747b253942bSMatthew G. Knepley   PetscFunctionBegin;
7489566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7499566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "cohesive"));
7509566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &label));
751b253942bSMatthew G. Knepley   for (h = 0; h <= dim; ++h) {
7529566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(dm, h, NULL, &pMax));
7539566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &pStart, &pEnd));
7549566063dSJacob Faibussowitsch     for (p = pMax; p < pEnd; ++p) PetscCall(DMLabelSetValue(label, p, 1));
755b253942bSMatthew G. Knepley   }
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757b253942bSMatthew G. Knepley }
758b253942bSMatthew G. Knepley 
759989fa639SBrad Aagaard /*
760989fa639SBrad Aagaard   Create a displacement field, and some number of vector fault fields
761989fa639SBrad Aagaard */
CreateDiscretization(DM dm,AppCtx * user)762d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateDiscretization(DM dm, AppCtx *user)
763d71ae5a4SJacob Faibussowitsch {
764b253942bSMatthew G. Knepley   PetscFE  fe;
765b253942bSMatthew G. Knepley   DMLabel  fault;
766b7519becSMatthew G. Knepley   PetscInt dim, Ncf = user->cohesiveFields, f;
767b253942bSMatthew G. Knepley 
768b253942bSMatthew G. Knepley   PetscFunctionBegin;
7699566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7709566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &fault));
7719566063dSJacob Faibussowitsch   PetscCall(DMLabelView(fault, PETSC_VIEWER_STDOUT_WORLD));
772b253942bSMatthew G. Knepley 
7739566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, user->cellSimplex, "displacement_", PETSC_DETERMINE, &fe));
7749566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(fe, "displacement"));
7759566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
7769566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
777b253942bSMatthew G. Knepley 
778b7519becSMatthew G. Knepley   if (Ncf > 0) {
7799566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim - 1, dim, user->cellSimplex, "faulttraction_", PETSC_DETERMINE, &fe));
7809566063dSJacob Faibussowitsch     PetscCall(PetscFESetName(fe, "fault traction"));
7819566063dSJacob Faibussowitsch     PetscCall(DMAddField(dm, fault, (PetscObject)fe));
7829566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
783b7519becSMatthew G. Knepley   }
784b7519becSMatthew G. Knepley   for (f = 1; f < Ncf; ++f) {
785b7519becSMatthew G. Knepley     char name[256], opt[256];
786b7519becSMatthew G. Knepley 
78763a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, 256, "fault field %" PetscInt_FMT, f));
78863a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(opt, 256, "faultfield_%" PetscInt_FMT "_", f));
7899566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim - 1, dim, user->cellSimplex, opt, PETSC_DETERMINE, &fe));
7909566063dSJacob Faibussowitsch     PetscCall(PetscFESetName(fe, name));
7919566063dSJacob Faibussowitsch     PetscCall(DMAddField(dm, fault, (PetscObject)fe));
7929566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
793b7519becSMatthew G. Knepley   }
794b253942bSMatthew G. Knepley 
7959566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
7963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
797b253942bSMatthew G. Knepley }
798b253942bSMatthew G. Knepley 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)799d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
800d71ae5a4SJacob Faibussowitsch {
801c4762a1bSJed Brown   PetscInt    dim         = user->dim;
802c4762a1bSJed Brown   PetscBool   cellSimplex = user->cellSimplex, hasFault, hasFault2, hasParallelFault;
803c4762a1bSJed Brown   PetscMPIInt rank, size;
804dd96b1f9SToby Isaac   DMLabel     matLabel;
805c4762a1bSJed Brown 
806c4762a1bSJed Brown   PetscFunctionBegin;
8079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
8099566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
8109566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
8119566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
812c4762a1bSJed Brown   switch (dim) {
813c4762a1bSJed Brown   case 2:
814c4762a1bSJed Brown     if (cellSimplex) {
8159566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_2D(comm, user->testNum, dm));
816c4762a1bSJed Brown     } else {
8179566063dSJacob Faibussowitsch       PetscCall(CreateQuad_2D(comm, user->testNum, dm));
818c4762a1bSJed Brown     }
819c4762a1bSJed Brown     break;
820c4762a1bSJed Brown   case 3:
821c4762a1bSJed Brown     if (cellSimplex) {
8229566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_3D(comm, user, *dm));
823c4762a1bSJed Brown     } else {
8249566063dSJacob Faibussowitsch       PetscCall(CreateHex_3D(comm, user->testNum, dm));
825c4762a1bSJed Brown     }
826c4762a1bSJed Brown     break;
827d71ae5a4SJacob Faibussowitsch   default:
828d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %" PetscInt_FMT, dim);
829c4762a1bSJed Brown   }
8309566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "orig_"));
8319566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
8329566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
8339566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
8349566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "fault", &hasFault));
835c4762a1bSJed Brown   if (hasFault) {
836b253942bSMatthew G. Knepley     DM      dmHybrid = NULL, dmInterface = NULL;
837b253942bSMatthew G. Knepley     DMLabel faultLabel, faultBdLabel, hybridLabel, splitLabel;
838c4762a1bSJed Brown 
8399566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
8409566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "faultBd", &faultBdLabel));
841caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, &splitLabel, &dmInterface, &dmHybrid));
8429566063dSJacob Faibussowitsch     PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD));
8439566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
8449566063dSJacob Faibussowitsch     PetscCall(DMLabelView(splitLabel, PETSC_VIEWER_STDOUT_WORLD));
8459566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&splitLabel));
8469566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_interface_view"));
8479566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmInterface));
8489566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
849c4762a1bSJed Brown     *dm = dmHybrid;
850c4762a1bSJed Brown   }
8519566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "fault2", &hasFault2));
852c4762a1bSJed Brown   if (hasFault2) {
853c4762a1bSJed Brown     DM      dmHybrid = NULL;
854c4762a1bSJed Brown     DMLabel faultLabel, faultBdLabel, hybridLabel;
855c4762a1bSJed Brown 
8569566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "faulted_"));
8579566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
8589566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
8599566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
8609566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
8619566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault2", &faultLabel));
8629566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault2Bd", &faultBdLabel));
863caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, NULL, &dmHybrid));
8649566063dSJacob Faibussowitsch     PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD));
8659566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
8669566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
867c4762a1bSJed Brown     *dm = dmHybrid;
868c4762a1bSJed Brown   }
869c4762a1bSJed Brown   if (user->testPartition && size > 1) {
870c4762a1bSJed Brown     PetscPartitioner part;
871c4762a1bSJed Brown     PetscInt        *sizes  = NULL;
872c4762a1bSJed Brown     PetscInt        *points = NULL;
873c4762a1bSJed Brown 
874dd400576SPatrick Sanan     if (rank == 0) {
875c4762a1bSJed Brown       if (dim == 2 && cellSimplex && size == 2) {
876c4762a1bSJed Brown         switch (user->testNum) {
877c4762a1bSJed Brown         case 0: {
878c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {1, 2};
879c4762a1bSJed Brown           PetscInt triPoints_p2[3] = {0, 1, 2};
880c4762a1bSJed Brown 
8819566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
8829566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, triSizes_p2, 2));
8839371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, triPoints_p2, 3));
8849371c9d4SSatish Balay           break;
8859371c9d4SSatish Balay         }
8865b9dfbb6SMatthew G. Knepley         case 3: {
8875b9dfbb6SMatthew G. Knepley           PetscInt triSizes_p2[2]  = {3, 3};
8885b9dfbb6SMatthew G. Knepley           PetscInt triPoints_p2[6] = {0, 1, 2, 3, 4, 5};
8895b9dfbb6SMatthew G. Knepley 
8905b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(2, &sizes, 6, &points));
8915b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes, triSizes_p2, 2));
8929371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, triPoints_p2, 6));
8939371c9d4SSatish Balay           break;
8949371c9d4SSatish Balay         }
895d71ae5a4SJacob Faibussowitsch         default:
896d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 2 procs", user->testNum);
8975b9dfbb6SMatthew G. Knepley         }
8985b9dfbb6SMatthew G. Knepley       } else if (dim == 2 && cellSimplex && size == 6) {
8995b9dfbb6SMatthew G. Knepley         switch (user->testNum) {
9005b9dfbb6SMatthew G. Knepley         case 4: {
9015b9dfbb6SMatthew G. Knepley           PetscInt triSizes_p6[6]  = {1, 1, 1, 1, 1, 1};
9025b9dfbb6SMatthew G. Knepley           PetscInt triPoints_p6[6] = {0, 1, 2, 3, 4, 5};
9035b9dfbb6SMatthew G. Knepley 
9045b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(6, &sizes, 6, &points));
9055b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes, triSizes_p6, 6));
9069371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, triPoints_p6, 6));
9079371c9d4SSatish Balay           break;
9089371c9d4SSatish Balay         }
909d71ae5a4SJacob Faibussowitsch         default:
910d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 6 procs", user->testNum);
911c4762a1bSJed Brown         }
912c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && size == 2) {
913c4762a1bSJed Brown         switch (user->testNum) {
914c4762a1bSJed Brown         case 0: {
915c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 2};
916c4762a1bSJed Brown           PetscInt quadPoints_p2[3] = {0, 1, 2};
917c4762a1bSJed Brown 
9189566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9199566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2));
9209371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, quadPoints_p2, 3));
9219371c9d4SSatish Balay           break;
9229371c9d4SSatish Balay         }
923c4762a1bSJed Brown         case 2: {
924c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 1};
925c4762a1bSJed Brown           PetscInt quadPoints_p2[2] = {0, 1};
926c4762a1bSJed Brown 
9279566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
9289566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2));
9299371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));
9309371c9d4SSatish Balay           break;
9319371c9d4SSatish Balay         }
9325b9dfbb6SMatthew G. Knepley         case 3: {
9335b9dfbb6SMatthew G. Knepley           PetscInt quadSizes_p2[2]  = {1, 1};
9345b9dfbb6SMatthew G. Knepley           PetscInt quadPoints_p2[2] = {1, 0};
9355b9dfbb6SMatthew G. Knepley 
9365b9dfbb6SMatthew G. Knepley           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
9375b9dfbb6SMatthew G. Knepley           PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2));
9389371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));
9399371c9d4SSatish Balay           break;
9409371c9d4SSatish Balay         }
941d71ae5a4SJacob Faibussowitsch         default:
942d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral mesh on 2 procs", user->testNum);
943c4762a1bSJed Brown         }
944c4762a1bSJed Brown       } else if (dim == 3 && cellSimplex && size == 2) {
945c4762a1bSJed Brown         switch (user->testNum) {
946c4762a1bSJed Brown         case 0: {
947c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 2};
948c4762a1bSJed Brown           PetscInt tetPoints_p2[3] = {0, 1, 2};
949c4762a1bSJed Brown 
9509566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9519566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
9529371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, tetPoints_p2, 3));
9539371c9d4SSatish Balay           break;
9549371c9d4SSatish Balay         }
955d71ae5a4SJacob Faibussowitsch         default:
956d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for teterehedral mesh on 2 procs", user->testNum);
957c4762a1bSJed Brown         }
958c4762a1bSJed Brown       } else if (dim == 3 && !cellSimplex && size == 2) {
959c4762a1bSJed Brown         switch (user->testNum) {
960c4762a1bSJed Brown         case 0: {
961c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {1, 2};
962c4762a1bSJed Brown           PetscInt hexPoints_p2[3] = {0, 1, 2};
963c4762a1bSJed Brown 
9649566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
9659566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2));
9669371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, hexPoints_p2, 3));
9679371c9d4SSatish Balay           break;
9689371c9d4SSatish Balay         }
969d71ae5a4SJacob Faibussowitsch         default:
970d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral mesh on 2 procs", user->testNum);
971c4762a1bSJed Brown         }
972c4762a1bSJed Brown       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
973c4762a1bSJed Brown     }
9749566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
9759566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
9769566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
9779566063dSJacob Faibussowitsch     PetscCall(PetscFree2(sizes, points));
978c4762a1bSJed Brown   }
97967de2c8aSMatthew G. Knepley   PetscCall(DMGetLabel(*dm, "material", &matLabel));
98067de2c8aSMatthew G. Knepley   if (matLabel) PetscCall(DMPlexLabelComplete(*dm, matLabel));
981c4762a1bSJed Brown   {
982c4762a1bSJed Brown     DM pdm = NULL;
983c4762a1bSJed Brown 
984c4762a1bSJed Brown     /* Distribute mesh over processes */
9859566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
986c4762a1bSJed Brown     if (pdm) {
9879566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
9889566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
989c4762a1bSJed Brown       *dm = pdm;
990c4762a1bSJed Brown     }
991c4762a1bSJed Brown   }
9929566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "pfault", &hasParallelFault));
993c4762a1bSJed Brown   if (hasParallelFault) {
9945b9dfbb6SMatthew G. Knepley     DM      dmHybrid = NULL, dmInterface;
995c4762a1bSJed Brown     DMLabel faultLabel, faultBdLabel, hybridLabel;
996c4762a1bSJed Brown 
9979566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "pfault", &faultLabel));
9989566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "pfaultBd", &faultBdLabel));
999caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, &dmInterface, &dmHybrid));
10005b9dfbb6SMatthew G. Knepley     PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_fault_view"));
10015b9dfbb6SMatthew G. Knepley     {
10025b9dfbb6SMatthew G. Knepley       PetscViewer viewer;
10035b9dfbb6SMatthew G. Knepley       PetscMPIInt rank;
10045b9dfbb6SMatthew G. Knepley 
10055b9dfbb6SMatthew G. Knepley       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)*dm), &rank));
10065b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
10075b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "Rank %d\n", rank));
10085b9dfbb6SMatthew G. Knepley       PetscCall(DMLabelView(hybridLabel, viewer));
10095b9dfbb6SMatthew G. Knepley       PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
10105b9dfbb6SMatthew G. Knepley     }
10119566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
10125b9dfbb6SMatthew G. Knepley     PetscCall(DMDestroy(&dmInterface));
10139566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
1014c4762a1bSJed Brown     *dm = dmHybrid;
1015c4762a1bSJed Brown   }
10169566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, "Hybrid Mesh"));
10179566063dSJacob Faibussowitsch   PetscCall(CreateFaultLabel(*dm));
10189566063dSJacob Faibussowitsch   PetscCall(CreateDiscretization(*dm, user));
10199566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
10209566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
10219566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
10229566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
10233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1024c4762a1bSJed Brown }
1025c4762a1bSJed Brown 
TestMesh(DM dm,AppCtx * user)1026d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestMesh(DM dm, AppCtx *user)
1027d71ae5a4SJacob Faibussowitsch {
1028b253942bSMatthew G. Knepley   PetscFunctionBegin;
10299566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckSymmetry(dm));
10309566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckSkeleton(dm, 0));
10319566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckFaces(dm, 0));
10323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1033b253942bSMatthew G. Knepley }
1034b253942bSMatthew G. Knepley 
TestDiscretization(DM dm,AppCtx * user)1035d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestDiscretization(DM dm, AppCtx *user)
1036d71ae5a4SJacob Faibussowitsch {
1037b253942bSMatthew G. Knepley   PetscSection s;
1038b253942bSMatthew G. Knepley 
1039b253942bSMatthew G. Knepley   PetscFunctionBegin;
1040882cb04dSMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &s));
10419566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)s, NULL, "-local_section_view"));
10423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1043b253942bSMatthew G. Knepley }
1044b253942bSMatthew G. Knepley 
r(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)1045*2a8381b2SBarry Smith static PetscErrorCode r(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
1046d71ae5a4SJacob Faibussowitsch {
1047b253942bSMatthew G. Knepley   PetscInt d;
1048b253942bSMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = x[d];
10493ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1050b253942bSMatthew G. Knepley }
1051b253942bSMatthew G. Knepley 
rp1(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)1052*2a8381b2SBarry Smith static PetscErrorCode rp1(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
1053d71ae5a4SJacob Faibussowitsch {
1054b253942bSMatthew G. Knepley   PetscInt d;
1055b253942bSMatthew G. Knepley   for (d = 0; d < dim; ++d) u[d] = x[d] + (d > 0 ? 1.0 : 0.0);
10563ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1057b253942bSMatthew G. Knepley }
1058b253942bSMatthew G. Knepley 
phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)1059*2a8381b2SBarry Smith static PetscErrorCode phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
1060d71ae5a4SJacob Faibussowitsch {
1061b253942bSMatthew G. Knepley   PetscInt d;
1062b253942bSMatthew G. Knepley   u[0] = -x[1];
1063b253942bSMatthew G. Knepley   u[1] = x[0];
1064b253942bSMatthew G. Knepley   for (d = 2; d < dim; ++d) u[d] = x[d];
10653ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1066b253942bSMatthew G. Knepley }
1067b253942bSMatthew G. Knepley 
add_fields(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f[])106867de2c8aSMatthew G. Knepley static void add_fields(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
106967de2c8aSMatthew G. Knepley {
107067de2c8aSMatthew G. Knepley   PetscInt       d;
107167de2c8aSMatthew G. Knepley   const PetscInt offN = 0;
107267de2c8aSMatthew G. Knepley   const PetscInt offP = dim;
107367de2c8aSMatthew G. Knepley   for (d = 0; d < dim; ++d) f[d] = u[offN + d] + u[offP + d];
107467de2c8aSMatthew G. Knepley }
107567de2c8aSMatthew G. Knepley 
normal_field(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f[])107667de2c8aSMatthew G. Knepley static void normal_field(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
107767de2c8aSMatthew G. Knepley {
107867de2c8aSMatthew G. Knepley   PetscInt d;
107967de2c8aSMatthew G. Knepley   for (d = 0; d < dim; ++d) f[d] = n[d];
108067de2c8aSMatthew G. Knepley }
108167de2c8aSMatthew G. Knepley 
1082b253942bSMatthew G. Knepley /* \lambda \cdot (\psi_u^- - \psi_u^+) */
f0_bd_u_neg(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])108381363e6fSMatthew G. Knepley static void f0_bd_u_neg(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1084d71ae5a4SJacob Faibussowitsch {
108581363e6fSMatthew G. Knepley   const PetscInt Nc = dim + 1;
108681363e6fSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] = -u[uOff[1] + c];
1087b253942bSMatthew G. Knepley }
108881363e6fSMatthew G. Knepley 
f0_bd_u_pos(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])108981363e6fSMatthew G. Knepley static void f0_bd_u_pos(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
109081363e6fSMatthew G. Knepley {
109181363e6fSMatthew G. Knepley   const PetscInt Nc = dim + 1;
109281363e6fSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] = u[uOff[1] + c];
1093b253942bSMatthew G. Knepley }
1094b253942bSMatthew G. Knepley 
1095b253942bSMatthew G. Knepley /* (d - u^+ + u^-) \cdot \psi_\lambda */
f0_bd_l(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])1096d71ae5a4SJacob Faibussowitsch static void f0_bd_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
1097d71ae5a4SJacob Faibussowitsch {
1098b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[2] - uOff[1];
1099b253942bSMatthew G. Knepley 
110081363e6fSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] = (c > 0 ? 1.0 : 0.0) + u[c] - u[Nc + c];
1101b253942bSMatthew G. Knepley }
1102b253942bSMatthew G. Knepley 
1103b253942bSMatthew G. Knepley /* \psi_lambda \cdot (\psi_u^- - \psi_u^+) */
g0_bd_ul_neg(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g0[])110481363e6fSMatthew G. Knepley static void g0_bd_ul_neg(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1105d71ae5a4SJacob Faibussowitsch {
110681363e6fSMatthew G. Knepley   const PetscInt Nc = dim + 1;
110781363e6fSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) g0[c * Nc + c] = -1.0;
1108b253942bSMatthew G. Knepley }
110981363e6fSMatthew G. Knepley 
g0_bd_ul_pos(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g0[])111081363e6fSMatthew G. Knepley static void g0_bd_ul_pos(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
111181363e6fSMatthew G. Knepley {
111281363e6fSMatthew G. Knepley   const PetscInt Nc = dim + 1;
111381363e6fSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0;
1114b253942bSMatthew G. Knepley }
1115b253942bSMatthew G. Knepley 
1116b253942bSMatthew G. Knepley /* (-\psi_u^+ + \psi_u^-) \cdot \psi_\lambda */
g0_bd_lu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g0[])1117d71ae5a4SJacob Faibussowitsch static void g0_bd_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1118d71ae5a4SJacob Faibussowitsch {
1119b253942bSMatthew G. Knepley   const PetscInt Nc = uOff[2] - uOff[1];
1120b253942bSMatthew G. Knepley 
112181363e6fSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) {
112281363e6fSMatthew G. Knepley     g0[c * Nc + c]           = -1.0;
112381363e6fSMatthew G. Knepley     g0[Nc * Nc + c * Nc + c] = 1.0;
1124b253942bSMatthew G. Knepley   }
1125b253942bSMatthew G. Knepley }
1126b253942bSMatthew G. Knepley 
TestAssembly(DM dm,AppCtx * user)1127d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestAssembly(DM dm, AppCtx *user)
1128d71ae5a4SJacob Faibussowitsch {
1129b253942bSMatthew G. Knepley   Mat           J;
113067de2c8aSMatthew G. Knepley   Vec           locX, locF, locW;
1131b253942bSMatthew G. Knepley   PetscDS       probh;
1132b253942bSMatthew G. Knepley   DMLabel       fault, material;
113367de2c8aSMatthew G. Knepley   DM            dmFault;
1134b253942bSMatthew G. Knepley   IS            cohesiveCells;
113567de2c8aSMatthew G. Knepley   PetscFE       fe;
1136b7519becSMatthew G. Knepley   PetscWeakForm wf;
113706ad1575SMatthew G. Knepley   PetscFormKey  keys[3];
1138*2a8381b2SBarry Smith   PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], PetscCtx ctx);
1139b253942bSMatthew G. Knepley   PetscInt    dim, Nf, cMax, cEnd, id;
1140b7519becSMatthew G. Knepley   PetscMPIInt rank;
1141b253942bSMatthew G. Knepley 
1142b253942bSMatthew G. Knepley   PetscFunctionBegin;
114367de2c8aSMatthew G. Knepley   if (!user->testAssembly) PetscFunctionReturn(PETSC_SUCCESS);
11449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
11459566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
11469566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, NULL, &cMax));
11479566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
11489566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, cEnd - cMax, cMax, 1, &cohesiveCells));
11499566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "cohesive", &fault));
11509566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
11519566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)locX, "Local Solution"));
11529566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
11539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)locF, "Local Residual"));
11549566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
11559566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
1156b253942bSMatthew G. Knepley 
1157b253942bSMatthew G. Knepley   /* The initial guess has displacement shifted by one unit in each fault parallel direction across the fault */
11589566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "material", &material));
1159b253942bSMatthew G. Knepley   id              = 1;
1160b253942bSMatthew G. Knepley   initialGuess[0] = r;
1161b253942bSMatthew G. Knepley   initialGuess[1] = NULL;
11629566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
1163b253942bSMatthew G. Knepley   id              = 2;
1164b253942bSMatthew G. Knepley   initialGuess[0] = rp1;
1165b253942bSMatthew G. Knepley   initialGuess[1] = NULL;
11669566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
1167b253942bSMatthew G. Knepley   id              = 1;
1168b253942bSMatthew G. Knepley   initialGuess[0] = NULL;
1169b253942bSMatthew G. Knepley   initialGuess[1] = phi;
11709566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, fault, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX));
11719566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(locX, NULL, "-local_solution_view"));
1172b253942bSMatthew G. Knepley 
117367de2c8aSMatthew G. Knepley   /* Test projection to fault mesh */
117467de2c8aSMatthew G. Knepley   PetscCall(DMPlexCreateCohesiveSubmesh(dm, PETSC_FALSE, NULL, 0, &dmFault));
117567de2c8aSMatthew G. Knepley   PetscCall(DMViewFromOptions(dmFault, NULL, "-fault_view"));
117667de2c8aSMatthew G. Knepley   PetscCall(DMPlexOrient(dmFault));
117767de2c8aSMatthew G. Knepley   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim - 1, dim, user->cellSimplex, "fault_field_", PETSC_DETERMINE, &fe));
117867de2c8aSMatthew G. Knepley   PetscCall(PetscFESetName(fe, "fault_field"));
117967de2c8aSMatthew G. Knepley   PetscCall(DMAddField(dmFault, NULL, (PetscObject)fe));
118067de2c8aSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe));
118167de2c8aSMatthew G. Knepley   PetscCall(DMCreateDS(dmFault));
118267de2c8aSMatthew G. Knepley   PetscCall(DMGetLocalVector(dmFault, &locW));
118367de2c8aSMatthew G. Knepley   PetscCall(DMViewFromOptions(dmFault, NULL, "-cohesive_view"));
118467de2c8aSMatthew G. Knepley   void (*faultFuncs[1])(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]);
118567de2c8aSMatthew G. Knepley 
118667de2c8aSMatthew G. Knepley   DMLabel  depthLabel;
118767de2c8aSMatthew G. Knepley   PetscInt depth;
118867de2c8aSMatthew G. Knepley   PetscCall(DMPlexGetDepthLabel(dmFault, &depthLabel));
118967de2c8aSMatthew G. Knepley   PetscCall(DMPlexGetDepth(dmFault, &depth));
119067de2c8aSMatthew G. Knepley   id = depth - 1;
119167de2c8aSMatthew G. Knepley   /* w = r + rp1 */
119267de2c8aSMatthew G. Knepley   faultFuncs[0] = add_fields;
119367de2c8aSMatthew G. Knepley   PetscCall(DMProjectBdFieldLabelLocal(dmFault, 0.0, depthLabel, 1, &id, PETSC_DETERMINE, NULL, locX, faultFuncs, INSERT_VALUES, locW));
119467de2c8aSMatthew G. Knepley   PetscCall(VecViewFromOptions(locW, NULL, "-local_projection_view"));
119567de2c8aSMatthew G. Knepley 
119667de2c8aSMatthew G. Knepley   /* w = fault_normal */
119767de2c8aSMatthew G. Knepley   faultFuncs[0] = normal_field;
119867de2c8aSMatthew G. Knepley   PetscCall(DMProjectBdFieldLabelLocal(dmFault, 0.0, depthLabel, 1, &id, PETSC_DETERMINE, NULL, locX, faultFuncs, INSERT_VALUES, locW));
119967de2c8aSMatthew G. Knepley   PetscCall(VecViewFromOptions(locW, NULL, "-local_projection_view"));
120067de2c8aSMatthew G. Knepley   PetscCall(DMRestoreLocalVector(dmFault, &locW));
120167de2c8aSMatthew G. Knepley   PetscCall(DMDestroy(&dmFault));
120267de2c8aSMatthew G. Knepley 
1203e5939c1dSMatthew G. Knepley   PetscCall(DMGetCellDS(dm, cMax, &probh, NULL));
12049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(probh, &wf));
12059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(probh, &Nf));
120681363e6fSMatthew G. Knepley   PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 1, 0, 0, 0, f0_bd_u_neg, 0, NULL));
120781363e6fSMatthew G. Knepley   PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 2, 0, 0, 0, f0_bd_u_pos, 0, NULL));
120881363e6fSMatthew G. Knepley   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 1, 0, 1, 0, 0, g0_bd_ul_neg, 0, NULL, 0, NULL, 0, NULL));
120981363e6fSMatthew G. Knepley   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 2, 0, 1, 0, 0, g0_bd_ul_pos, 0, NULL, 0, NULL, 0, NULL));
1210b7519becSMatthew G. Knepley   if (Nf > 1) {
12119566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdResidual(wf, fault, 1, 1, 0, 0, f0_bd_l, 0, NULL));
12129566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdJacobian(wf, fault, 1, 1, 0, 0, 0, g0_bd_lu, 0, NULL, 0, NULL, 0, NULL));
1213b7519becSMatthew G. Knepley   }
1214c5853193SPierre Jolivet   if (rank == 0) PetscCall(PetscDSView(probh, NULL));
1215b253942bSMatthew G. Knepley 
121681363e6fSMatthew G. Knepley   keys[0].label = material;
121781363e6fSMatthew G. Knepley   keys[0].value = 1;
12186528b96dSMatthew G. Knepley   keys[0].field = 0;
1219dd96b1f9SToby Isaac   keys[0].part  = 0;
1220b7519becSMatthew G. Knepley   keys[1].label = material;
1221b7519becSMatthew G. Knepley   keys[1].value = 2;
12226528b96dSMatthew G. Knepley   keys[1].field = 0;
1223dd96b1f9SToby Isaac   keys[1].part  = 0;
1224b7519becSMatthew G. Knepley   keys[2].label = fault;
1225b7519becSMatthew G. Knepley   keys[2].value = 1;
1226b7519becSMatthew G. Knepley   keys[2].field = 1;
1227dd96b1f9SToby Isaac   keys[2].part  = 0;
12289566063dSJacob Faibussowitsch   PetscCall(VecSet(locF, 0.));
1229754e4fbaSMatthew G. Knepley   PetscCall(DMPlexComputeResidualHybridByKey(dm, keys, cohesiveCells, 0.0, locX, NULL, 0.0, locF, user));
12309566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(locF, NULL, "-local_residual_view"));
12319566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(J));
1232754e4fbaSMatthew G. Knepley   PetscCall(DMPlexComputeJacobianHybridByKey(dm, keys, cohesiveCells, 0.0, 0.0, locX, NULL, J, J, user));
123349664dceSMatthew G. Knepley   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
123449664dceSMatthew G. Knepley   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
12359566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-local_jacobian_view"));
1236b253942bSMatthew G. Knepley 
12379566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
12389566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
12399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
12409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cohesiveCells));
1241e5939c1dSMatthew G. Knepley 
1242e5939c1dSMatthew G. Knepley   if (cMax < cEnd) {
1243e5939c1dSMatthew G. Knepley     PetscDS         ds;
1244e5939c1dSMatthew G. Knepley     PetscFE         fe;
1245e5939c1dSMatthew G. Knepley     PetscQuadrature quad;
1246e5939c1dSMatthew G. Knepley     IS             *perm;
1247e5939c1dSMatthew G. Knepley     const PetscInt *cone;
1248e5939c1dSMatthew G. Knepley     PetscInt        Na, a;
1249e5939c1dSMatthew G. Knepley 
1250e5939c1dSMatthew G. Knepley     PetscCall(DMPlexGetCone(dm, cMax, &cone));
1251e5939c1dSMatthew G. Knepley     PetscCall(DMGetCellDS(dm, cMax, &ds, NULL));
1252e5939c1dSMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
1253e5939c1dSMatthew G. Knepley     PetscCall(PetscFEGetQuadrature(fe, &quad));
1254e5939c1dSMatthew G. Knepley     PetscCall(PetscQuadratureComputePermutations(quad, &Na, &perm));
1255e5939c1dSMatthew G. Knepley     for (a = 0; a < Na; ++a) PetscCall(ISDestroy(&perm[a]));
1256e5939c1dSMatthew G. Knepley     PetscCall(PetscFree(perm));
1257e5939c1dSMatthew G. Knepley   }
12583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1259b253942bSMatthew G. Knepley }
1260b253942bSMatthew G. Knepley 
main(int argc,char ** argv)1261d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
1262d71ae5a4SJacob Faibussowitsch {
1263c4762a1bSJed Brown   DM     dm;
1264c4762a1bSJed Brown   AppCtx user; /* user-defined work context */
1265c4762a1bSJed Brown 
1266327415f7SBarry Smith   PetscFunctionBeginUser;
12679566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
12689566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
12699566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
12709566063dSJacob Faibussowitsch   PetscCall(TestMesh(dm, &user));
12719566063dSJacob Faibussowitsch   PetscCall(TestDiscretization(dm, &user));
12729566063dSJacob Faibussowitsch   PetscCall(TestAssembly(dm, &user));
12739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
12749566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1275b122ec5aSJacob Faibussowitsch   return 0;
1276c4762a1bSJed Brown }
1277c4762a1bSJed Brown 
1278c4762a1bSJed Brown /*TEST
1279c4762a1bSJed Brown   testset:
1280b253942bSMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1281b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1 -local_section_view \
1282b253942bSMatthew G. Knepley           -local_solution_view -local_residual_view -local_jacobian_view
1283c4762a1bSJed Brown     test:
1284c4762a1bSJed Brown       suffix: tri_0
1285c4762a1bSJed Brown       args: -dim 2
128681363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1287c4762a1bSJed Brown     test:
1288882cb04dSMatthew G. Knepley       suffix: tri_0_perm
1289adc21957SMatthew G. Knepley       args: -dim 2 -dm_reorder_section -dm_reorder_section_type cohesive
1290882cb04dSMatthew 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" -e "s/_neg//g" -e "s/_pos//g" -e "s~_ZL.*~~g"
1291882cb04dSMatthew G. Knepley     test:
1292c4762a1bSJed Brown       suffix: tri_t1_0
1293c4762a1bSJed Brown       args: -dim 2 -test_num 1
129481363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1295c4762a1bSJed Brown     test:
1296882cb04dSMatthew G. Knepley       suffix: tri_t1_0_perm
1297adc21957SMatthew G. Knepley       args: -dim 2 -test_num 1 -dm_reorder_section -dm_reorder_section_type cohesive
1298882cb04dSMatthew 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" -e "s/_neg//g" -e "s/_pos//g" -e "s~_ZL.*~~g"
1299882cb04dSMatthew G. Knepley     test:
13005b9dfbb6SMatthew G. Knepley       suffix: tri_t2_0
13015b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 2
130281363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
130381363e6fSMatthew G. Knepley 
13045b9dfbb6SMatthew G. Knepley     test:
13053bdf08b3SMatthew G. Knepley       suffix: tri_t5_0
13063bdf08b3SMatthew G. Knepley       args: -dim 2 -test_num 5
130781363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
130881363e6fSMatthew G. Knepley 
13093bdf08b3SMatthew G. Knepley     test:
1310c4762a1bSJed Brown       suffix: tet_0
1311c4762a1bSJed Brown       args: -dim 3
131281363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
131381363e6fSMatthew G. Knepley 
1314c4762a1bSJed Brown     test:
1315b253942bSMatthew G. Knepley       suffix: tet_t1_0
1316b253942bSMatthew G. Knepley       args: -dim 3 -test_num 1
131781363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1318b253942bSMatthew G. Knepley 
1319b253942bSMatthew G. Knepley   testset:
1320b253942bSMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1321b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1
1322b253942bSMatthew G. Knepley     test:
1323c4762a1bSJed Brown       suffix: tet_1
1324c4762a1bSJed Brown       nsize: 2
1325c4762a1bSJed Brown       args: -dim 3
132681363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
132781363e6fSMatthew G. Knepley 
1328c4762a1bSJed Brown     test:
1329b253942bSMatthew G. Knepley       suffix: tri_1
1330b253942bSMatthew G. Knepley       nsize: 2
1331b253942bSMatthew G. Knepley       args: -dim 2
133281363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
133381363e6fSMatthew G. Knepley 
13345b9dfbb6SMatthew G. Knepley     test:
13355b9dfbb6SMatthew G. Knepley       suffix: tri_t3_0
13365b9dfbb6SMatthew G. Knepley       nsize: 2
13375b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 3
133881363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
133981363e6fSMatthew G. Knepley 
13405b9dfbb6SMatthew G. Knepley     test:
13415b9dfbb6SMatthew G. Knepley       suffix: tri_t4_0
13425b9dfbb6SMatthew G. Knepley       nsize: 6
13435b9dfbb6SMatthew G. Knepley       args: -dim 2 -test_num 4
134481363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1345c4762a1bSJed Brown 
1346c4762a1bSJed Brown   testset:
1347ecfb78b5SMatthew G. Knepley     args: -orig_dm_plex_check_all -dm_plex_check_all \
1348b7519becSMatthew G. Knepley           -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1
1349c4762a1bSJed Brown     # 2D Quads
1350c4762a1bSJed Brown     test:
1351c4762a1bSJed Brown       suffix: quad_0
1352c4762a1bSJed Brown       args: -dim 2 -cell_simplex 0
135381363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
135481363e6fSMatthew G. Knepley 
1355c4762a1bSJed Brown     test:
1356c4762a1bSJed Brown       suffix: quad_1
1357c4762a1bSJed Brown       nsize: 2
1358c4762a1bSJed Brown       args: -dim 2 -cell_simplex 0
135981363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
136081363e6fSMatthew G. Knepley 
136167de2c8aSMatthew G. Knepley     # Caanot run the assembly test because we cannot orient a fault mesh in a T-shape
1362c4762a1bSJed Brown     test:
1363c4762a1bSJed Brown       suffix: quad_t1_0
136467de2c8aSMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 1 -faulted_dm_plex_check_all -test_assembly 0
136581363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
136681363e6fSMatthew G. Knepley 
13675b9dfbb6SMatthew G. Knepley     test:
13685b9dfbb6SMatthew G. Knepley       suffix: quad_t2_0
13695b9dfbb6SMatthew G. Knepley       nsize: 2
13705b9dfbb6SMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 2
137181363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
137281363e6fSMatthew G. Knepley 
13735b9dfbb6SMatthew G. Knepley     test:
13745b9dfbb6SMatthew G. Knepley       # TODO: The PetscSF is wrong here (connects to wrong side of split)
13755b9dfbb6SMatthew G. Knepley       suffix: quad_t3_0
13765b9dfbb6SMatthew G. Knepley       nsize: 2
13775b9dfbb6SMatthew G. Knepley       args: -dim 2 -cell_simplex 0 -test_num 3
137881363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
137981363e6fSMatthew G. Knepley 
1380c4762a1bSJed Brown     # 3D Hex
1381c4762a1bSJed Brown     test:
1382c4762a1bSJed Brown       suffix: hex_0
1383c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0
138481363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1385c4762a1bSJed Brown     test:
1386c4762a1bSJed Brown       suffix: hex_1
1387c4762a1bSJed Brown       nsize: 2
1388c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0
138981363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1390c4762a1bSJed Brown     test:
1391c4762a1bSJed Brown       suffix: hex_t1_0
1392c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0 -test_num 1
139381363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1394c4762a1bSJed Brown     test:
1395c4762a1bSJed Brown       suffix: hex_t2_0
1396c4762a1bSJed Brown       args: -dim 3 -cell_simplex 0 -test_num 2
139781363e6fSMatthew G. Knepley       filter: sed -e "s/_start//g" -e "s/f0_bd_u_neg//g" -e "s/f0_bd_u_pos//g" -e "s/f0_bd_l//g" -e "s/g0_bd_ul_neg//g" -e "s/g0_bd_ul_pos//g" -e "s/g0_bd_lu//g" -e "s~_ZL.*~~g"
1398c4762a1bSJed Brown 
1399c4762a1bSJed Brown TEST*/
1400