1c4762a1bSJed Brown static char help[] = "Tests for creation of hybrid meshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* TODO 4c4762a1bSJed Brown - Propagate hybridSize with distribution 5c4762a1bSJed Brown - Test with multiple fault segments 6c4762a1bSJed Brown - Test with embedded fault 7c4762a1bSJed Brown - Test with multiple faults 8c4762a1bSJed Brown - Move over all PyLith tests? 9c4762a1bSJed Brown */ 10c4762a1bSJed Brown 11c4762a1bSJed Brown #include <petscdmplex.h> 12b253942bSMatthew G. Knepley #include <petscds.h> 13b253942bSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> 14c4762a1bSJed Brown /* List of test meshes 15c4762a1bSJed Brown 16c4762a1bSJed Brown Triangle 17c4762a1bSJed Brown -------- 18c4762a1bSJed Brown Test 0: 19c4762a1bSJed Brown Two triangles sharing a face 20c4762a1bSJed Brown 21c4762a1bSJed Brown 4 22c4762a1bSJed Brown / | \ 23c4762a1bSJed Brown 8 | 9 24c4762a1bSJed Brown / | \ 25c4762a1bSJed Brown 2 0 7 1 5 26c4762a1bSJed Brown \ | / 27c4762a1bSJed Brown 6 | 10 28c4762a1bSJed Brown \ | / 29c4762a1bSJed Brown 3 30c4762a1bSJed Brown 31c4762a1bSJed Brown should become two triangles separated by a zero-volume cell with 4 vertices 32c4762a1bSJed Brown 33c4762a1bSJed Brown 5--16--8 4--12--6 3 34c4762a1bSJed Brown / | | \ / | | | \ 35c4762a1bSJed Brown 11 | | 12 9 | | | 4 36c4762a1bSJed Brown / | | \ / | | | \ 37c4762a1bSJed Brown 3 0 10 2 14 1 6 2 0 8 1 10 6 0 1 38c4762a1bSJed Brown \ | | / \ | | | / 39c4762a1bSJed Brown 9 | | 13 7 | | | 5 40c4762a1bSJed Brown \ | | / \ | | | / 41c4762a1bSJed Brown 4--15--7 3--11--5 2 42c4762a1bSJed Brown 43c4762a1bSJed Brown Test 1: 44c4762a1bSJed Brown Four triangles sharing two faces which are oriented against each other 45c4762a1bSJed Brown 46c4762a1bSJed Brown 9 47c4762a1bSJed Brown / \ 48c4762a1bSJed Brown / \ 49c4762a1bSJed Brown 17 2 16 50c4762a1bSJed Brown / \ 51c4762a1bSJed Brown / \ 52c4762a1bSJed Brown 8-----15----5 53c4762a1bSJed Brown \ /|\ 54c4762a1bSJed Brown \ / | \ 55c4762a1bSJed Brown 18 3 12 | 14 56c4762a1bSJed Brown \ / | \ 57c4762a1bSJed Brown \ / | \ 58c4762a1bSJed Brown 4 0 11 1 7 59c4762a1bSJed Brown \ | / 60c4762a1bSJed Brown \ | / 61c4762a1bSJed Brown 10 | 13 62c4762a1bSJed Brown \ | / 63c4762a1bSJed Brown \|/ 64c4762a1bSJed Brown 6 65c4762a1bSJed Brown 66c4762a1bSJed Brown Fault mesh 67c4762a1bSJed Brown 68c4762a1bSJed Brown 0 --> 0 69c4762a1bSJed Brown 1 --> 1 70c4762a1bSJed Brown 2 --> 2 71c4762a1bSJed Brown 3 --> 3 72c4762a1bSJed Brown 4 --> 5 73c4762a1bSJed Brown 5 --> 6 74c4762a1bSJed Brown 6 --> 8 75c4762a1bSJed Brown 7 --> 11 76c4762a1bSJed Brown 8 --> 15 77c4762a1bSJed Brown 78c4762a1bSJed Brown 2 79c4762a1bSJed Brown | 80c4762a1bSJed Brown 6----8----4 81c4762a1bSJed Brown | | 82c4762a1bSJed Brown 3 | 83c4762a1bSJed Brown 0-7-1 84c4762a1bSJed Brown | 85c4762a1bSJed Brown | 86c4762a1bSJed Brown 5 87c4762a1bSJed Brown 88c4762a1bSJed Brown should become four triangles separated by two zero-volume cells with 4 vertices 89c4762a1bSJed Brown 90c4762a1bSJed Brown 11 91c4762a1bSJed Brown / \ 92c4762a1bSJed Brown / \ 93c4762a1bSJed Brown / \ 94c4762a1bSJed Brown 22 2 21 95c4762a1bSJed Brown / \ 96c4762a1bSJed Brown / \ 97c4762a1bSJed Brown 10-----20------7 98c4762a1bSJed Brown 28 | 5 26/ \ 99c4762a1bSJed Brown 14----25----12 \ 100c4762a1bSJed Brown \ /| |\ 101c4762a1bSJed Brown \ / | | \ 102c4762a1bSJed Brown 23 3 17 | | 19 103c4762a1bSJed Brown \ / | | \ 104c4762a1bSJed Brown \ / | | \ 105c4762a1bSJed Brown 6 0 24 4 16 1 9 106c4762a1bSJed Brown \ | | / 107c4762a1bSJed Brown \ | | / 108c4762a1bSJed Brown 15 | | 18 109c4762a1bSJed Brown \ | | / 110c4762a1bSJed Brown \| |/ 111c4762a1bSJed Brown 13---8 112c4762a1bSJed Brown 27 113c4762a1bSJed Brown 1145b9dfbb6SMatthew G. Knepley Test 2: 1155b9dfbb6SMatthew G. Knepley Six triangles sharing one face 1165b9dfbb6SMatthew G. Knepley 1175b9dfbb6SMatthew G. Knepley 11-----12------13 1185b9dfbb6SMatthew G. Knepley | /|\ | 1195b9dfbb6SMatthew G. Knepley | 1 / | \ 4 | 1205b9dfbb6SMatthew G. Knepley | / | \ | 1215b9dfbb6SMatthew G. Knepley | / | \ | 1225b9dfbb6SMatthew G. Knepley | / | \ | 1235b9dfbb6SMatthew G. Knepley |/ | \| 1245b9dfbb6SMatthew G. Knepley 9 2 | 5 10 1255b9dfbb6SMatthew G. Knepley |\ | /| 1265b9dfbb6SMatthew G. Knepley | \ | / | 1275b9dfbb6SMatthew G. Knepley | \ | / | 1285b9dfbb6SMatthew G. Knepley | \ | / | 1295b9dfbb6SMatthew G. Knepley | 0 \ | / 3 | 1305b9dfbb6SMatthew G. Knepley | \|/ | 1315b9dfbb6SMatthew G. Knepley 6------7------8 1325b9dfbb6SMatthew G. Knepley 1335b9dfbb6SMatthew G. Knepley Test 3: 1345b9dfbb6SMatthew G. Knepley This is Test 2 on two processes. After the fault, we have 1355b9dfbb6SMatthew G. Knepley 1365b9dfbb6SMatthew G. Knepley 6--12--7 7--20-10--16-8 1375b9dfbb6SMatthew G. Knepley | /| | |\ | 1385b9dfbb6SMatthew G. Knepley | 1 / | | | \ 1 | 1395b9dfbb6SMatthew G. Knepley 13 11 | | | 17 15 1405b9dfbb6SMatthew G. Knepley | / | | | \ | 1415b9dfbb6SMatthew G. Knepley | / | | | \ | 1425b9dfbb6SMatthew G. Knepley |/ | | | \| 1435b9dfbb6SMatthew G. Knepley 5 2 14 11 3 18 2 6 1445b9dfbb6SMatthew G. Knepley |\ | | | /| 1455b9dfbb6SMatthew G. Knepley | \ | | | / | 1465b9dfbb6SMatthew G. Knepley | \ | | | / | 1475b9dfbb6SMatthew G. Knepley 10 9 | | | 14 13 1485b9dfbb6SMatthew G. Knepley | 0 \ | | | / 0 | 1495b9dfbb6SMatthew G. Knepley | \| | |/ | 1505b9dfbb6SMatthew G. Knepley 3---8--4 4--19-9--12--5 1515b9dfbb6SMatthew G. Knepley 1525b9dfbb6SMatthew G. Knepley Test 4: 1535b9dfbb6SMatthew G. Knepley This is Test 2 on six processes. After the fault, we have 1545b9dfbb6SMatthew G. Knepley 1553bdf08b3SMatthew G. Knepley Test 5: 1563bdf08b3SMatthew G. Knepley 1573bdf08b3SMatthew G. Knepley Fault only on points 2 and 5: 1583bdf08b3SMatthew G. Knepley 1593bdf08b3SMatthew G. Knepley 6 1603bdf08b3SMatthew G. Knepley / | \ 1613bdf08b3SMatthew G. Knepley 13 | 17 1623bdf08b3SMatthew G. Knepley / 15 \ 1633bdf08b3SMatthew G. Knepley 7 0 | 1 9 1643bdf08b3SMatthew G. Knepley |\ | /| 1653bdf08b3SMatthew G. Knepley | 14 | 16 | 1663bdf08b3SMatthew G. Knepley | \ | / | 1673bdf08b3SMatthew G. Knepley 18| 2 8 3 |21 1683bdf08b3SMatthew G. Knepley | / | \ | 1693bdf08b3SMatthew G. Knepley | 19 | 20 | 1703bdf08b3SMatthew G. Knepley |/ | \| 1713bdf08b3SMatthew G. Knepley 10 4 | 5 12 1723bdf08b3SMatthew G. Knepley \ 23 / 1733bdf08b3SMatthew G. Knepley 22 | 24 1743bdf08b3SMatthew G. Knepley \ | / 1753bdf08b3SMatthew G. Knepley 11 1763bdf08b3SMatthew G. Knepley 177c4762a1bSJed Brown Tetrahedron 178c4762a1bSJed Brown ----------- 179c4762a1bSJed Brown Test 0: 180c4762a1bSJed Brown Two tets sharing a face 181c4762a1bSJed Brown 182c4762a1bSJed Brown cell 5 _______ cell 183c4762a1bSJed Brown 0 / | \ \ 1 184c4762a1bSJed Brown 16 | 18 22 185c4762a1bSJed Brown /8 19 10\ \ 186c4762a1bSJed Brown 2-15-|----4--21--6 187c4762a1bSJed Brown \ 9| 7 / / 188c4762a1bSJed Brown 14 | 17 20 189c4762a1bSJed Brown \ | / / 190c4762a1bSJed Brown 3------- 191c4762a1bSJed Brown 192c4762a1bSJed Brown should become two tetrahedrons separated by a zero-volume cell with 3 faces/3 edges/6 vertices 193c4762a1bSJed Brown 194c4762a1bSJed Brown cell 6 ___36___10______ cell 195c4762a1bSJed Brown 0 / | \ |\ \ 1 196c4762a1bSJed Brown 24 | 26 | 32 30 197c4762a1bSJed Brown /12 27 14\ 33 \ \ 198c4762a1bSJed Brown 3-23-|----5--35-|---9--29--7 199c4762a1bSJed Brown \ 13| 11/ |18 / / 200c4762a1bSJed Brown 22 | 25 | 31 28 201c4762a1bSJed Brown \ | / |/ / 202c4762a1bSJed Brown 4----34----8------ 203c4762a1bSJed Brown cell 2 204c4762a1bSJed Brown 205c4762a1bSJed Brown In parallel, 206c4762a1bSJed Brown 207c4762a1bSJed Brown cell 5 ___28____8 4______ cell 208c4762a1bSJed Brown 0 / | \ |\ |\ \ 0 209c4762a1bSJed Brown 19 | 21 | 24 | 13 6 11 210c4762a1bSJed Brown /10 22 12\ 25 \ |8 \ \ 211c4762a1bSJed Brown 2-18-|----4--27-|---7 14 3--10--1 212c4762a1bSJed Brown \ 11| 9 / |13 / | / / 213c4762a1bSJed Brown 17 | 20 | 23 | 12 5 9 214c4762a1bSJed Brown \ | / |/ |/ / 215c4762a1bSJed Brown 3----26----6 2------ 216c4762a1bSJed Brown cell 1 217c4762a1bSJed Brown 218c4762a1bSJed Brown Test 1: 219c4762a1bSJed Brown Four tets sharing two faces 220c4762a1bSJed Brown 221c4762a1bSJed Brown Cells: 0-3,4-5 222c4762a1bSJed Brown Vertices: 6-15 223c4762a1bSJed Brown Faces: 16-29,30-34 224c4762a1bSJed Brown Edges: 35-52,53-56 225c4762a1bSJed Brown 226c4762a1bSJed Brown Quadrilateral 227c4762a1bSJed Brown ------------- 228c4762a1bSJed Brown Test 0: 229c4762a1bSJed Brown Two quads sharing a face 230c4762a1bSJed Brown 231c4762a1bSJed Brown 5--10---4--14---7 232c4762a1bSJed Brown | | | 233c4762a1bSJed Brown 11 0 9 1 13 234c4762a1bSJed Brown | | | 235c4762a1bSJed Brown 2---8---3--12---6 236c4762a1bSJed Brown 237c4762a1bSJed Brown should become two quads separated by a zero-volume cell with 4 vertices 238c4762a1bSJed Brown 239c4762a1bSJed Brown 6--13---5-20-10--17---8 5--10---4-14--7 4---7---2 240c4762a1bSJed Brown | | | | | | | | | 241c4762a1bSJed Brown 14 0 12 2 18 1 16 11 0 9 1 12 8 0 6 242c4762a1bSJed Brown | | | | | | | | | 243c4762a1bSJed Brown 3--11---4-19--9--15---7 2---8---3-13--6 3---5---1 244c4762a1bSJed Brown 245c4762a1bSJed Brown Test 1: 246c4762a1bSJed Brown 247c4762a1bSJed Brown Original mesh with 9 cells, 248c4762a1bSJed Brown 2495b9dfbb6SMatthew G. Knepley 9-----10-----11-----12 2505b9dfbb6SMatthew G. Knepley | | || | 2515b9dfbb6SMatthew G. Knepley | | || | 2525b9dfbb6SMatthew G. Knepley | 0 | 1 || 2 | 2535b9dfbb6SMatthew G. Knepley | | || | 2545b9dfbb6SMatthew G. Knepley 13-----14-----15-----16 2555b9dfbb6SMatthew G. Knepley | | || | 2565b9dfbb6SMatthew G. Knepley | | || | 2575b9dfbb6SMatthew G. Knepley | 3 | 4 || 5 | 2585b9dfbb6SMatthew G. Knepley | | || | 2595b9dfbb6SMatthew G. Knepley 17-----18-----19=====20 260c4762a1bSJed Brown | | | | 261c4762a1bSJed Brown | | | | 2625b9dfbb6SMatthew G. Knepley | 6 | 7 | 8 | 263c4762a1bSJed Brown | | | | 2645b9dfbb6SMatthew G. Knepley 21-----22-----23-----24 265c4762a1bSJed Brown 266c4762a1bSJed Brown After first fault, 267c4762a1bSJed Brown 268c4762a1bSJed Brown 12 ----13 ----14-28 ----15 269c4762a1bSJed Brown | | | | | 270c4762a1bSJed Brown | 0 | 1 | 9| 2 | 271c4762a1bSJed Brown | | | | | 272c4762a1bSJed Brown | | | | | 273c4762a1bSJed Brown 16 ----17 ----18-29 ----19 274c4762a1bSJed Brown | | | | | 275c4762a1bSJed Brown | 3 | 4 |10| 5 | 276c4762a1bSJed Brown | | | | | 277c4762a1bSJed Brown | | | | | 278c4762a1bSJed Brown 20 ----21-----22-30 ----23 279c4762a1bSJed Brown | | | \--11- | 280c4762a1bSJed Brown | 6 | 7 | 8 | 281c4762a1bSJed Brown | | | | 282c4762a1bSJed Brown | | | | 283c4762a1bSJed Brown 24 ----25 ----26--------27 284c4762a1bSJed Brown 285c4762a1bSJed Brown After second fault, 286c4762a1bSJed Brown 287c4762a1bSJed Brown 14 ----15 ----16-30 ----17 288c4762a1bSJed Brown | | | | | 289c4762a1bSJed Brown | 0 | 1 | 9| 2 | 290c4762a1bSJed Brown | | | | | 291c4762a1bSJed Brown | | | | | 292c4762a1bSJed Brown 18 ----19 ----20-31 ----21 293c4762a1bSJed Brown | | | | | 294c4762a1bSJed Brown | 3 | 4 |10| 5 | 295c4762a1bSJed Brown | | | | | 296c4762a1bSJed Brown | | | | | 297c4762a1bSJed Brown 33 ----34-----24-32 ----25 298c4762a1bSJed Brown | 12 | 13 / | \-11-- | 299c4762a1bSJed Brown 22 ----23---/ | | 3005b9dfbb6SMatthew G. Knepley | | | | 3015b9dfbb6SMatthew G. Knepley | 6 | 7 | 8 | 302c4762a1bSJed Brown | | | | 303c4762a1bSJed Brown | | | | 304c4762a1bSJed Brown 26 ----27 ----28--------29 305c4762a1bSJed Brown 3065b9dfbb6SMatthew G. Knepley Test 2: 3075b9dfbb6SMatthew G. Knepley Two quads sharing a face in parallel 3085b9dfbb6SMatthew G. Knepley 3095b9dfbb6SMatthew G. Knepley 4---7---3 2---8---4 3105b9dfbb6SMatthew G. Knepley | | | | 3115b9dfbb6SMatthew G. Knepley 8 0 6 5 0 7 3125b9dfbb6SMatthew G. Knepley | | | | 3135b9dfbb6SMatthew G. Knepley 1---5---2 1---6---3 3145b9dfbb6SMatthew G. Knepley 3155b9dfbb6SMatthew G. Knepley should become two quads separated by a zero-volume cell with 4 vertices 3165b9dfbb6SMatthew G. Knepley 3175b9dfbb6SMatthew G. Knepley 4---7---3 3-14--7--11---5 3185b9dfbb6SMatthew G. Knepley | | | | | 3195b9dfbb6SMatthew G. Knepley 8 0 6 8 1 12 0 10 3205b9dfbb6SMatthew G. Knepley | | | | | 3215b9dfbb6SMatthew G. Knepley 1---5---2 2-13--6---9---4 3225b9dfbb6SMatthew G. Knepley 3235b9dfbb6SMatthew G. Knepley Test 3: 3245b9dfbb6SMatthew G. Knepley Like Test 2, but with different partition 3255b9dfbb6SMatthew G. Knepley 3265b9dfbb6SMatthew G. Knepley 5--10---4-14--7 2---8---4 3275b9dfbb6SMatthew G. Knepley | | | | | 3285b9dfbb6SMatthew G. Knepley 11 0 9 1 12 5 0 7 3295b9dfbb6SMatthew G. Knepley | | | | | 3305b9dfbb6SMatthew G. Knepley 2---8---3-13--6 1---6---3 3315b9dfbb6SMatthew G. Knepley 332c4762a1bSJed Brown Hexahedron 333c4762a1bSJed Brown ---------- 334c4762a1bSJed Brown Test 0: 335c4762a1bSJed Brown Two hexes sharing a face 336c4762a1bSJed Brown 337c4762a1bSJed Brown cell 9-----31------8-----42------13 cell 338c4762a1bSJed Brown 0 /| /| /| 1 339c4762a1bSJed Brown 32 | 15 30| 21 41| 340c4762a1bSJed Brown / | / | / | 341c4762a1bSJed Brown 6-----29------7-----40------12 | 342c4762a1bSJed Brown | | 18 | | 24 | | 343c4762a1bSJed Brown | 36 | 35 | 44 344c4762a1bSJed Brown |19 | |17 | |23 | 345c4762a1bSJed Brown 33 | 16 34 | 22 43 | 346c4762a1bSJed Brown | 5-----27--|---4-----39--|---11 347c4762a1bSJed Brown | / | / | / 348c4762a1bSJed Brown | 28 14 | 26 20 | 38 349c4762a1bSJed Brown |/ |/ |/ 350c4762a1bSJed Brown 2-----25------3-----37------10 351c4762a1bSJed Brown 352c4762a1bSJed Brown should become two hexes separated by a zero-volume cell with 8 vertices 353c4762a1bSJed Brown 354c4762a1bSJed Brown cell 2 355c4762a1bSJed Brown cell 10-----41------9-----62------18----52------14 cell 356c4762a1bSJed Brown 0 /| /| /| /| 1 357c4762a1bSJed Brown 42 | 20 40| 32 56| 26 51| 358c4762a1bSJed Brown / | / | / | / | 359c4762a1bSJed Brown 7-----39------8-----61------17--|-50------13 | 360c4762a1bSJed Brown | | 23 | | | | 29 | | 361c4762a1bSJed Brown | 46 | 45 | 58 | 54 362c4762a1bSJed Brown |24 | |22 | |30 | |28 | 363c4762a1bSJed Brown 43 | 21 44 | 57 | 27 53 | 364c4762a1bSJed Brown | 6-----37--|---5-----60--|---16----49--|---12 365c4762a1bSJed Brown | / | / | / | / 366c4762a1bSJed Brown | 38 19 | 36 31 | 55 25 | 48 367c4762a1bSJed Brown |/ |/ |/ |/ 368c4762a1bSJed Brown 3-----35------4-----59------15----47------11 369c4762a1bSJed Brown 370c4762a1bSJed Brown In parallel, 371c4762a1bSJed Brown 372c4762a1bSJed Brown cell 2 373c4762a1bSJed Brown cell 9-----31------8-----44------13 8----20------4 cell 374c4762a1bSJed Brown 0 /| /| /| /| /| 1 375c4762a1bSJed Brown 32 | 15 30| 22 38| 24 | 10 19| 376c4762a1bSJed Brown / | / | / | / | / | 377c4762a1bSJed Brown 6-----29------7-----43------12 | 7----18------3 | 378c4762a1bSJed Brown | | 18 | | | | | | 13 | | 379c4762a1bSJed Brown | 36 | 35 | 40 | 26 | 22 380c4762a1bSJed Brown |19 | |17 | |20 | |14 | |12 | 381c4762a1bSJed Brown 33 | 16 34 | 39 | 25 | 11 21 | 382c4762a1bSJed Brown | 5-----27--|---4-----42--|---11 | 6----17--|---2 383c4762a1bSJed Brown | / | / | / | / | / 384c4762a1bSJed Brown | 28 14 | 26 21 | 37 |23 9 | 16 385c4762a1bSJed Brown |/ |/ |/ |/ |/ 386c4762a1bSJed Brown 2-----25------3-----41------10 5----15------1 387c4762a1bSJed Brown 388c4762a1bSJed Brown Test 1: 389c4762a1bSJed Brown 390c4762a1bSJed Brown */ 391c4762a1bSJed Brown 392c4762a1bSJed Brown typedef struct { 393c4762a1bSJed Brown PetscInt debug; /* The debugging level */ 394c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 395c4762a1bSJed Brown PetscBool cellSimplex; /* Use simplices or hexes */ 396c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 397c4762a1bSJed Brown PetscInt testNum; /* The particular mesh to test */ 398b7519becSMatthew G. Knepley PetscInt cohesiveFields; /* The number of cohesive fields */ 399c4762a1bSJed Brown } AppCtx; 400c4762a1bSJed Brown 401d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 402d71ae5a4SJacob Faibussowitsch { 403c4762a1bSJed Brown PetscFunctionBegin; 404c4762a1bSJed Brown options->debug = 0; 405c4762a1bSJed Brown options->dim = 2; 406c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 407c4762a1bSJed Brown options->testPartition = PETSC_TRUE; 408c4762a1bSJed Brown options->testNum = 0; 409b7519becSMatthew G. Knepley options->cohesiveFields = 1; 410c4762a1bSJed Brown 411d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 4129566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex5.c", options->debug, &options->debug, NULL, 0)); 4139566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex5.c", options->dim, &options->dim, NULL, 1, 3)); 4149566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex5.c", options->cellSimplex, &options->cellSimplex, NULL)); 4159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex5.c", options->testPartition, &options->testPartition, NULL)); 4169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex5.c", options->testNum, &options->testNum, NULL, 0)); 4179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-cohesive_fields", "The number of cohesive fields", "ex5.c", options->cohesiveFields, &options->cohesiveFields, NULL, 0)); 418d0609cedSBarry Smith PetscOptionsEnd(); 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 420c4762a1bSJed Brown } 421c4762a1bSJed Brown 422d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 423d71ae5a4SJacob Faibussowitsch { 424c4762a1bSJed Brown DM idm; 425c4762a1bSJed Brown PetscInt p; 426c4762a1bSJed Brown PetscMPIInt rank; 427c4762a1bSJed Brown 428c4762a1bSJed Brown PetscFunctionBegin; 4299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 430dd400576SPatrick Sanan if (rank == 0) { 431c4762a1bSJed Brown switch (testNum) { 4329371c9d4SSatish Balay case 0: { 433c4762a1bSJed Brown PetscInt numPoints[2] = {4, 2}; 434c4762a1bSJed Brown PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 435c4762a1bSJed Brown PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 436c4762a1bSJed Brown PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 437c4762a1bSJed Brown PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 438c4762a1bSJed Brown PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 439c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 440c4762a1bSJed Brown 4419566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 4429566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4439566063dSJacob Faibussowitsch for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 4449566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 4459566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 1, 2)); 4469371c9d4SSatish Balay } break; 4479371c9d4SSatish Balay case 1: { 448c4762a1bSJed Brown PetscInt numPoints[2] = {6, 4}; 449c4762a1bSJed Brown PetscInt coneSize[10] = {3, 3, 3, 3, 0, 0, 0, 0, 0, 0}; 450c4762a1bSJed Brown PetscInt cones[12] = {4, 6, 5, 5, 6, 7, 8, 5, 9, 8, 4, 5}; 451c4762a1bSJed Brown PetscInt coneOrientations[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 452c4762a1bSJed 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}; 453c4762a1bSJed Brown PetscInt markerPoints[6] = {4, 1, 6, 1, 8, 1}; 454c4762a1bSJed Brown PetscInt faultPoints[3] = {5, 6, 8}; 455c4762a1bSJed Brown 4569566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 4579566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4589566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 4599371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 4609371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 3, 1)); 4619371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 1, 2)); 4629371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 2, 2)); 4639371c9d4SSatish Balay } break; 4645b9dfbb6SMatthew G. Knepley case 2: 4655b9dfbb6SMatthew G. Knepley case 3: 4669371c9d4SSatish Balay case 4: { 4675b9dfbb6SMatthew G. Knepley PetscInt numPoints[2] = {8, 6}; 4685b9dfbb6SMatthew G. Knepley PetscInt coneSize[14] = {3, 3, 3, 3, 3, 3, 0, 0, 0, 0, 0, 0, 0, 0}; 4699371c9d4SSatish Balay PetscInt cones[18] = {6, 7, 9, 9, 12, 11, 7, 12, 9, 7, 8, 10, 10, 13, 12, 7, 10, 12}; 4705b9dfbb6SMatthew G. Knepley PetscInt coneOrientations[18] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 4719371c9d4SSatish Balay PetscScalar vertexCoords[16] = { 4729371c9d4SSatish Balay -1., -1., 0., -1., 1., -1., -1., 0., 1., 0., -1., 1., 0., 1., 1., 1., 4739371c9d4SSatish Balay }; 4745b9dfbb6SMatthew G. Knepley PetscInt markerPoints[16] = {6, 1, 7, 1, 8, 1, 9, 1, 10, 1, 11, 1, 12, 1, 13, 1}; 4755b9dfbb6SMatthew G. Knepley PetscInt faultPoints[2] = {7, 12}; 4765b9dfbb6SMatthew G. Knepley 4775b9dfbb6SMatthew G. Knepley PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 4785b9dfbb6SMatthew G. Knepley PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 4795b9dfbb6SMatthew G. Knepley PetscCall(DMSetLabelValue(*dm, "material", 1, 1)); 4805b9dfbb6SMatthew G. Knepley PetscCall(DMSetLabelValue(*dm, "material", 2, 1)); 4815b9dfbb6SMatthew G. Knepley PetscCall(DMSetLabelValue(*dm, "material", 3, 2)); 4825b9dfbb6SMatthew G. Knepley PetscCall(DMSetLabelValue(*dm, "material", 4, 2)); 4835b9dfbb6SMatthew G. Knepley PetscCall(DMSetLabelValue(*dm, "material", 5, 2)); 4845b9dfbb6SMatthew G. Knepley for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4859371c9d4SSatish Balay if (testNum == 2) 4869371c9d4SSatish Balay for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 4879371c9d4SSatish Balay if (testNum == 3 || testNum == 4) 4889371c9d4SSatish Balay for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1)); 4899371c9d4SSatish Balay } break; 4909371c9d4SSatish Balay case 5: { 4913bdf08b3SMatthew G. Knepley PetscInt numPoints[2] = {7, 6}; 4923bdf08b3SMatthew G. Knepley PetscInt coneSize[13] = {3, 3, 3, 3, 3, 3, 0, 0, 0, 0, 0, 0, 0}; 4933bdf08b3SMatthew G. Knepley PetscInt cones[18] = {6, 7, 8, 8, 9, 6, 7, 10, 8, 9, 8, 12, 8, 10, 11, 11, 12, 8}; 4943bdf08b3SMatthew G. Knepley PetscInt coneOrientations[18] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 4953bdf08b3SMatthew 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}; 4963bdf08b3SMatthew G. Knepley PetscInt faultPoints[2] = {8, 11}; 4973bdf08b3SMatthew G. Knepley PetscInt faultBdPoints[1] = {8}; 4983bdf08b3SMatthew G. Knepley 4993bdf08b3SMatthew G. Knepley PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 5003bdf08b3SMatthew G. Knepley for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 5013bdf08b3SMatthew G. Knepley for (p = 0; p < 1; ++p) PetscCall(DMSetLabelValue(*dm, "faultBd", faultBdPoints[p], 1)); 5029371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 0, 0)); 5039371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 2, 0)); 5049371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 4, 0)); 5059371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 1, 2)); 5069371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 3, 2)); 5079371c9d4SSatish Balay PetscCall(DMSetLabelValue(*dm, "material", 5, 2)); 5089371c9d4SSatish Balay } break; 509d71ae5a4SJacob Faibussowitsch default: 510d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 511c4762a1bSJed Brown } 512c4762a1bSJed Brown } else { 513c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 514c4762a1bSJed Brown 5159566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 5165b9dfbb6SMatthew G. Knepley if (testNum == 3 || testNum == 4) PetscCall(DMCreateLabel(*dm, "pfault")); 5175b9dfbb6SMatthew G. Knepley else PetscCall(DMCreateLabel(*dm, "fault")); 518c4762a1bSJed Brown } 5199566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 5209566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view")); 5219566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 522c4762a1bSJed Brown *dm = idm; 5233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 524c4762a1bSJed Brown } 525c4762a1bSJed Brown 526d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, AppCtx *user, DM dm) 527d71ae5a4SJacob Faibussowitsch { 528c4762a1bSJed Brown PetscInt depth = 3, testNum = user->testNum, p; 529c4762a1bSJed Brown PetscMPIInt rank; 530c4762a1bSJed Brown 531c4762a1bSJed Brown PetscFunctionBegin; 5329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 533dd400576SPatrick Sanan if (rank == 0) { 534c4762a1bSJed Brown switch (testNum) { 5359371c9d4SSatish Balay case 0: { 536c4762a1bSJed Brown PetscInt numPoints[4] = {5, 7, 9, 2}; 537c4762a1bSJed 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}; 538c4762a1bSJed 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}; 539b5a892a1SMatthew 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}; 540c4762a1bSJed 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}; 541c4762a1bSJed Brown PetscInt markerPoints[20] = {2, 1, 3, 1, 4, 1, 5, 1, 14, 1, 15, 1, 16, 1, 17, 1, 18, 1, 19, 1}; 542c4762a1bSJed Brown PetscInt faultPoints[3] = {3, 4, 5}; 543c4762a1bSJed Brown 5449566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 54548a46eb9SPierre Jolivet for (p = 0; p < 10; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 54648a46eb9SPierre Jolivet for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1)); 5479566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "material", 0, 1)); 5489566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "material", 1, 2)); 5499371c9d4SSatish Balay } break; 5509371c9d4SSatish Balay case 1: { 551c4762a1bSJed Brown PetscInt numPoints[4] = {6, 13, 12, 4}; 552c4762a1bSJed 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}; 5539371c9d4SSatish 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, 5549371c9d4SSatish 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}; 5559371c9d4SSatish 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, 5569371c9d4SSatish 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}; 557c4762a1bSJed 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}; 558c4762a1bSJed Brown PetscInt markerPoints[14] = {5, 1, 6, 1, 7, 1, 10, 1, 22, 1, 23, 1, 24, 1}; 559c4762a1bSJed Brown PetscInt faultPoints[4] = {5, 6, 7, 8}; 560c4762a1bSJed Brown 5619566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 56248a46eb9SPierre Jolivet for (p = 0; p < 7; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 56348a46eb9SPierre Jolivet for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "fault", faultPoints[p], 1)); 5649566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "material", 0, 1)); 5659566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "material", 1, 2)); 5669371c9d4SSatish Balay } break; 567d71ae5a4SJacob Faibussowitsch default: 568d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 569c4762a1bSJed Brown } 570c4762a1bSJed Brown } else { 571c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 572c4762a1bSJed Brown 5739566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL)); 5749566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "fault")); 575c4762a1bSJed Brown } 5763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 577c4762a1bSJed Brown } 578c4762a1bSJed Brown 579d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM *dm) 580d71ae5a4SJacob Faibussowitsch { 581c4762a1bSJed Brown DM idm; 582c4762a1bSJed Brown PetscInt p; 583c4762a1bSJed Brown PetscMPIInt rank; 584c4762a1bSJed Brown 585c4762a1bSJed Brown PetscFunctionBegin; 5869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 587dd400576SPatrick Sanan if (rank == 0) { 588c4762a1bSJed Brown switch (testNum) { 589c4762a1bSJed Brown case 0: 590c4762a1bSJed Brown case 2: 5919371c9d4SSatish Balay case 3: { 592c4762a1bSJed Brown PetscInt numPoints[2] = {6, 2}; 593c4762a1bSJed Brown PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 594c4762a1bSJed Brown PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 595c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 596c4762a1bSJed 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}; 597c4762a1bSJed Brown PetscInt markerPoints[12] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1}; 598c4762a1bSJed Brown PetscInt faultPoints[2] = {3, 4}; 599c4762a1bSJed Brown 6009566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 6019566063dSJacob Faibussowitsch for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 6029371c9d4SSatish Balay if (testNum == 0) 6039371c9d4SSatish Balay for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 6049371c9d4SSatish Balay if (testNum == 2 || testNum == 3) 6059371c9d4SSatish Balay for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "pfault", faultPoints[p], 1)); 6069566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 6079566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 1, 2)); 6089371c9d4SSatish Balay } break; 6099371c9d4SSatish Balay case 1: { 610c4762a1bSJed Brown PetscInt numPoints[2] = {16, 9}; 611c4762a1bSJed 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}; 6129371c9d4SSatish 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}; 613c4762a1bSJed 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}; 6149371c9d4SSatish 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}; 6155b9dfbb6SMatthew G. Knepley PetscInt faultPoints[4] = {11, 15, 19, 20}; 616c4762a1bSJed Brown PetscInt fault2Points[2] = {17, 18}; 617c4762a1bSJed Brown 6189566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 6195b9dfbb6SMatthew G. Knepley for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1)); 6205b9dfbb6SMatthew G. Knepley for (p = 3; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "faultBd", faultPoints[p], 1)); 6219566063dSJacob Faibussowitsch for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault2", fault2Points[p], 1)); 6229566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 6239566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 1, 1)); 6249566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 2, 1)); 6259566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 3, 1)); 6269566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 4, 1)); 6279566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 5, 2)); 6289566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 6, 2)); 6299566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 7, 2)); 6309566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 8, 2)); 6319371c9d4SSatish Balay } break; 632d71ae5a4SJacob Faibussowitsch default: 633d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 634c4762a1bSJed Brown } 635c4762a1bSJed Brown } else { 636c4762a1bSJed Brown PetscInt numPoints[3] = {0, 0, 0}; 637c4762a1bSJed Brown 6389566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL)); 6395b9dfbb6SMatthew G. Knepley if (testNum == 2 || testNum == 3) PetscCall(DMCreateLabel(*dm, "pfault")); 6409566063dSJacob Faibussowitsch else PetscCall(DMCreateLabel(*dm, "fault")); 641c4762a1bSJed Brown } 6429566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 6439566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view")); 6449566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 645c4762a1bSJed Brown *dm = idm; 6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 647c4762a1bSJed Brown } 648c4762a1bSJed Brown 649d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscInt testNum, DM *dm) 650d71ae5a4SJacob Faibussowitsch { 651c4762a1bSJed Brown DM idm; 652c4762a1bSJed Brown PetscInt depth = 3, p; 653c4762a1bSJed Brown PetscMPIInt rank; 654c4762a1bSJed Brown 655c4762a1bSJed Brown PetscFunctionBegin; 6569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 657dd400576SPatrick Sanan if (rank == 0) { 658c4762a1bSJed Brown switch (testNum) { 6599371c9d4SSatish Balay case 0: { 660c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 661c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 662c4762a1bSJed Brown PetscInt cones[16] = {2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8}; 663c4762a1bSJed Brown PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 6649371c9d4SSatish 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}; 665c4762a1bSJed Brown PetscInt markerPoints[52] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1}; 666c4762a1bSJed Brown PetscInt faultPoints[4] = {3, 4, 7, 8}; 667c4762a1bSJed Brown 6689566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 6699566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 6709566063dSJacob Faibussowitsch for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(idm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 6719566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1)); 6729566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 6739566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 1, 2)); 6749371c9d4SSatish Balay } break; 6759371c9d4SSatish Balay case 1: { 676c4762a1bSJed Brown /* Cell Adjacency Graph: 677c4762a1bSJed Brown 0 -- { 8, 13, 21, 24} --> 1 678c4762a1bSJed Brown 0 -- {20, 21, 23, 24} --> 5 F 679c4762a1bSJed Brown 1 -- {10, 15, 21, 24} --> 2 680c4762a1bSJed Brown 1 -- {13, 14, 15, 24} --> 6 681c4762a1bSJed Brown 2 -- {21, 22, 24, 25} --> 4 F 682c4762a1bSJed Brown 3 -- {21, 24, 30, 35} --> 4 683c4762a1bSJed Brown 3 -- {21, 24, 28, 33} --> 5 684c4762a1bSJed Brown */ 685c4762a1bSJed Brown PetscInt numPoints[2] = {30, 7}; 686c4762a1bSJed 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}; 6879371c9d4SSatish 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}; 688c4762a1bSJed 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}; 6899371c9d4SSatish 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, 6909371c9d4SSatish 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, 6919371c9d4SSatish 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}; 692c4762a1bSJed Brown PetscInt faultPoints[6] = {20, 21, 22, 23, 24, 25}; 693c4762a1bSJed Brown 6949566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 6959566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 6969566063dSJacob Faibussowitsch for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1)); 6979566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 6989566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 1, 1)); 6999566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 2, 1)); 7009566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 3, 2)); 7019566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 4, 2)); 7029566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 5, 2)); 7039566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 6, 2)); 7049371c9d4SSatish Balay } break; 7059371c9d4SSatish Balay case 2: { 706c4762a1bSJed Brown /* Buried fault edge */ 707c4762a1bSJed Brown PetscInt numPoints[2] = {18, 4}; 708c4762a1bSJed 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}; 7099371c9d4SSatish 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}; 710c4762a1bSJed Brown PetscInt coneOrientations[32] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 7119371c9d4SSatish 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, 7129371c9d4SSatish 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}; 713c4762a1bSJed Brown PetscInt faultPoints[4] = {7, 8, 16, 17}; 714c4762a1bSJed Brown 7159566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 7169566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 7179566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(idm, "fault", faultPoints[p], 1)); 7189566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 0, 1)); 7199566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 1, 1)); 7209566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 2, 2)); 7219566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "material", 3, 2)); 7229371c9d4SSatish Balay } break; 723d71ae5a4SJacob Faibussowitsch default: 724d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 725c4762a1bSJed Brown } 726c4762a1bSJed Brown } else { 727c4762a1bSJed Brown PetscInt numPoints[4] = {0, 0, 0, 0}; 728c4762a1bSJed Brown 7299566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL)); 7309566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 7319566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(idm, "fault")); 732c4762a1bSJed Brown } 7339566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view")); 7349566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 735c4762a1bSJed Brown *dm = idm; 7363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 737c4762a1bSJed Brown } 738c4762a1bSJed Brown 739d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFaultLabel(DM dm) 740d71ae5a4SJacob Faibussowitsch { 741b253942bSMatthew G. Knepley DMLabel label; 742b253942bSMatthew G. Knepley PetscInt dim, h, pStart, pEnd, pMax, p; 743b253942bSMatthew G. Knepley 744b253942bSMatthew G. Knepley PetscFunctionBegin; 7459566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7469566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "cohesive")); 7479566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "cohesive", &label)); 748b253942bSMatthew G. Knepley for (h = 0; h <= dim; ++h) { 7499566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, h, NULL, &pMax)); 7509566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h, &pStart, &pEnd)); 7519566063dSJacob Faibussowitsch for (p = pMax; p < pEnd; ++p) PetscCall(DMLabelSetValue(label, p, 1)); 752b253942bSMatthew G. Knepley } 7533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 754b253942bSMatthew G. Knepley } 755b253942bSMatthew G. Knepley 756d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateDiscretization(DM dm, AppCtx *user) 757d71ae5a4SJacob Faibussowitsch { 758b253942bSMatthew G. Knepley PetscFE fe; 759b253942bSMatthew G. Knepley DMLabel fault; 760b7519becSMatthew G. Knepley PetscInt dim, Ncf = user->cohesiveFields, f; 761b253942bSMatthew G. Knepley 762b253942bSMatthew G. Knepley PetscFunctionBegin; 7639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7649566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "cohesive", &fault)); 7659566063dSJacob Faibussowitsch PetscCall(DMLabelView(fault, PETSC_VIEWER_STDOUT_WORLD)); 766b253942bSMatthew G. Knepley 7679566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, user->cellSimplex, "displacement_", PETSC_DETERMINE, &fe)); 7689566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "displacement")); 7699566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 7709566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 771b253942bSMatthew G. Knepley 772b7519becSMatthew G. Knepley if (Ncf > 0) { 7739566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim - 1, dim, user->cellSimplex, "faulttraction_", PETSC_DETERMINE, &fe)); 7749566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "fault traction")); 7759566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, fault, (PetscObject)fe)); 7769566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 777b7519becSMatthew G. Knepley } 778b7519becSMatthew G. Knepley for (f = 1; f < Ncf; ++f) { 779b7519becSMatthew G. Knepley char name[256], opt[256]; 780b7519becSMatthew G. Knepley 78163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, 256, "fault field %" PetscInt_FMT, f)); 78263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(opt, 256, "faultfield_%" PetscInt_FMT "_", f)); 7839566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim - 1, dim, user->cellSimplex, opt, PETSC_DETERMINE, &fe)); 7849566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, name)); 7859566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, fault, (PetscObject)fe)); 7869566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 787b7519becSMatthew G. Knepley } 788b253942bSMatthew G. Knepley 7899566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 7903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 791b253942bSMatthew G. Knepley } 792b253942bSMatthew G. Knepley 793d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 794d71ae5a4SJacob Faibussowitsch { 795c4762a1bSJed Brown PetscInt dim = user->dim; 796c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex, hasFault, hasFault2, hasParallelFault; 797c4762a1bSJed Brown PetscMPIInt rank, size; 798dd96b1f9SToby Isaac DMLabel matLabel; 799c4762a1bSJed Brown 800c4762a1bSJed Brown PetscFunctionBegin; 8019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 8039566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 8049566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 8059566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 806c4762a1bSJed Brown switch (dim) { 807c4762a1bSJed Brown case 2: 808c4762a1bSJed Brown if (cellSimplex) { 8099566063dSJacob Faibussowitsch PetscCall(CreateSimplex_2D(comm, user->testNum, dm)); 810c4762a1bSJed Brown } else { 8119566063dSJacob Faibussowitsch PetscCall(CreateQuad_2D(comm, user->testNum, dm)); 812c4762a1bSJed Brown } 813c4762a1bSJed Brown break; 814c4762a1bSJed Brown case 3: 815c4762a1bSJed Brown if (cellSimplex) { 8169566063dSJacob Faibussowitsch PetscCall(CreateSimplex_3D(comm, user, *dm)); 817c4762a1bSJed Brown } else { 8189566063dSJacob Faibussowitsch PetscCall(CreateHex_3D(comm, user->testNum, dm)); 819c4762a1bSJed Brown } 820c4762a1bSJed Brown break; 821d71ae5a4SJacob Faibussowitsch default: 822d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %" PetscInt_FMT, dim); 823c4762a1bSJed Brown } 8249566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "orig_")); 8259566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 8269566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 8279566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "material", &matLabel)); 8281baa6e33SBarry Smith if (matLabel) PetscCall(DMPlexLabelComplete(*dm, matLabel)); 8299566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 8309566063dSJacob Faibussowitsch PetscCall(DMHasLabel(*dm, "fault", &hasFault)); 831c4762a1bSJed Brown if (hasFault) { 832b253942bSMatthew G. Knepley DM dmHybrid = NULL, dmInterface = NULL; 833b253942bSMatthew G. Knepley DMLabel faultLabel, faultBdLabel, hybridLabel, splitLabel; 834c4762a1bSJed Brown 8359566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "fault", &faultLabel)); 8369566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "faultBd", &faultBdLabel)); 837caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, &splitLabel, &dmInterface, &dmHybrid)); 8389566063dSJacob Faibussowitsch PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD)); 8399566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&hybridLabel)); 8409566063dSJacob Faibussowitsch PetscCall(DMLabelView(splitLabel, PETSC_VIEWER_STDOUT_WORLD)); 8419566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&splitLabel)); 8429566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_interface_view")); 8439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmInterface)); 8449566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 845c4762a1bSJed Brown *dm = dmHybrid; 846c4762a1bSJed Brown } 8479566063dSJacob Faibussowitsch PetscCall(DMHasLabel(*dm, "fault2", &hasFault2)); 848c4762a1bSJed Brown if (hasFault2) { 849c4762a1bSJed Brown DM dmHybrid = NULL; 850c4762a1bSJed Brown DMLabel faultLabel, faultBdLabel, hybridLabel; 851c4762a1bSJed Brown 8529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "faulted_")); 8539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre")); 8549566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 8559566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 8569566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 8579566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "fault2", &faultLabel)); 8589566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "fault2Bd", &faultBdLabel)); 859caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, NULL, &dmHybrid)); 8609566063dSJacob Faibussowitsch PetscCall(DMLabelView(hybridLabel, PETSC_VIEWER_STDOUT_WORLD)); 8619566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&hybridLabel)); 8629566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 863c4762a1bSJed Brown *dm = dmHybrid; 864c4762a1bSJed Brown } 865c4762a1bSJed Brown if (user->testPartition && size > 1) { 866c4762a1bSJed Brown PetscPartitioner part; 867c4762a1bSJed Brown PetscInt *sizes = NULL; 868c4762a1bSJed Brown PetscInt *points = NULL; 869c4762a1bSJed Brown 870dd400576SPatrick Sanan if (rank == 0) { 871c4762a1bSJed Brown if (dim == 2 && cellSimplex && size == 2) { 872c4762a1bSJed Brown switch (user->testNum) { 873c4762a1bSJed Brown case 0: { 874c4762a1bSJed Brown PetscInt triSizes_p2[2] = {1, 2}; 875c4762a1bSJed Brown PetscInt triPoints_p2[3] = {0, 1, 2}; 876c4762a1bSJed Brown 8779566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 3, &points)); 8789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, triSizes_p2, 2)); 8799371c9d4SSatish Balay PetscCall(PetscArraycpy(points, triPoints_p2, 3)); 8809371c9d4SSatish Balay break; 8819371c9d4SSatish Balay } 8825b9dfbb6SMatthew G. Knepley case 3: { 8835b9dfbb6SMatthew G. Knepley PetscInt triSizes_p2[2] = {3, 3}; 8845b9dfbb6SMatthew G. Knepley PetscInt triPoints_p2[6] = {0, 1, 2, 3, 4, 5}; 8855b9dfbb6SMatthew G. Knepley 8865b9dfbb6SMatthew G. Knepley PetscCall(PetscMalloc2(2, &sizes, 6, &points)); 8875b9dfbb6SMatthew G. Knepley PetscCall(PetscArraycpy(sizes, triSizes_p2, 2)); 8889371c9d4SSatish Balay PetscCall(PetscArraycpy(points, triPoints_p2, 6)); 8899371c9d4SSatish Balay break; 8909371c9d4SSatish Balay } 891d71ae5a4SJacob Faibussowitsch default: 892d71ae5a4SJacob 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); 8935b9dfbb6SMatthew G. Knepley } 8945b9dfbb6SMatthew G. Knepley } else if (dim == 2 && cellSimplex && size == 6) { 8955b9dfbb6SMatthew G. Knepley switch (user->testNum) { 8965b9dfbb6SMatthew G. Knepley case 4: { 8975b9dfbb6SMatthew G. Knepley PetscInt triSizes_p6[6] = {1, 1, 1, 1, 1, 1}; 8985b9dfbb6SMatthew G. Knepley PetscInt triPoints_p6[6] = {0, 1, 2, 3, 4, 5}; 8995b9dfbb6SMatthew G. Knepley 9005b9dfbb6SMatthew G. Knepley PetscCall(PetscMalloc2(6, &sizes, 6, &points)); 9015b9dfbb6SMatthew G. Knepley PetscCall(PetscArraycpy(sizes, triSizes_p6, 6)); 9029371c9d4SSatish Balay PetscCall(PetscArraycpy(points, triPoints_p6, 6)); 9039371c9d4SSatish Balay break; 9049371c9d4SSatish Balay } 905d71ae5a4SJacob Faibussowitsch default: 906d71ae5a4SJacob 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); 907c4762a1bSJed Brown } 908c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && size == 2) { 909c4762a1bSJed Brown switch (user->testNum) { 910c4762a1bSJed Brown case 0: { 911c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 2}; 912c4762a1bSJed Brown PetscInt quadPoints_p2[3] = {0, 1, 2}; 913c4762a1bSJed Brown 9149566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 3, &points)); 9159566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2)); 9169371c9d4SSatish Balay PetscCall(PetscArraycpy(points, quadPoints_p2, 3)); 9179371c9d4SSatish Balay break; 9189371c9d4SSatish Balay } 919c4762a1bSJed Brown case 2: { 920c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {1, 1}; 921c4762a1bSJed Brown PetscInt quadPoints_p2[2] = {0, 1}; 922c4762a1bSJed Brown 9239566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 9249566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2)); 9259371c9d4SSatish Balay PetscCall(PetscArraycpy(points, quadPoints_p2, 2)); 9269371c9d4SSatish Balay break; 9279371c9d4SSatish Balay } 9285b9dfbb6SMatthew G. Knepley case 3: { 9295b9dfbb6SMatthew G. Knepley PetscInt quadSizes_p2[2] = {1, 1}; 9305b9dfbb6SMatthew G. Knepley PetscInt quadPoints_p2[2] = {1, 0}; 9315b9dfbb6SMatthew G. Knepley 9325b9dfbb6SMatthew G. Knepley PetscCall(PetscMalloc2(2, &sizes, 2, &points)); 9335b9dfbb6SMatthew G. Knepley PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2)); 9349371c9d4SSatish Balay PetscCall(PetscArraycpy(points, quadPoints_p2, 2)); 9359371c9d4SSatish Balay break; 9369371c9d4SSatish Balay } 937d71ae5a4SJacob Faibussowitsch default: 938d71ae5a4SJacob 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); 939c4762a1bSJed Brown } 940c4762a1bSJed Brown } else if (dim == 3 && cellSimplex && size == 2) { 941c4762a1bSJed Brown switch (user->testNum) { 942c4762a1bSJed Brown case 0: { 943c4762a1bSJed Brown PetscInt tetSizes_p2[2] = {1, 2}; 944c4762a1bSJed Brown PetscInt tetPoints_p2[3] = {0, 1, 2}; 945c4762a1bSJed Brown 9469566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 3, &points)); 9479566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2)); 9489371c9d4SSatish Balay PetscCall(PetscArraycpy(points, tetPoints_p2, 3)); 9499371c9d4SSatish Balay break; 9509371c9d4SSatish Balay } 951d71ae5a4SJacob Faibussowitsch default: 952d71ae5a4SJacob 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); 953c4762a1bSJed Brown } 954c4762a1bSJed Brown } else if (dim == 3 && !cellSimplex && size == 2) { 955c4762a1bSJed Brown switch (user->testNum) { 956c4762a1bSJed Brown case 0: { 957c4762a1bSJed Brown PetscInt hexSizes_p2[2] = {1, 2}; 958c4762a1bSJed Brown PetscInt hexPoints_p2[3] = {0, 1, 2}; 959c4762a1bSJed Brown 9609566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 3, &points)); 9619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2)); 9629371c9d4SSatish Balay PetscCall(PetscArraycpy(points, hexPoints_p2, 3)); 9639371c9d4SSatish Balay break; 9649371c9d4SSatish Balay } 965d71ae5a4SJacob Faibussowitsch default: 966d71ae5a4SJacob 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); 967c4762a1bSJed Brown } 968c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition"); 969c4762a1bSJed Brown } 9709566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 9719566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 9729566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points)); 9739566063dSJacob Faibussowitsch PetscCall(PetscFree2(sizes, points)); 974c4762a1bSJed Brown } 975c4762a1bSJed Brown { 976c4762a1bSJed Brown DM pdm = NULL; 977c4762a1bSJed Brown 978c4762a1bSJed Brown /* Distribute mesh over processes */ 9799566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm)); 980c4762a1bSJed Brown if (pdm) { 9819566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view")); 9829566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 983c4762a1bSJed Brown *dm = pdm; 984c4762a1bSJed Brown } 985c4762a1bSJed Brown } 9869566063dSJacob Faibussowitsch PetscCall(DMHasLabel(*dm, "pfault", &hasParallelFault)); 987c4762a1bSJed Brown if (hasParallelFault) { 9885b9dfbb6SMatthew G. Knepley DM dmHybrid = NULL, dmInterface; 989c4762a1bSJed Brown DMLabel faultLabel, faultBdLabel, hybridLabel; 990c4762a1bSJed Brown 9919566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "pfault", &faultLabel)); 9929566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "pfaultBd", &faultBdLabel)); 993caf9e14dSMatthew G. Knepley PetscCall(DMPlexCreateHybridMesh(*dm, faultLabel, faultBdLabel, 1, &hybridLabel, NULL, &dmInterface, &dmHybrid)); 9945b9dfbb6SMatthew G. Knepley PetscCall(DMViewFromOptions(dmInterface, NULL, "-dm_fault_view")); 9955b9dfbb6SMatthew G. Knepley { 9965b9dfbb6SMatthew G. Knepley PetscViewer viewer; 9975b9dfbb6SMatthew G. Knepley PetscMPIInt rank; 9985b9dfbb6SMatthew G. Knepley 9995b9dfbb6SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)*dm), &rank)); 10005b9dfbb6SMatthew G. Knepley PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer)); 10015b9dfbb6SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "Rank %d\n", rank)); 10025b9dfbb6SMatthew G. Knepley PetscCall(DMLabelView(hybridLabel, viewer)); 10035b9dfbb6SMatthew G. Knepley PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer)); 10045b9dfbb6SMatthew G. Knepley PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 10055b9dfbb6SMatthew G. Knepley } 10069566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&hybridLabel)); 10075b9dfbb6SMatthew G. Knepley PetscCall(DMDestroy(&dmInterface)); 10089566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 1009c4762a1bSJed Brown *dm = dmHybrid; 1010c4762a1bSJed Brown } 10119566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Hybrid Mesh")); 10129566063dSJacob Faibussowitsch PetscCall(CreateFaultLabel(*dm)); 10139566063dSJacob Faibussowitsch PetscCall(CreateDiscretization(*dm, user)); 10149566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre")); 10159566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 10169566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 10179566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 10183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1019c4762a1bSJed Brown } 1020c4762a1bSJed Brown 1021d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestMesh(DM dm, AppCtx *user) 1022d71ae5a4SJacob Faibussowitsch { 1023b253942bSMatthew G. Knepley PetscFunctionBegin; 10249566063dSJacob Faibussowitsch PetscCall(DMPlexCheckSymmetry(dm)); 10259566063dSJacob Faibussowitsch PetscCall(DMPlexCheckSkeleton(dm, 0)); 10269566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(dm, 0)); 10273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1028b253942bSMatthew G. Knepley } 1029b253942bSMatthew G. Knepley 1030d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestDiscretization(DM dm, AppCtx *user) 1031d71ae5a4SJacob Faibussowitsch { 1032b253942bSMatthew G. Knepley PetscSection s; 1033b253942bSMatthew G. Knepley 1034b253942bSMatthew G. Knepley PetscFunctionBegin; 10359566063dSJacob Faibussowitsch PetscCall(DMGetSection(dm, &s)); 10369566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)s, NULL, "-local_section_view")); 10373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1038b253942bSMatthew G. Knepley } 1039b253942bSMatthew G. Knepley 1040d71ae5a4SJacob Faibussowitsch static PetscErrorCode r(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 1041d71ae5a4SJacob Faibussowitsch { 1042b253942bSMatthew G. Knepley PetscInt d; 1043b253942bSMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = x[d]; 10443ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 1045b253942bSMatthew G. Knepley } 1046b253942bSMatthew G. Knepley 1047d71ae5a4SJacob Faibussowitsch static PetscErrorCode rp1(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 1048d71ae5a4SJacob Faibussowitsch { 1049b253942bSMatthew G. Knepley PetscInt d; 1050b253942bSMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = x[d] + (d > 0 ? 1.0 : 0.0); 10513ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 1052b253942bSMatthew G. Knepley } 1053b253942bSMatthew G. Knepley 1054d71ae5a4SJacob Faibussowitsch static PetscErrorCode phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 1055d71ae5a4SJacob Faibussowitsch { 1056b253942bSMatthew G. Knepley PetscInt d; 1057b253942bSMatthew G. Knepley u[0] = -x[1]; 1058b253942bSMatthew G. Knepley u[1] = x[0]; 1059b253942bSMatthew G. Knepley for (d = 2; d < dim; ++d) u[d] = x[d]; 10603ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 1061b253942bSMatthew G. Knepley } 1062b253942bSMatthew G. Knepley 1063b253942bSMatthew G. Knepley /* \lambda \cdot (\psi_u^- - \psi_u^+) */ 1064*81363e6fSMatthew 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[]) 1065d71ae5a4SJacob Faibussowitsch { 1066*81363e6fSMatthew G. Knepley const PetscInt Nc = dim + 1; 1067*81363e6fSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] = -u[uOff[1] + c]; 1068b253942bSMatthew G. Knepley } 1069*81363e6fSMatthew G. Knepley 1070*81363e6fSMatthew 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[]) 1071*81363e6fSMatthew G. Knepley { 1072*81363e6fSMatthew G. Knepley const PetscInt Nc = dim + 1; 1073*81363e6fSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] = u[uOff[1] + c]; 1074b253942bSMatthew G. Knepley } 1075b253942bSMatthew G. Knepley 1076b253942bSMatthew G. Knepley /* (d - u^+ + u^-) \cdot \psi_\lambda */ 1077d71ae5a4SJacob 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[]) 1078d71ae5a4SJacob Faibussowitsch { 1079b253942bSMatthew G. Knepley const PetscInt Nc = uOff[2] - uOff[1]; 1080b253942bSMatthew G. Knepley 1081*81363e6fSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] = (c > 0 ? 1.0 : 0.0) + u[c] - u[Nc + c]; 1082b253942bSMatthew G. Knepley } 1083b253942bSMatthew G. Knepley 1084b253942bSMatthew G. Knepley /* \psi_lambda \cdot (\psi_u^- - \psi_u^+) */ 1085*81363e6fSMatthew 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[]) 1086d71ae5a4SJacob Faibussowitsch { 1087*81363e6fSMatthew G. Knepley const PetscInt Nc = dim + 1; 1088*81363e6fSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) g0[c * Nc + c] = -1.0; 1089b253942bSMatthew G. Knepley } 1090*81363e6fSMatthew G. Knepley 1091*81363e6fSMatthew 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[]) 1092*81363e6fSMatthew G. Knepley { 1093*81363e6fSMatthew G. Knepley const PetscInt Nc = dim + 1; 1094*81363e6fSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0; 1095b253942bSMatthew G. Knepley } 1096b253942bSMatthew G. Knepley 1097b253942bSMatthew G. Knepley /* (-\psi_u^+ + \psi_u^-) \cdot \psi_\lambda */ 1098d71ae5a4SJacob 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[]) 1099d71ae5a4SJacob Faibussowitsch { 1100b253942bSMatthew G. Knepley const PetscInt Nc = uOff[2] - uOff[1]; 1101b253942bSMatthew G. Knepley 1102*81363e6fSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) { 1103*81363e6fSMatthew G. Knepley g0[c * Nc + c] = -1.0; 1104*81363e6fSMatthew G. Knepley g0[Nc * Nc + c * Nc + c] = 1.0; 1105b253942bSMatthew G. Knepley } 1106b253942bSMatthew G. Knepley } 1107b253942bSMatthew G. Knepley 1108d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestAssembly(DM dm, AppCtx *user) 1109d71ae5a4SJacob Faibussowitsch { 1110b253942bSMatthew G. Knepley Mat J; 1111b253942bSMatthew G. Knepley Vec locX, locF; 1112b253942bSMatthew G. Knepley PetscDS probh; 1113b253942bSMatthew G. Knepley DMLabel fault, material; 1114b253942bSMatthew G. Knepley IS cohesiveCells; 1115b7519becSMatthew G. Knepley PetscWeakForm wf; 111606ad1575SMatthew G. Knepley PetscFormKey keys[3]; 1117b253942bSMatthew G. Knepley PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx); 1118b253942bSMatthew G. Knepley PetscInt dim, Nf, cMax, cEnd, id; 1119b7519becSMatthew G. Knepley PetscMPIInt rank; 1120b253942bSMatthew G. Knepley 1121b253942bSMatthew G. Knepley PetscFunctionBegin; 11229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 11239566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 11249566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, NULL, &cMax)); 11259566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd)); 11269566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, cEnd - cMax, cMax, 1, &cohesiveCells)); 11279566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "cohesive", &fault)); 11289566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX)); 11299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)locX, "Local Solution")); 11309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locF)); 11319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)locF, "Local Residual")); 11329566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J)); 11339566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian")); 1134b253942bSMatthew G. Knepley 1135b253942bSMatthew G. Knepley /* The initial guess has displacement shifted by one unit in each fault parallel direction across the fault */ 11369566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "material", &material)); 1137b253942bSMatthew G. Knepley id = 1; 1138b253942bSMatthew G. Knepley initialGuess[0] = r; 1139b253942bSMatthew G. Knepley initialGuess[1] = NULL; 11409566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX)); 1141b253942bSMatthew G. Knepley id = 2; 1142b253942bSMatthew G. Knepley initialGuess[0] = rp1; 1143b253942bSMatthew G. Knepley initialGuess[1] = NULL; 11449566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, material, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX)); 1145b253942bSMatthew G. Knepley id = 1; 1146b253942bSMatthew G. Knepley initialGuess[0] = NULL; 1147b253942bSMatthew G. Knepley initialGuess[1] = phi; 11489566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, fault, 1, &id, PETSC_DETERMINE, NULL, initialGuess, NULL, INSERT_VALUES, locX)); 11499566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(locX, NULL, "-local_solution_view")); 1150b253942bSMatthew G. Knepley 1151e5939c1dSMatthew G. Knepley PetscCall(DMGetCellDS(dm, cMax, &probh, NULL)); 11529566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(probh, &wf)); 11539566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(probh, &Nf)); 1154*81363e6fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 1, 0, 0, 0, f0_bd_u_neg, 0, NULL)); 1155*81363e6fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, material, 2, 0, 0, 0, f0_bd_u_pos, 0, NULL)); 1156*81363e6fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 1, 0, 1, 0, 0, g0_bd_ul_neg, 0, NULL, 0, NULL, 0, NULL)); 1157*81363e6fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdJacobian(wf, material, 2, 0, 1, 0, 0, g0_bd_ul_pos, 0, NULL, 0, NULL, 0, NULL)); 1158b7519becSMatthew G. Knepley if (Nf > 1) { 11599566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, fault, 1, 1, 0, 0, f0_bd_l, 0, NULL)); 11609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdJacobian(wf, fault, 1, 1, 0, 0, 0, g0_bd_lu, 0, NULL, 0, NULL, 0, NULL)); 1161b7519becSMatthew G. Knepley } 1162c5853193SPierre Jolivet if (rank == 0) PetscCall(PetscDSView(probh, NULL)); 1163b253942bSMatthew G. Knepley 1164*81363e6fSMatthew G. Knepley keys[0].label = material; 1165*81363e6fSMatthew G. Knepley keys[0].value = 1; 11666528b96dSMatthew G. Knepley keys[0].field = 0; 1167dd96b1f9SToby Isaac keys[0].part = 0; 1168b7519becSMatthew G. Knepley keys[1].label = material; 1169b7519becSMatthew G. Knepley keys[1].value = 2; 11706528b96dSMatthew G. Knepley keys[1].field = 0; 1171dd96b1f9SToby Isaac keys[1].part = 0; 1172b7519becSMatthew G. Knepley keys[2].label = fault; 1173b7519becSMatthew G. Knepley keys[2].value = 1; 1174b7519becSMatthew G. Knepley keys[2].field = 1; 1175dd96b1f9SToby Isaac keys[2].part = 0; 11769566063dSJacob Faibussowitsch PetscCall(VecSet(locF, 0.)); 11779566063dSJacob Faibussowitsch PetscCall(DMPlexComputeResidual_Hybrid_Internal(dm, keys, cohesiveCells, 0.0, locX, NULL, 0.0, locF, user)); 11789566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(locF, NULL, "-local_residual_view")); 11799566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J)); 11809566063dSJacob Faibussowitsch PetscCall(DMPlexComputeJacobian_Hybrid_Internal(dm, keys, cohesiveCells, 0.0, 0.0, locX, NULL, J, J, user)); 118149664dceSMatthew G. Knepley PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 118249664dceSMatthew G. Knepley PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 11839566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(J, NULL, "-local_jacobian_view")); 1184b253942bSMatthew G. Knepley 11859566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX)); 11869566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locF)); 11879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 11889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cohesiveCells)); 1189e5939c1dSMatthew G. Knepley 1190e5939c1dSMatthew G. Knepley if (cMax < cEnd) { 1191e5939c1dSMatthew G. Knepley PetscDS ds; 1192e5939c1dSMatthew G. Knepley PetscFE fe; 1193e5939c1dSMatthew G. Knepley PetscQuadrature quad; 1194e5939c1dSMatthew G. Knepley IS *perm; 1195e5939c1dSMatthew G. Knepley const PetscInt *cone; 1196e5939c1dSMatthew G. Knepley PetscInt Na, a; 1197e5939c1dSMatthew G. Knepley 1198e5939c1dSMatthew G. Knepley PetscCall(DMPlexGetCone(dm, cMax, &cone)); 1199e5939c1dSMatthew G. Knepley PetscCall(DMGetCellDS(dm, cMax, &ds, NULL)); 1200e5939c1dSMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe)); 1201e5939c1dSMatthew G. Knepley PetscCall(PetscFEGetQuadrature(fe, &quad)); 1202e5939c1dSMatthew G. Knepley PetscCall(PetscQuadratureComputePermutations(quad, &Na, &perm)); 1203e5939c1dSMatthew G. Knepley for (a = 0; a < Na; ++a) PetscCall(ISDestroy(&perm[a])); 1204e5939c1dSMatthew G. Knepley PetscCall(PetscFree(perm)); 1205e5939c1dSMatthew G. Knepley } 12063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1207b253942bSMatthew G. Knepley } 1208b253942bSMatthew G. Knepley 1209d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 1210d71ae5a4SJacob Faibussowitsch { 1211c4762a1bSJed Brown DM dm; 1212c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 1213c4762a1bSJed Brown 1214327415f7SBarry Smith PetscFunctionBeginUser; 12159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 12169566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 12179566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 12189566063dSJacob Faibussowitsch PetscCall(TestMesh(dm, &user)); 12199566063dSJacob Faibussowitsch PetscCall(TestDiscretization(dm, &user)); 12209566063dSJacob Faibussowitsch PetscCall(TestAssembly(dm, &user)); 12219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 12229566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1223b122ec5aSJacob Faibussowitsch return 0; 1224c4762a1bSJed Brown } 1225c4762a1bSJed Brown 1226c4762a1bSJed Brown /*TEST 1227c4762a1bSJed Brown testset: 1228b253942bSMatthew G. Knepley args: -orig_dm_plex_check_all -dm_plex_check_all \ 1229b7519becSMatthew G. Knepley -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1 -local_section_view \ 1230b253942bSMatthew G. Knepley -local_solution_view -local_residual_view -local_jacobian_view 1231c4762a1bSJed Brown test: 1232c4762a1bSJed Brown suffix: tri_0 1233c4762a1bSJed Brown args: -dim 2 1234*81363e6fSMatthew 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" 1235c4762a1bSJed Brown test: 1236c4762a1bSJed Brown suffix: tri_t1_0 1237c4762a1bSJed Brown args: -dim 2 -test_num 1 1238*81363e6fSMatthew 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" 1239c4762a1bSJed Brown test: 12405b9dfbb6SMatthew G. Knepley suffix: tri_t2_0 12415b9dfbb6SMatthew G. Knepley args: -dim 2 -test_num 2 1242*81363e6fSMatthew 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" 1243*81363e6fSMatthew G. Knepley 12445b9dfbb6SMatthew G. Knepley test: 12453bdf08b3SMatthew G. Knepley suffix: tri_t5_0 12463bdf08b3SMatthew G. Knepley args: -dim 2 -test_num 5 1247*81363e6fSMatthew 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" 1248*81363e6fSMatthew G. Knepley 12493bdf08b3SMatthew G. Knepley test: 1250c4762a1bSJed Brown suffix: tet_0 1251c4762a1bSJed Brown args: -dim 3 1252*81363e6fSMatthew 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" 1253*81363e6fSMatthew G. Knepley 1254c4762a1bSJed Brown test: 1255b253942bSMatthew G. Knepley suffix: tet_t1_0 1256b253942bSMatthew G. Knepley args: -dim 3 -test_num 1 1257*81363e6fSMatthew 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" 1258b253942bSMatthew G. Knepley 1259b253942bSMatthew G. Knepley testset: 1260b253942bSMatthew G. Knepley args: -orig_dm_plex_check_all -dm_plex_check_all \ 1261b7519becSMatthew G. Knepley -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1 1262b253942bSMatthew G. Knepley test: 1263c4762a1bSJed Brown suffix: tet_1 1264c4762a1bSJed Brown nsize: 2 1265c4762a1bSJed Brown args: -dim 3 1266*81363e6fSMatthew 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" 1267*81363e6fSMatthew G. Knepley 1268c4762a1bSJed Brown test: 1269b253942bSMatthew G. Knepley suffix: tri_1 1270b253942bSMatthew G. Knepley nsize: 2 1271b253942bSMatthew G. Knepley args: -dim 2 1272*81363e6fSMatthew 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" 1273*81363e6fSMatthew G. Knepley 12745b9dfbb6SMatthew G. Knepley test: 12755b9dfbb6SMatthew G. Knepley suffix: tri_t3_0 12765b9dfbb6SMatthew G. Knepley nsize: 2 12775b9dfbb6SMatthew G. Knepley args: -dim 2 -test_num 3 1278*81363e6fSMatthew 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" 1279*81363e6fSMatthew G. Knepley 12805b9dfbb6SMatthew G. Knepley test: 12815b9dfbb6SMatthew G. Knepley suffix: tri_t4_0 12825b9dfbb6SMatthew G. Knepley nsize: 6 12835b9dfbb6SMatthew G. Knepley args: -dim 2 -test_num 4 1284*81363e6fSMatthew 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" 1285c4762a1bSJed Brown 1286c4762a1bSJed Brown testset: 1287ecfb78b5SMatthew G. Knepley args: -orig_dm_plex_check_all -dm_plex_check_all \ 1288b7519becSMatthew G. Knepley -displacement_petscspace_degree 1 -faulttraction_petscspace_degree 1 1289c4762a1bSJed Brown # 2D Quads 1290c4762a1bSJed Brown test: 1291c4762a1bSJed Brown suffix: quad_0 1292c4762a1bSJed Brown args: -dim 2 -cell_simplex 0 1293*81363e6fSMatthew 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" 1294*81363e6fSMatthew G. Knepley 1295c4762a1bSJed Brown test: 1296c4762a1bSJed Brown suffix: quad_1 1297c4762a1bSJed Brown nsize: 2 1298c4762a1bSJed Brown args: -dim 2 -cell_simplex 0 1299*81363e6fSMatthew 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" 1300*81363e6fSMatthew G. Knepley 1301c4762a1bSJed Brown test: 1302c4762a1bSJed Brown suffix: quad_t1_0 130354fcfd0cSMatthew G. Knepley args: -dim 2 -cell_simplex 0 -test_num 1 -faulted_dm_plex_check_all 1304*81363e6fSMatthew 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" 1305*81363e6fSMatthew G. Knepley 13065b9dfbb6SMatthew G. Knepley test: 13075b9dfbb6SMatthew G. Knepley suffix: quad_t2_0 13085b9dfbb6SMatthew G. Knepley nsize: 2 13095b9dfbb6SMatthew G. Knepley args: -dim 2 -cell_simplex 0 -test_num 2 1310*81363e6fSMatthew 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" 1311*81363e6fSMatthew G. Knepley 13125b9dfbb6SMatthew G. Knepley test: 13135b9dfbb6SMatthew G. Knepley # TODO: The PetscSF is wrong here (connects to wrong side of split) 13145b9dfbb6SMatthew G. Knepley suffix: quad_t3_0 13155b9dfbb6SMatthew G. Knepley nsize: 2 13165b9dfbb6SMatthew G. Knepley args: -dim 2 -cell_simplex 0 -test_num 3 1317*81363e6fSMatthew 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" 1318*81363e6fSMatthew G. Knepley 1319c4762a1bSJed Brown # 3D Hex 1320c4762a1bSJed Brown test: 1321c4762a1bSJed Brown suffix: hex_0 1322c4762a1bSJed Brown args: -dim 3 -cell_simplex 0 1323*81363e6fSMatthew 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" 1324c4762a1bSJed Brown test: 1325c4762a1bSJed Brown suffix: hex_1 1326c4762a1bSJed Brown nsize: 2 1327c4762a1bSJed Brown args: -dim 3 -cell_simplex 0 1328*81363e6fSMatthew 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" 1329c4762a1bSJed Brown test: 1330c4762a1bSJed Brown suffix: hex_t1_0 1331c4762a1bSJed Brown args: -dim 3 -cell_simplex 0 -test_num 1 1332*81363e6fSMatthew 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" 1333c4762a1bSJed Brown test: 1334c4762a1bSJed Brown suffix: hex_t2_0 1335c4762a1bSJed Brown args: -dim 3 -cell_simplex 0 -test_num 2 1336*81363e6fSMatthew 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" 1337c4762a1bSJed Brown 1338c4762a1bSJed Brown TEST*/ 1339