1231de6a2SMatthew G. Knepley static char help[] = "The main goal of this code is to retrieve the original element numbers as found in the " 2231de6a2SMatthew G. Knepley "initial partitions (sInitialPartition)... but after the call to DMPlexDistribute"; 3231de6a2SMatthew G. Knepley 4231de6a2SMatthew G. Knepley #include <petsc.h> 5231de6a2SMatthew G. Knepley 69f046c2dSEric Chamberland /* Coordinates of a 2x5 rectangular mesh of quads : */ 7231de6a2SMatthew G. Knepley PetscReal sCoords2x5Mesh[18][2] = { 8231de6a2SMatthew G. Knepley {0.00000000000000000e+00, 0.00000000000000000e+00}, 9231de6a2SMatthew G. Knepley {2.00000000000000000e+00, 0.00000000000000000e+00}, 10231de6a2SMatthew G. Knepley {0.00000000000000000e+00, 1.00000000000000000e+00}, 11231de6a2SMatthew G. Knepley {2.00000000000000000e+00, 1.00000000000000000e+00}, 12231de6a2SMatthew G. Knepley {9.99999999997387978e-01, 0.00000000000000000e+00}, 13231de6a2SMatthew G. Knepley {9.99999999997387978e-01, 1.00000000000000000e+00}, 14231de6a2SMatthew G. Knepley {0.00000000000000000e+00, 2.00000000000000011e-01}, 15231de6a2SMatthew G. Knepley {0.00000000000000000e+00, 4.00000000000000022e-01}, 16231de6a2SMatthew G. Knepley {0.00000000000000000e+00, 5.99999999999999978e-01}, 17231de6a2SMatthew G. Knepley {0.00000000000000000e+00, 8.00000000000000044e-01}, 18231de6a2SMatthew G. Knepley {2.00000000000000000e+00, 2.00000000000000011e-01}, 19231de6a2SMatthew G. Knepley {2.00000000000000000e+00, 4.00000000000000022e-01}, 20231de6a2SMatthew G. Knepley {2.00000000000000000e+00, 5.99999999999999978e-01}, 21231de6a2SMatthew G. Knepley {2.00000000000000000e+00, 8.00000000000000044e-01}, 22231de6a2SMatthew G. Knepley {9.99999999997387756e-01, 2.00000000000000011e-01}, 23231de6a2SMatthew G. Knepley {9.99999999997387978e-01, 4.00000000000000022e-01}, 24231de6a2SMatthew G. Knepley {9.99999999997387978e-01, 6.00000000000000089e-01}, 259371c9d4SSatish Balay {9.99999999997388089e-01, 8.00000000000000044e-01} 269371c9d4SSatish Balay }; 27231de6a2SMatthew G. Knepley 289f046c2dSEric Chamberland /* Connectivity of a 2x5 rectangular mesh of quads : */ 29231de6a2SMatthew G. Knepley const PetscInt sConnectivity2x5Mesh[10][4] = { 30231de6a2SMatthew G. Knepley {0, 4, 14, 6 }, 31231de6a2SMatthew G. Knepley {6, 14, 15, 7 }, 32231de6a2SMatthew G. Knepley {7, 15, 16, 8 }, 33231de6a2SMatthew G. Knepley {8, 16, 17, 9 }, 34231de6a2SMatthew G. Knepley {9, 17, 5, 2 }, 35231de6a2SMatthew G. Knepley {4, 1, 10, 14}, 36231de6a2SMatthew G. Knepley {14, 10, 11, 15}, 37231de6a2SMatthew G. Knepley {15, 11, 12, 16}, 38231de6a2SMatthew G. Knepley {16, 12, 13, 17}, 399371c9d4SSatish Balay {17, 13, 3, 5 } 409371c9d4SSatish Balay }; 41231de6a2SMatthew G. Knepley 429f046c2dSEric Chamberland /* Partitions of a 2x5 rectangular mesh of quads : */ 43231de6a2SMatthew G. Knepley const PetscInt sInitialPartition2x5Mesh[2][5] = { 44231de6a2SMatthew G. Knepley {0, 2, 4, 6, 8}, 45231de6a2SMatthew G. Knepley {1, 3, 5, 7, 9} 46231de6a2SMatthew G. Knepley }; 47231de6a2SMatthew G. Knepley 48231de6a2SMatthew G. Knepley const PetscInt sNLoclCells2x5Mesh = 5; 49231de6a2SMatthew G. Knepley const PetscInt sNGlobVerts2x5Mesh = 18; 50231de6a2SMatthew G. Knepley 5172f9f102SEric Chamberland /* Mixed mesh : quads and triangles (4 first quads above divided into triangles*/ 5272f9f102SEric Chamberland /* Connectivity of a 2x5 rectangular mesh of quads : */ 5372f9f102SEric Chamberland const PetscInt sConnectivityMixedTQMesh[14][4] = { 5472f9f102SEric Chamberland {0, 4, 6, -1}, 5572f9f102SEric Chamberland {4, 14, 6, -1}, 5672f9f102SEric Chamberland {6, 14, 7, -1}, 5772f9f102SEric Chamberland {14, 15, 7, -1}, 5872f9f102SEric Chamberland {7, 15, 8, -1}, 5972f9f102SEric Chamberland {15, 16, 8, -1}, 6072f9f102SEric Chamberland {8, 16, 9, -1}, 6172f9f102SEric Chamberland {16, 17, 9, -1}, 6272f9f102SEric Chamberland {9, 17, 5, 2 }, 6372f9f102SEric Chamberland {4, 1, 10, 14}, 6472f9f102SEric Chamberland {14, 10, 11, 15}, 6572f9f102SEric Chamberland {15, 11, 12, 16}, 6672f9f102SEric Chamberland {16, 12, 13, 17}, 6772f9f102SEric Chamberland {17, 13, 3, 5 } 6872f9f102SEric Chamberland }; 6972f9f102SEric Chamberland 7072f9f102SEric Chamberland /* Partitions for the rectangular mesh of quads and triangles: */ 7172f9f102SEric Chamberland const PetscInt sInitialPartitionMixedTQMesh[2][7] = { 7272f9f102SEric Chamberland {0, 1, 4, 5, 8, 10, 12}, 7372f9f102SEric Chamberland {2, 3, 6, 7, 9, 11, 13} 7472f9f102SEric Chamberland }; 7572f9f102SEric Chamberland 7672f9f102SEric Chamberland const PetscInt sNLoclCellsMixedTQMesh = 7; 7772f9f102SEric Chamberland const PetscInt sNGlobVertsMixedTQMesh = 18; 7872f9f102SEric Chamberland 799f046c2dSEric Chamberland /* Prisms mesh */ 809f046c2dSEric Chamberland PetscReal sCoordsPrismsMesh[125][3] = { 819f046c2dSEric Chamberland {2.24250931694056355e-01, 0.00000000000000000e+00, 0.00000000000000000e+00}, 829f046c2dSEric Chamberland {2.20660660151932697e-01, 2.87419338850266937e-01, 0.00000000000000000e+00}, 839f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 2.70243537720639027e-01}, 849f046c2dSEric Chamberland {2.32445727460992402e-01, 0.00000000000000000e+00, 2.60591845015572310e-01}, 859f046c2dSEric Chamberland {2.41619971105419079e-01, 2.69894910706158231e-01, 2.42844781736072490e-01}, 869f046c2dSEric Chamberland {0.00000000000000000e+00, 2.46523339883120779e-01, 2.69072907562752262e-01}, 879f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 0.00000000000000000e+00}, 889f046c2dSEric Chamberland {1.00000000000000000e+00, 2.75433417601945563e-01, 0.00000000000000000e+00}, 899f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 2.33748605950385602e-01}, 909f046c2dSEric Chamberland {7.32445727460992457e-01, 0.00000000000000000e+00, 2.42344379130445597e-01}, 919f046c2dSEric Chamberland {1.00000000000000000e+00, 2.78258478013028610e-01, 2.57379172987105553e-01}, 929f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 0.00000000000000000e+00}, 939f046c2dSEric Chamberland {7.49586880891153995e-01, 1.00000000000000000e+00, 0.00000000000000000e+00}, 949f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 2.51949651266657582e-01}, 959f046c2dSEric Chamberland {7.41619971105419107e-01, 7.69894910706158120e-01, 2.33697838509081768e-01}, 969f046c2dSEric Chamberland {1.00000000000000000e+00, 7.78258478013028610e-01, 2.66479695645241543e-01}, 979f046c2dSEric Chamberland {7.55042653233710115e-01, 1.00000000000000000e+00, 2.58019637386860512e-01}, 989f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 0.00000000000000000e+00}, 999f046c2dSEric Chamberland {0.00000000000000000e+00, 7.59235710423095789e-01, 0.00000000000000000e+00}, 1009f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 2.17232187874490473e-01}, 1019f046c2dSEric Chamberland {0.00000000000000000e+00, 7.46523339883120807e-01, 2.42567232639677999e-01}, 1029f046c2dSEric Chamberland {2.55042653233710115e-01, 1.00000000000000000e+00, 2.40660905690776916e-01}, 1039f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 0.00000000000000000e+00}, 1049f046c2dSEric Chamberland {2.38934376044866809e-01, 0.00000000000000000e+00, 1.00000000000000000e+00}, 1059f046c2dSEric Chamberland {2.18954188589218168e-01, 2.26916038449581692e-01, 1.00000000000000000e+00}, 1069f046c2dSEric Chamberland {2.39787449636397643e-01, 0.00000000000000000e+00, 7.60591845015572310e-01}, 1079f046c2dSEric Chamberland {2.40766735324061815e-01, 2.39643260505815608e-01, 7.42844781736072490e-01}, 1089f046c2dSEric Chamberland {0.00000000000000000e+00, 2.57448248192627016e-01, 7.69072907562752262e-01}, 1099f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 1.00000000000000000e+00}, 1109f046c2dSEric Chamberland {1.00000000000000000e+00, 2.38666970143638080e-01, 1.00000000000000000e+00}, 1119f046c2dSEric Chamberland {7.39787449636397643e-01, 0.00000000000000000e+00, 7.42344379130445597e-01}, 1129f046c2dSEric Chamberland {1.00000000000000000e+00, 2.59875254283874868e-01, 7.57379172987105553e-01}, 1139f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 1.00000000000000000e+00}, 1149f046c2dSEric Chamberland {7.76318984007844159e-01, 1.00000000000000000e+00, 1.00000000000000000e+00}, 1159f046c2dSEric Chamberland {7.40766735324061787e-01, 7.39643260505815636e-01, 7.33697838509081768e-01}, 1169f046c2dSEric Chamberland {1.00000000000000000e+00, 7.59875254283874924e-01, 7.66479695645241543e-01}, 1179f046c2dSEric Chamberland {7.68408704792055142e-01, 1.00000000000000000e+00, 7.58019637386860512e-01}, 1189f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 1.00000000000000000e+00}, 1199f046c2dSEric Chamberland {0.00000000000000000e+00, 7.81085527042108207e-01, 1.00000000000000000e+00}, 1209f046c2dSEric Chamberland {0.00000000000000000e+00, 7.57448248192627016e-01, 7.42567232639678054e-01}, 1219f046c2dSEric Chamberland {2.68408704792055197e-01, 1.00000000000000000e+00, 7.40660905690776916e-01}, 1229f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 1.00000000000000000e+00}, 1239f046c2dSEric Chamberland {7.24250931694056410e-01, 0.00000000000000000e+00, 0.00000000000000000e+00}, 1249f046c2dSEric Chamberland {7.24250931694056410e-01, 2.75433417601945563e-01, 0.00000000000000000e+00}, 1259f046c2dSEric Chamberland {4.44911591845989052e-01, 2.87419338850266937e-01, 0.00000000000000000e+00}, 1269f046c2dSEric Chamberland {4.64891454921984804e-01, 0.00000000000000000e+00, 2.50940152310505593e-01}, 1279f046c2dSEric Chamberland {4.74065698566411453e-01, 2.69894910706158231e-01, 2.33193089031005774e-01}, 1289f046c2dSEric Chamberland {4.48501863388112709e-01, 0.00000000000000000e+00, 0.00000000000000000e+00}, 1299f046c2dSEric Chamberland {7.20660660151932753e-01, 7.87419338850266937e-01, 0.00000000000000000e+00}, 1309f046c2dSEric Chamberland {7.20660660151932753e-01, 5.62852756452212555e-01, 0.00000000000000000e+00}, 1319f046c2dSEric Chamberland {2.20660660151932697e-01, 5.46655049273362614e-01, 0.00000000000000000e+00}, 1329f046c2dSEric Chamberland {4.83239942210838158e-01, 5.39789821412316462e-01, 2.15446025751505982e-01}, 1339f046c2dSEric Chamberland {7.41619971105419107e-01, 5.48153388719186951e-01, 2.48227882887665785e-01}, 1349f046c2dSEric Chamberland {2.41619971105419079e-01, 5.16418250589278927e-01, 2.41674151578185781e-01}, 1359f046c2dSEric Chamberland {4.41321320303865394e-01, 5.74838677700533873e-01, 0.00000000000000000e+00}, 1369f046c2dSEric Chamberland {1.00000000000000000e+00, 7.75433417601945507e-01, 0.00000000000000000e+00}, 1379f046c2dSEric Chamberland {1.00000000000000000e+00, 5.56516956026057219e-01, 2.81009740023825560e-01}, 1389f046c2dSEric Chamberland {7.32445727460992457e-01, 2.78258478013028610e-01, 2.65974946167165549e-01}, 1399f046c2dSEric Chamberland {1.00000000000000000e+00, 5.50866835203891125e-01, 0.00000000000000000e+00}, 1409f046c2dSEric Chamberland {0.00000000000000000e+00, 2.59235710423095733e-01, 0.00000000000000000e+00}, 1419f046c2dSEric Chamberland {0.00000000000000000e+00, 4.93046679766241558e-01, 2.67902277404865552e-01}, 1429f046c2dSEric Chamberland {2.55042653233710115e-01, 7.46523339883120807e-01, 2.65995950455964469e-01}, 1439f046c2dSEric Chamberland {0.00000000000000000e+00, 5.18471420846191466e-01, 0.00000000000000000e+00}, 1449f046c2dSEric Chamberland {2.49586880891154023e-01, 1.00000000000000000e+00, 0.00000000000000000e+00}, 1459f046c2dSEric Chamberland {2.49586880891154023e-01, 7.59235710423095789e-01, 0.00000000000000000e+00}, 1469f046c2dSEric Chamberland {4.70247541043086748e-01, 7.87419338850266937e-01, 0.00000000000000000e+00}, 1479f046c2dSEric Chamberland {5.10085306467420230e-01, 1.00000000000000000e+00, 2.64089623507063387e-01}, 1489f046c2dSEric Chamberland {4.96662624339129222e-01, 7.69894910706158231e-01, 2.39767824629284698e-01}, 1499f046c2dSEric Chamberland {4.99173761782308045e-01, 1.00000000000000000e+00, 0.00000000000000000e+00}, 1509f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 7.70243537720639027e-01}, 1519f046c2dSEric Chamberland {2.40640523227928449e-01, 0.00000000000000000e+00, 5.21183690031144620e-01}, 1529f046c2dSEric Chamberland {2.62579282058905461e-01, 2.52370482562049525e-01, 4.85689563472144981e-01}, 1539f046c2dSEric Chamberland {0.00000000000000000e+00, 2.33810969343145825e-01, 5.38145815125504523e-01}, 1549f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 5.40487075441278053e-01}, 1559f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 7.33748605950385602e-01}, 1569f046c2dSEric Chamberland {7.40640523227928504e-01, 0.00000000000000000e+00, 4.84688758260891195e-01}, 1579f046c2dSEric Chamberland {1.00000000000000000e+00, 2.81083538424111656e-01, 5.14758345974211107e-01}, 1589f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 4.67497211900771203e-01}, 1599f046c2dSEric Chamberland {7.38934376044866781e-01, 0.00000000000000000e+00, 1.00000000000000000e+00}, 1609f046c2dSEric Chamberland {4.79574899272795285e-01, 0.00000000000000000e+00, 7.50940152310505593e-01}, 1619f046c2dSEric Chamberland {4.77868752089733617e-01, 0.00000000000000000e+00, 1.00000000000000000e+00}, 1629f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 7.51949651266657582e-01}, 1639f046c2dSEric Chamberland {7.62579282058905461e-01, 7.52370482562049525e-01, 4.67395677018163536e-01}, 1649f046c2dSEric Chamberland {1.00000000000000000e+00, 7.81083538424111712e-01, 5.32959391290483087e-01}, 1659f046c2dSEric Chamberland {7.60498425576266124e-01, 1.00000000000000000e+00, 5.16039274773721024e-01}, 1669f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 5.03899302533315163e-01}, 1679f046c2dSEric Chamberland {7.18954188589218113e-01, 7.26916038449581636e-01, 1.00000000000000000e+00}, 1689f046c2dSEric Chamberland {4.81533470648123629e-01, 4.79286521011631217e-01, 7.15446025751505954e-01}, 1699f046c2dSEric Chamberland {4.57888564634085005e-01, 2.26916038449581692e-01, 1.00000000000000000e+00}, 1709f046c2dSEric Chamberland {4.95273172597062383e-01, 7.26916038449581636e-01, 1.00000000000000000e+00}, 1719f046c2dSEric Chamberland {4.37908377178436337e-01, 4.53832076899163384e-01, 1.00000000000000000e+00}, 1729f046c2dSEric Chamberland {1.00000000000000000e+00, 7.38666970143638135e-01, 1.00000000000000000e+00}, 1739f046c2dSEric Chamberland {1.00000000000000000e+00, 5.19750508567749736e-01, 7.81009740023825616e-01}, 1749f046c2dSEric Chamberland {7.38934376044866781e-01, 2.38666970143638080e-01, 1.00000000000000000e+00}, 1759f046c2dSEric Chamberland {7.18954188589218113e-01, 4.65583008593219771e-01, 1.00000000000000000e+00}, 1769f046c2dSEric Chamberland {1.00000000000000000e+00, 4.77333940287276159e-01, 1.00000000000000000e+00}, 1779f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 7.17232187874490501e-01}, 1789f046c2dSEric Chamberland {0.00000000000000000e+00, 7.33810969343145825e-01, 4.85134465279355998e-01}, 1799f046c2dSEric Chamberland {2.60498425576266179e-01, 1.00000000000000000e+00, 4.81321811381553832e-01}, 1809f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 4.34464375748980947e-01}, 1819f046c2dSEric Chamberland {0.00000000000000000e+00, 2.81085527042108152e-01, 1.00000000000000000e+00}, 1829f046c2dSEric Chamberland {0.00000000000000000e+00, 5.14896496385254032e-01, 7.67902277404865607e-01}, 1839f046c2dSEric Chamberland {2.76318984007844215e-01, 7.81085527042108207e-01, 1.00000000000000000e+00}, 1849f046c2dSEric Chamberland {2.18954188589218168e-01, 5.08001565491689844e-01, 1.00000000000000000e+00}, 1859f046c2dSEric Chamberland {0.00000000000000000e+00, 5.62171054084216304e-01, 1.00000000000000000e+00}, 1869f046c2dSEric Chamberland {2.76318984007844215e-01, 1.00000000000000000e+00, 1.00000000000000000e+00}, 1879f046c2dSEric Chamberland {5.36817409584110394e-01, 1.00000000000000000e+00, 7.64089623507063331e-01}, 1889f046c2dSEric Chamberland {5.52637968015688430e-01, 1.00000000000000000e+00, 1.00000000000000000e+00}, 1899f046c2dSEric Chamberland {5.03219805286833965e-01, 2.52370482562049525e-01, 4.66386178062011547e-01}, 1909f046c2dSEric Chamberland {4.80554184960459430e-01, 2.39643260505815608e-01, 7.33193089031005774e-01}, 1919f046c2dSEric Chamberland {4.81281046455856898e-01, 0.00000000000000000e+00, 5.01880304621011186e-01}, 1929f046c2dSEric Chamberland {7.62579282058905461e-01, 5.33454020986161126e-01, 4.96455765775331570e-01}, 1939f046c2dSEric Chamberland {2.62579282058905461e-01, 4.86181451905195350e-01, 4.83348303156371562e-01}, 1949f046c2dSEric Chamberland {7.40766735324061787e-01, 4.99518514789690449e-01, 7.48227882887665841e-01}, 1959f046c2dSEric Chamberland {2.40766735324061815e-01, 4.97091508698442541e-01, 7.41674151578185725e-01}, 1969f046c2dSEric Chamberland {5.25158564117810922e-01, 5.04740965124099050e-01, 4.30892051503011964e-01}, 1979f046c2dSEric Chamberland {7.40640523227928504e-01, 2.81083538424111656e-01, 5.31949892334331098e-01}, 1989f046c2dSEric Chamberland {7.39787449636397643e-01, 2.59875254283874868e-01, 7.65974946167165549e-01}, 1999f046c2dSEric Chamberland {1.00000000000000000e+00, 5.62167076848223313e-01, 5.62019480047651121e-01}, 2009f046c2dSEric Chamberland {2.60498425576266179e-01, 7.33810969343145825e-01, 5.31991900911928939e-01}, 2019f046c2dSEric Chamberland {2.68408704792055197e-01, 7.57448248192627016e-01, 7.65995950455964469e-01}, 2029f046c2dSEric Chamberland {0.00000000000000000e+00, 4.67621938686291649e-01, 5.35804554809731104e-01}, 2039f046c2dSEric Chamberland {5.23077707635171585e-01, 7.52370482562049525e-01, 4.79535649258569396e-01}, 2049f046c2dSEric Chamberland {5.09175440116116929e-01, 7.39643260505815636e-01, 7.39767824629284698e-01}, 2059f046c2dSEric Chamberland {5.20996851152532359e-01, 1.00000000000000000e+00, 5.28179247014126774e-01} 2069f046c2dSEric Chamberland }; 2079f046c2dSEric Chamberland 2089f046c2dSEric Chamberland const PetscInt sConnectivityPrismsMesh[128][6] = { 2099f046c2dSEric Chamberland /* rank 0 */ 2109f046c2dSEric Chamberland {11, 7, 42, 8, 10, 9 }, 2119f046c2dSEric Chamberland {47, 42, 43, 45, 9, 57 }, 2129f046c2dSEric Chamberland {8, 10, 9, 77, 76, 75 }, 2139f046c2dSEric Chamberland {45, 9, 57, 110, 75, 116}, 2149f046c2dSEric Chamberland {17, 48, 55, 13, 14, 15 }, 2159f046c2dSEric Chamberland {58, 55, 49, 56, 15, 52 }, 2169f046c2dSEric Chamberland {13, 14, 15, 85, 82, 83 }, 2179f046c2dSEric Chamberland {56, 15, 52, 118, 83, 111}, 2189f046c2dSEric Chamberland {6, 0, 1, 2, 3, 4 }, 2199f046c2dSEric Chamberland {54, 1, 44, 51, 4, 46 }, 2209f046c2dSEric Chamberland {2, 3, 4, 73, 70, 71 }, 2219f046c2dSEric Chamberland {51, 4, 46, 115, 71, 108}, 2229f046c2dSEric Chamberland {58, 49, 43, 56, 52, 57 }, 2239f046c2dSEric Chamberland {47, 43, 44, 45, 57, 46 }, 2249f046c2dSEric Chamberland {56, 52, 57, 118, 111, 116}, 2259f046c2dSEric Chamberland {45, 57, 46, 110, 116, 108}, 2269f046c2dSEric Chamberland {77, 76, 75, 74, 31, 30 }, 2279f046c2dSEric Chamberland {110, 75, 116, 79, 30, 117}, 2289f046c2dSEric Chamberland {74, 31, 30, 32, 29, 78 }, 2299f046c2dSEric Chamberland {79, 30, 117, 80, 78, 93 }, 2309f046c2dSEric Chamberland {85, 82, 83, 81, 34, 35 }, 2319f046c2dSEric Chamberland {118, 83, 111, 92, 35, 113}, 2329f046c2dSEric Chamberland {81, 34, 35, 37, 86, 91 }, 2339f046c2dSEric Chamberland {92, 35, 113, 95, 91, 94 }, 2349f046c2dSEric Chamberland {73, 70, 71, 69, 25, 26 }, 2359f046c2dSEric Chamberland {115, 71, 108, 87, 26, 109}, 2369f046c2dSEric Chamberland {69, 25, 26, 28, 23, 24 }, 2379f046c2dSEric Chamberland {87, 26, 109, 90, 24, 88 }, 2389f046c2dSEric Chamberland {118, 111, 116, 92, 113, 117}, 2399f046c2dSEric Chamberland {110, 116, 108, 79, 117, 109}, 2409f046c2dSEric Chamberland {92, 113, 117, 95, 94, 93 }, 2419f046c2dSEric Chamberland {79, 117, 109, 80, 93, 88 }, 2429f046c2dSEric Chamberland {22, 18, 63, 19, 20, 21 }, 2439f046c2dSEric Chamberland {68, 63, 64, 66, 21, 61 }, 2449f046c2dSEric Chamberland {19, 20, 21, 99, 97, 98 }, 2459f046c2dSEric Chamberland {66, 21, 61, 124, 98, 119}, 2469f046c2dSEric Chamberland {6, 1, 59, 2, 4, 5 }, 2479f046c2dSEric Chamberland {62, 59, 50, 60, 5, 53 }, 2489f046c2dSEric Chamberland {2, 4, 5, 73, 71, 72 }, 2499f046c2dSEric Chamberland {60, 5, 53, 121, 72, 112}, 2509f046c2dSEric Chamberland {17, 12, 48, 13, 16, 14 }, 2519f046c2dSEric Chamberland {54, 48, 65, 51, 14, 67 }, 2529f046c2dSEric Chamberland {13, 16, 14, 85, 84, 82 }, 2539f046c2dSEric Chamberland {51, 14, 67, 115, 82, 122}, 2549f046c2dSEric Chamberland {62, 50, 64, 60, 53, 61 }, 2559f046c2dSEric Chamberland {68, 64, 65, 66, 61, 67 }, 2569f046c2dSEric Chamberland {60, 53, 61, 121, 112, 119}, 2579f046c2dSEric Chamberland {66, 61, 67, 124, 119, 122}, 2589f046c2dSEric Chamberland {99, 97, 98, 96, 39, 40 }, 2599f046c2dSEric Chamberland {124, 98, 119, 106, 40, 120}, 2609f046c2dSEric Chamberland {96, 39, 40, 41, 38, 105}, 2619f046c2dSEric Chamberland {106, 40, 120, 107, 105, 102}, 2629f046c2dSEric Chamberland {73, 71, 72, 69, 26, 27 }, 2639f046c2dSEric Chamberland {121, 72, 112, 101, 27, 114}, 2649f046c2dSEric Chamberland {69, 26, 27, 28, 24, 100}, 2659f046c2dSEric Chamberland {101, 27, 114, 104, 100, 103}, 2669f046c2dSEric Chamberland {85, 84, 82, 81, 36, 34 }, 2679f046c2dSEric Chamberland {115, 82, 122, 87, 34, 123}, 2689f046c2dSEric Chamberland {81, 36, 34, 37, 33, 86 }, 2699f046c2dSEric Chamberland {87, 34, 123, 90, 86, 89 }, 2709f046c2dSEric Chamberland {121, 112, 119, 101, 114, 120}, 2719f046c2dSEric Chamberland {124, 119, 122, 106, 120, 123}, 2729f046c2dSEric Chamberland {101, 114, 120, 104, 103, 102}, 2739f046c2dSEric Chamberland {106, 120, 123, 107, 102, 89 }, 2749f046c2dSEric Chamberland /* rank 1 */ 2759f046c2dSEric Chamberland {58, 43, 7, 56, 57, 10 }, 2769f046c2dSEric Chamberland {7, 43, 42, 10, 57, 9 }, 2779f046c2dSEric Chamberland {56, 57, 10, 118, 116, 76 }, 2789f046c2dSEric Chamberland {10, 57, 9, 76, 116, 75 }, 2799f046c2dSEric Chamberland {54, 49, 48, 51, 52, 14 }, 2809f046c2dSEric Chamberland {48, 49, 55, 14, 52, 15 }, 2819f046c2dSEric Chamberland {51, 52, 14, 115, 111, 82 }, 2829f046c2dSEric Chamberland {14, 52, 15, 82, 111, 83 }, 2839f046c2dSEric Chamberland {47, 44, 0, 45, 46, 3 }, 2849f046c2dSEric Chamberland {0, 44, 1, 3, 46, 4 }, 2859f046c2dSEric Chamberland {45, 46, 3, 110, 108, 70 }, 2869f046c2dSEric Chamberland {3, 46, 4, 70, 108, 71 }, 2879f046c2dSEric Chamberland {54, 44, 49, 51, 46, 52 }, 2889f046c2dSEric Chamberland {49, 44, 43, 52, 46, 57 }, 2899f046c2dSEric Chamberland {51, 46, 52, 115, 108, 111}, 2909f046c2dSEric Chamberland {52, 46, 57, 111, 108, 116}, 2919f046c2dSEric Chamberland {118, 116, 76, 92, 117, 31 }, 2929f046c2dSEric Chamberland {76, 116, 75, 31, 117, 30 }, 2939f046c2dSEric Chamberland {92, 117, 31, 95, 93, 29 }, 2949f046c2dSEric Chamberland {31, 117, 30, 29, 93, 78 }, 2959f046c2dSEric Chamberland {115, 111, 82, 87, 113, 34 }, 2969f046c2dSEric Chamberland {82, 111, 83, 34, 113, 35 }, 2979f046c2dSEric Chamberland {87, 113, 34, 90, 94, 86 }, 2989f046c2dSEric Chamberland {34, 113, 35, 86, 94, 91 }, 2999f046c2dSEric Chamberland {110, 108, 70, 79, 109, 25 }, 3009f046c2dSEric Chamberland {70, 108, 71, 25, 109, 26 }, 3019f046c2dSEric Chamberland {79, 109, 25, 80, 88, 23 }, 3029f046c2dSEric Chamberland {25, 109, 26, 23, 88, 24 }, 3039f046c2dSEric Chamberland {115, 108, 111, 87, 109, 113}, 3049f046c2dSEric Chamberland {111, 108, 116, 113, 109, 117}, 3059f046c2dSEric Chamberland {87, 109, 113, 90, 88, 94 }, 3069f046c2dSEric Chamberland {113, 109, 117, 94, 88, 93 }, 3079f046c2dSEric Chamberland {62, 64, 18, 60, 61, 20 }, 3089f046c2dSEric Chamberland {18, 64, 63, 20, 61, 21 }, 3099f046c2dSEric Chamberland {60, 61, 20, 121, 119, 97 }, 3109f046c2dSEric Chamberland {20, 61, 21, 97, 119, 98 }, 3119f046c2dSEric Chamberland {54, 50, 1, 51, 53, 4 }, 3129f046c2dSEric Chamberland {1, 50, 59, 4, 53, 5 }, 3139f046c2dSEric Chamberland {51, 53, 4, 115, 112, 71 }, 3149f046c2dSEric Chamberland {4, 53, 5, 71, 112, 72 }, 3159f046c2dSEric Chamberland {68, 65, 12, 66, 67, 16 }, 3169f046c2dSEric Chamberland {12, 65, 48, 16, 67, 14 }, 3179f046c2dSEric Chamberland {66, 67, 16, 124, 122, 84 }, 3189f046c2dSEric Chamberland {16, 67, 14, 84, 122, 82 }, 3199f046c2dSEric Chamberland {54, 65, 50, 51, 67, 53 }, 3209f046c2dSEric Chamberland {50, 65, 64, 53, 67, 61 }, 3219f046c2dSEric Chamberland {51, 67, 53, 115, 122, 112}, 3229f046c2dSEric Chamberland {53, 67, 61, 112, 122, 119}, 3239f046c2dSEric Chamberland {121, 119, 97, 101, 120, 39 }, 3249f046c2dSEric Chamberland {97, 119, 98, 39, 120, 40 }, 3259f046c2dSEric Chamberland {101, 120, 39, 104, 102, 38 }, 3269f046c2dSEric Chamberland {39, 120, 40, 38, 102, 105}, 3279f046c2dSEric Chamberland {115, 112, 71, 87, 114, 26 }, 3289f046c2dSEric Chamberland {71, 112, 72, 26, 114, 27 }, 3299f046c2dSEric Chamberland {87, 114, 26, 90, 103, 24 }, 3309f046c2dSEric Chamberland {26, 114, 27, 24, 103, 100}, 3319f046c2dSEric Chamberland {124, 122, 84, 106, 123, 36 }, 3329f046c2dSEric Chamberland {84, 122, 82, 36, 123, 34 }, 3339f046c2dSEric Chamberland {106, 123, 36, 107, 89, 33 }, 3349f046c2dSEric Chamberland {36, 123, 34, 33, 89, 86 }, 3359f046c2dSEric Chamberland {115, 122, 112, 87, 123, 114}, 3369f046c2dSEric Chamberland {112, 122, 119, 114, 123, 120}, 3379f046c2dSEric Chamberland {87, 123, 114, 90, 89, 103}, 3389f046c2dSEric Chamberland {114, 123, 120, 103, 89, 102} 3399f046c2dSEric Chamberland }; 3409f046c2dSEric Chamberland 3419f046c2dSEric Chamberland /* Partitions of prisms mesh : */ 3429f046c2dSEric Chamberland const PetscInt sInitialPartitionPrismsMesh[2][64] = { 3439f046c2dSEric Chamberland {0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63 }, 3449f046c2dSEric Chamberland {64, 65, 66, 67, 68, 69, 70, 71, 72, 73, 74, 75, 76, 77, 78, 79, 80, 81, 82, 83, 84, 85, 86, 87, 88, 89, 90, 91, 92, 93, 94, 95, 3459f046c2dSEric Chamberland 96, 97, 98, 99, 100, 101, 102, 103, 104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127} 3469f046c2dSEric Chamberland }; 3479f046c2dSEric Chamberland 3489f046c2dSEric Chamberland const PetscInt sNLoclCellsPrismsMesh = 64; 3499f046c2dSEric Chamberland const PetscInt sNGlobVertsPrismsMesh = 125; 3509f046c2dSEric Chamberland 351d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 352d71ae5a4SJacob Faibussowitsch { 3539f046c2dSEric Chamberland PetscInt Nc = 0; 3549f046c2dSEric Chamberland const PetscInt *InitPartForRank[2]; 355231de6a2SMatthew G. Knepley DM dm, idm, ddm; 356231de6a2SMatthew G. Knepley PetscSF sfVert, sfMig, sfPart; 357231de6a2SMatthew G. Knepley PetscPartitioner part; 358231de6a2SMatthew G. Knepley PetscSection s; 359231de6a2SMatthew G. Knepley PetscInt *cells, c; 360231de6a2SMatthew G. Knepley PetscMPIInt size, rank; 36172f9f102SEric Chamberland PetscBool box = PETSC_FALSE, field = PETSC_FALSE, quadsmesh = PETSC_FALSE, trisquadsmesh = PETSC_FALSE, prismsmesh = PETSC_FALSE; 362231de6a2SMatthew G. Knepley 363327415f7SBarry Smith PetscFunctionBeginUser; 364231de6a2SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 365231de6a2SMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 366231de6a2SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 367231de6a2SMatthew G. Knepley PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a 2 processors example only"); 368231de6a2SMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-field", &field, NULL)); 3699f046c2dSEric Chamberland PetscCall(PetscOptionsGetBool(NULL, NULL, "-box", &box, NULL)); 3709f046c2dSEric Chamberland PetscCall(PetscOptionsGetBool(NULL, NULL, "-quadsmesh", &quadsmesh, NULL)); 37172f9f102SEric Chamberland PetscCall(PetscOptionsGetBool(NULL, NULL, "-trisquadsmesh", &trisquadsmesh, NULL)); 3729f046c2dSEric Chamberland PetscCall(PetscOptionsGetBool(NULL, NULL, "-prismsmesh", &prismsmesh, NULL)); 37372f9f102SEric Chamberland PetscCheck(1 == (box ? 1 : 0) + (quadsmesh ? 1 : 0) + (trisquadsmesh ? 1 : 0) + (prismsmesh ? 1 : 0), PETSC_COMM_WORLD, PETSC_ERR_SUP, "Specify one and only one of -box, -quadsmesh or -prismsmesh"); 374231de6a2SMatthew G. Knepley 375231de6a2SMatthew G. Knepley PetscCall(DMPlexCreate(PETSC_COMM_WORLD, &dm)); 376231de6a2SMatthew G. Knepley if (box) { 377231de6a2SMatthew G. Knepley PetscCall(DMSetType(dm, DMPLEX)); 378231de6a2SMatthew G. Knepley PetscCall(DMSetFromOptions(dm)); 379231de6a2SMatthew G. Knepley } else { 3809f046c2dSEric Chamberland if (quadsmesh) { 3819f046c2dSEric Chamberland Nc = sNLoclCells2x5Mesh; //Same on each rank for this example... 3829f046c2dSEric Chamberland PetscInt Nv = sNGlobVerts2x5Mesh; 3839f046c2dSEric Chamberland InitPartForRank[0] = &sInitialPartition2x5Mesh[0][0]; 3849f046c2dSEric Chamberland InitPartForRank[1] = &sInitialPartition2x5Mesh[1][0]; 3859f046c2dSEric Chamberland const PetscInt(*Conn)[4] = sConnectivity2x5Mesh; 3869f046c2dSEric Chamberland 3879f046c2dSEric Chamberland const PetscInt Ncor = 4; 3889f046c2dSEric Chamberland const PetscInt dim = 2; 3899f046c2dSEric Chamberland 390231de6a2SMatthew G. Knepley PetscCall(PetscMalloc1(Nc * Ncor, &cells)); 391231de6a2SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 392231de6a2SMatthew G. Knepley PetscInt cell = (InitPartForRank[rank])[c], cor; 393231de6a2SMatthew G. Knepley 394ad540459SPierre Jolivet for (cor = 0; cor < Ncor; ++cor) cells[c * Ncor + cor] = Conn[cell][cor]; 395231de6a2SMatthew G. Knepley } 396231de6a2SMatthew G. Knepley PetscCall(DMSetDimension(dm, dim)); 397231de6a2SMatthew G. Knepley PetscCall(DMPlexBuildFromCellListParallel(dm, Nc, PETSC_DECIDE, Nv, Ncor, cells, &sfVert, NULL)); 39872f9f102SEric Chamberland } else if (trisquadsmesh) { 39972f9f102SEric Chamberland Nc = sNLoclCellsMixedTQMesh; //Same on each rank for this example... 40072f9f102SEric Chamberland PetscInt Nv = sNGlobVertsMixedTQMesh; 40172f9f102SEric Chamberland InitPartForRank[0] = &sInitialPartitionMixedTQMesh[0][0]; 40272f9f102SEric Chamberland InitPartForRank[1] = &sInitialPartitionMixedTQMesh[1][0]; 40372f9f102SEric Chamberland const PetscInt(*Conn)[4] = sConnectivityMixedTQMesh; 40472f9f102SEric Chamberland 40572f9f102SEric Chamberland const PetscInt NcorMax = 4; 40672f9f102SEric Chamberland const PetscInt dim = 2; 40772f9f102SEric Chamberland 408d7c1f440SPierre Jolivet /* Create a PetscSection and taking care to exclude nodes with "-1" into element connectivity: */ 40972f9f102SEric Chamberland PetscSection s; 41072f9f102SEric Chamberland PetscInt vStart = 0, vEnd = Nc; 41172f9f102SEric Chamberland PetscCall(PetscSectionCreate(PETSC_COMM_WORLD, &s)); 41272f9f102SEric Chamberland PetscCall(PetscSectionSetNumFields(s, 1)); 41372f9f102SEric Chamberland PetscCall(PetscSectionSetFieldComponents(s, 0, 1)); 41472f9f102SEric Chamberland PetscCall(PetscSectionSetChart(s, vStart, vEnd)); 41572f9f102SEric Chamberland 41672f9f102SEric Chamberland PetscCall(PetscMalloc1(Nc * NcorMax, &cells)); 41772f9f102SEric Chamberland PetscInt count = 0; 41872f9f102SEric Chamberland for (c = 0; c < Nc; ++c) { 41972f9f102SEric Chamberland PetscInt cell = (InitPartForRank[rank])[c], cor; 42072f9f102SEric Chamberland PetscInt nbElemVertex = ((-1 == Conn[cell][NcorMax - 1]) ? 3 : 4); 42172f9f102SEric Chamberland for (cor = 0; cor < nbElemVertex; ++cor) { 42272f9f102SEric Chamberland cells[count] = Conn[cell][cor]; 42372f9f102SEric Chamberland ++count; 42472f9f102SEric Chamberland } 42572f9f102SEric Chamberland PetscCall(PetscSectionSetDof(s, c, nbElemVertex)); 42672f9f102SEric Chamberland PetscCall(PetscSectionSetFieldDof(s, c, 0, nbElemVertex)); 42772f9f102SEric Chamberland } 42872f9f102SEric Chamberland PetscCall(PetscSectionSetUp(s)); 42972f9f102SEric Chamberland PetscCall(DMSetDimension(dm, dim)); 43072f9f102SEric Chamberland PetscCall(DMPlexBuildFromCellSectionParallel(dm, Nc, PETSC_DECIDE, Nv, s, cells, &sfVert, NULL)); 43172f9f102SEric Chamberland PetscCall(PetscSectionDestroy(&s)); 4329f046c2dSEric Chamberland } else if (prismsmesh) { 4339f046c2dSEric Chamberland Nc = sNLoclCellsPrismsMesh; //Same on each rank for this example... 4349f046c2dSEric Chamberland PetscInt Nv = sNGlobVertsPrismsMesh; 4359f046c2dSEric Chamberland InitPartForRank[0] = &sInitialPartitionPrismsMesh[0][0]; 4369f046c2dSEric Chamberland InitPartForRank[1] = &sInitialPartitionPrismsMesh[1][0]; 4379f046c2dSEric Chamberland const PetscInt(*Conn)[6] = sConnectivityPrismsMesh; 4389f046c2dSEric Chamberland 4399f046c2dSEric Chamberland const PetscInt Ncor = 6; 4409f046c2dSEric Chamberland const PetscInt dim = 3; 4419f046c2dSEric Chamberland 4429f046c2dSEric Chamberland PetscCall(PetscMalloc1(Nc * Ncor, &cells)); 4439f046c2dSEric Chamberland for (c = 0; c < Nc; ++c) { 4449f046c2dSEric Chamberland PetscInt cell = (InitPartForRank[rank])[c], cor; 4459f046c2dSEric Chamberland 4469f046c2dSEric Chamberland for (cor = 0; cor < Ncor; ++cor) cells[c * Ncor + cor] = Conn[cell][cor]; 4479f046c2dSEric Chamberland } 4489f046c2dSEric Chamberland PetscCall(DMSetDimension(dm, dim)); 4499f046c2dSEric Chamberland PetscCall(DMPlexBuildFromCellListParallel(dm, Nc, PETSC_DECIDE, Nv, Ncor, cells, &sfVert, NULL)); 4509f046c2dSEric Chamberland } 451231de6a2SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfVert)); 452231de6a2SMatthew G. Knepley PetscCall(PetscFree(cells)); 453231de6a2SMatthew G. Knepley PetscCall(DMPlexInterpolate(dm, &idm)); 454231de6a2SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 455231de6a2SMatthew G. Knepley dm = idm; 456231de6a2SMatthew G. Knepley } 457231de6a2SMatthew G. Knepley PetscCall(DMSetUseNatural(dm, PETSC_TRUE)); 458231de6a2SMatthew G. Knepley PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 459231de6a2SMatthew G. Knepley 460231de6a2SMatthew G. Knepley if (field) { 461231de6a2SMatthew G. Knepley const PetscInt Nf = 1; 462231de6a2SMatthew G. Knepley const PetscInt numBC = 0; 4639f046c2dSEric Chamberland const PetscInt numComp[1] = {1}; 4649f046c2dSEric Chamberland PetscInt numDof[4] = {0, 0, 0, 0}; 4659f046c2dSEric Chamberland PetscInt dim; 4669f046c2dSEric Chamberland 4679f046c2dSEric Chamberland PetscCall(DMGetDimension(dm, &dim)); 4689f046c2dSEric Chamberland numDof[dim] = 1; 469231de6a2SMatthew G. Knepley 470231de6a2SMatthew G. Knepley PetscCall(DMSetNumFields(dm, Nf)); 471231de6a2SMatthew G. Knepley PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &s)); 472231de6a2SMatthew G. Knepley PetscCall(DMSetLocalSection(dm, s)); 473231de6a2SMatthew G. Knepley PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD)); 474231de6a2SMatthew G. Knepley PetscCall(PetscSectionDestroy(&s)); 475231de6a2SMatthew G. Knepley } 476231de6a2SMatthew G. Knepley 477231de6a2SMatthew G. Knepley PetscCall(DMPlexGetPartitioner(dm, &part)); 478231de6a2SMatthew G. Knepley PetscCall(PetscPartitionerSetFromOptions(part)); 479231de6a2SMatthew G. Knepley 480231de6a2SMatthew G. Knepley PetscCall(DMPlexDistribute(dm, 0, &sfMig, &ddm)); 481231de6a2SMatthew G. Knepley PetscCall(PetscSFView(sfMig, PETSC_VIEWER_STDOUT_WORLD)); 482231de6a2SMatthew G. Knepley PetscCall(PetscSFCreateInverseSF(sfMig, &sfPart)); 483231de6a2SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sfPart, "Inverse Migration SF")); 484231de6a2SMatthew G. Knepley PetscCall(PetscSFView(sfPart, PETSC_VIEWER_STDOUT_WORLD)); 485231de6a2SMatthew G. Knepley 486231de6a2SMatthew G. Knepley Vec lGlobalVec, lNatVec; 487231de6a2SMatthew G. Knepley PetscScalar *lNatVecArray; 488231de6a2SMatthew G. Knepley 489231de6a2SMatthew G. Knepley { 490231de6a2SMatthew G. Knepley PetscSection s; 491231de6a2SMatthew G. Knepley 492231de6a2SMatthew G. Knepley PetscCall(DMGetGlobalSection(dm, &s)); 493231de6a2SMatthew G. Knepley PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD)); 494231de6a2SMatthew G. Knepley } 495231de6a2SMatthew G. Knepley PetscCall(DMGetGlobalVector(dm, &lNatVec)); 496231de6a2SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)lNatVec, "Natural Vector (initial partition)")); 497231de6a2SMatthew G. Knepley 498231de6a2SMatthew G. Knepley //Copying the initial partition into the "natural" vector: 4999f046c2dSEric Chamberland PetscCall(VecZeroEntries(lNatVec)); 500231de6a2SMatthew G. Knepley PetscCall(VecGetArray(lNatVec, &lNatVecArray)); 501231de6a2SMatthew G. Knepley for (c = 0; c < Nc; ++c) lNatVecArray[c] = (InitPartForRank[rank])[c]; 502231de6a2SMatthew G. Knepley PetscCall(VecRestoreArray(lNatVec, &lNatVecArray)); 503231de6a2SMatthew G. Knepley 504231de6a2SMatthew G. Knepley PetscCall(DMGetGlobalVector(ddm, &lGlobalVec)); 505231de6a2SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)lGlobalVec, "Global Vector (reordered element numbers in the petsc distributed order)")); 506231de6a2SMatthew G. Knepley PetscCall(VecZeroEntries(lGlobalVec)); 507231de6a2SMatthew G. Knepley 508231de6a2SMatthew G. Knepley // The call to DMPlexNaturalToGlobalBegin/End does not produce our expected result... 509231de6a2SMatthew G. Knepley // In lGlobalVec, we expect to have: 510231de6a2SMatthew G. Knepley /* 511231de6a2SMatthew G. Knepley * Process [0] 512231de6a2SMatthew G. Knepley * 2. 513231de6a2SMatthew G. Knepley * 4. 514231de6a2SMatthew G. Knepley * 8. 515231de6a2SMatthew G. Knepley * 3. 516231de6a2SMatthew G. Knepley * 9. 517231de6a2SMatthew G. Knepley * Process [1] 518231de6a2SMatthew G. Knepley * 1. 519231de6a2SMatthew G. Knepley * 5. 520231de6a2SMatthew G. Knepley * 7. 521231de6a2SMatthew G. Knepley * 0. 522231de6a2SMatthew G. Knepley * 6. 523231de6a2SMatthew G. Knepley * 524231de6a2SMatthew G. Knepley * but we obtained: 525231de6a2SMatthew G. Knepley * 526231de6a2SMatthew G. Knepley * Process [0] 527231de6a2SMatthew G. Knepley * 2. 528231de6a2SMatthew G. Knepley * 4. 529231de6a2SMatthew G. Knepley * 8. 530231de6a2SMatthew G. Knepley * 0. 531231de6a2SMatthew G. Knepley * 0. 532231de6a2SMatthew G. Knepley * Process [1] 533231de6a2SMatthew G. Knepley * 0. 534231de6a2SMatthew G. Knepley * 0. 535231de6a2SMatthew G. Knepley * 0. 536231de6a2SMatthew G. Knepley * 0. 537231de6a2SMatthew G. Knepley * 0. 538231de6a2SMatthew G. Knepley */ 539231de6a2SMatthew G. Knepley 540231de6a2SMatthew G. Knepley { 541231de6a2SMatthew G. Knepley PetscSF nsf; 542231de6a2SMatthew G. Knepley 543*e535cce4SJames Wright PetscCall(DMGetNaturalSF(ddm, &nsf)); 544231de6a2SMatthew G. Knepley PetscCall(PetscSFView(nsf, NULL)); 545231de6a2SMatthew G. Knepley } 546231de6a2SMatthew G. Knepley PetscCall(DMPlexNaturalToGlobalBegin(ddm, lNatVec, lGlobalVec)); 547231de6a2SMatthew G. Knepley PetscCall(DMPlexNaturalToGlobalEnd(ddm, lNatVec, lGlobalVec)); 548231de6a2SMatthew G. Knepley 549231de6a2SMatthew G. Knepley PetscCall(VecView(lNatVec, PETSC_VIEWER_STDOUT_WORLD)); 550231de6a2SMatthew G. Knepley PetscCall(VecView(lGlobalVec, PETSC_VIEWER_STDOUT_WORLD)); 551231de6a2SMatthew G. Knepley 552231de6a2SMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dm, &lNatVec)); 553231de6a2SMatthew G. Knepley PetscCall(DMRestoreGlobalVector(ddm, &lGlobalVec)); 554231de6a2SMatthew G. Knepley 5559f046c2dSEric Chamberland const PetscBool lUseCone = PETSC_FALSE; 5569f046c2dSEric Chamberland const PetscBool lUseClosure = PETSC_TRUE; 5579f046c2dSEric Chamberland PetscCall(DMSetBasicAdjacency(ddm, lUseCone, lUseClosure)); 5589f046c2dSEric Chamberland const PetscInt lNbCellsInOverlap = 1; 5599f046c2dSEric Chamberland PetscSF lSFMigrationOvl; 5609f046c2dSEric Chamberland DM ddm_with_overlap; 5619f046c2dSEric Chamberland 5629f046c2dSEric Chamberland PetscCall(DMPlexDistributeOverlap(ddm, lNbCellsInOverlap, &lSFMigrationOvl, &ddm_with_overlap)); 5639f046c2dSEric Chamberland 5649f046c2dSEric Chamberland IS lISCellWithOvl = 0; 5659f046c2dSEric Chamberland /* This is the buggy call with prisms since commit 5ae96e2b862 */ 566e2b8d0fcSMatthew G. Knepley PetscCall(DMPlexCreateCellNumbering(ddm_with_overlap, PETSC_TRUE, &lISCellWithOvl)); 5679f046c2dSEric Chamberland /* Here, we can see the elements in the overlap within the IS: they are the ones with negative indices */ 5689f046c2dSEric Chamberland PetscCall(ISView(lISCellWithOvl, PETSC_VIEWER_STDOUT_WORLD)); 569e2b8d0fcSMatthew G. Knepley PetscCall(ISDestroy(&lISCellWithOvl)); 5709f046c2dSEric Chamberland 5719f046c2dSEric Chamberland PetscCall(PetscSFDestroy(&lSFMigrationOvl)); 5729f046c2dSEric Chamberland PetscCall(DMDestroy(&ddm_with_overlap)); 573231de6a2SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfMig)); 574231de6a2SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfPart)); 575231de6a2SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 576231de6a2SMatthew G. Knepley PetscCall(DMDestroy(&ddm)); 577231de6a2SMatthew G. Knepley PetscCall(PetscFinalize()); 578231de6a2SMatthew G. Knepley return 0; 579231de6a2SMatthew G. Knepley } 580231de6a2SMatthew G. Knepley 581231de6a2SMatthew G. Knepley /*TEST 582231de6a2SMatthew G. Knepley 583231de6a2SMatthew G. Knepley testset: 584231de6a2SMatthew G. Knepley args: -field -petscpartitioner_type simple 585231de6a2SMatthew G. Knepley nsize: 2 586231de6a2SMatthew G. Knepley 587231de6a2SMatthew G. Knepley test: 588231de6a2SMatthew G. Knepley suffix: 0 5899f046c2dSEric Chamberland args: -quadsmesh 5909f046c2dSEric Chamberland output_file: output/ex47_0.out 591231de6a2SMatthew G. Knepley 592231de6a2SMatthew G. Knepley test: 593231de6a2SMatthew G. Knepley suffix: 1 594231de6a2SMatthew G. Knepley args: -box -dm_plex_simplex 0 -dm_plex_box_faces 2,5 -dm_distribute 5959f046c2dSEric Chamberland output_file: output/ex47_1.out 5969f046c2dSEric Chamberland 5979f046c2dSEric Chamberland test: 5989f046c2dSEric Chamberland suffix: 2 5999f046c2dSEric Chamberland args: -prismsmesh 6009f046c2dSEric Chamberland output_file: output/ex47_2.out 601231de6a2SMatthew G. Knepley 60272f9f102SEric Chamberland test: 60372f9f102SEric Chamberland suffix: 3 60472f9f102SEric Chamberland args: -trisquadsmesh 60572f9f102SEric Chamberland output_file: output/ex47_3.out 60672f9f102SEric Chamberland 607231de6a2SMatthew G. Knepley TEST*/ 608