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 6*9f046c2dSEric 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 28*9f046c2dSEric 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 42*9f046c2dSEric 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 51*9f046c2dSEric Chamberland /* Prisms mesh */ 52*9f046c2dSEric Chamberland PetscReal sCoordsPrismsMesh[125][3] = { 53*9f046c2dSEric Chamberland {2.24250931694056355e-01, 0.00000000000000000e+00, 0.00000000000000000e+00}, 54*9f046c2dSEric Chamberland {2.20660660151932697e-01, 2.87419338850266937e-01, 0.00000000000000000e+00}, 55*9f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 2.70243537720639027e-01}, 56*9f046c2dSEric Chamberland {2.32445727460992402e-01, 0.00000000000000000e+00, 2.60591845015572310e-01}, 57*9f046c2dSEric Chamberland {2.41619971105419079e-01, 2.69894910706158231e-01, 2.42844781736072490e-01}, 58*9f046c2dSEric Chamberland {0.00000000000000000e+00, 2.46523339883120779e-01, 2.69072907562752262e-01}, 59*9f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 0.00000000000000000e+00}, 60*9f046c2dSEric Chamberland {1.00000000000000000e+00, 2.75433417601945563e-01, 0.00000000000000000e+00}, 61*9f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 2.33748605950385602e-01}, 62*9f046c2dSEric Chamberland {7.32445727460992457e-01, 0.00000000000000000e+00, 2.42344379130445597e-01}, 63*9f046c2dSEric Chamberland {1.00000000000000000e+00, 2.78258478013028610e-01, 2.57379172987105553e-01}, 64*9f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 0.00000000000000000e+00}, 65*9f046c2dSEric Chamberland {7.49586880891153995e-01, 1.00000000000000000e+00, 0.00000000000000000e+00}, 66*9f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 2.51949651266657582e-01}, 67*9f046c2dSEric Chamberland {7.41619971105419107e-01, 7.69894910706158120e-01, 2.33697838509081768e-01}, 68*9f046c2dSEric Chamberland {1.00000000000000000e+00, 7.78258478013028610e-01, 2.66479695645241543e-01}, 69*9f046c2dSEric Chamberland {7.55042653233710115e-01, 1.00000000000000000e+00, 2.58019637386860512e-01}, 70*9f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 0.00000000000000000e+00}, 71*9f046c2dSEric Chamberland {0.00000000000000000e+00, 7.59235710423095789e-01, 0.00000000000000000e+00}, 72*9f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 2.17232187874490473e-01}, 73*9f046c2dSEric Chamberland {0.00000000000000000e+00, 7.46523339883120807e-01, 2.42567232639677999e-01}, 74*9f046c2dSEric Chamberland {2.55042653233710115e-01, 1.00000000000000000e+00, 2.40660905690776916e-01}, 75*9f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 0.00000000000000000e+00}, 76*9f046c2dSEric Chamberland {2.38934376044866809e-01, 0.00000000000000000e+00, 1.00000000000000000e+00}, 77*9f046c2dSEric Chamberland {2.18954188589218168e-01, 2.26916038449581692e-01, 1.00000000000000000e+00}, 78*9f046c2dSEric Chamberland {2.39787449636397643e-01, 0.00000000000000000e+00, 7.60591845015572310e-01}, 79*9f046c2dSEric Chamberland {2.40766735324061815e-01, 2.39643260505815608e-01, 7.42844781736072490e-01}, 80*9f046c2dSEric Chamberland {0.00000000000000000e+00, 2.57448248192627016e-01, 7.69072907562752262e-01}, 81*9f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 1.00000000000000000e+00}, 82*9f046c2dSEric Chamberland {1.00000000000000000e+00, 2.38666970143638080e-01, 1.00000000000000000e+00}, 83*9f046c2dSEric Chamberland {7.39787449636397643e-01, 0.00000000000000000e+00, 7.42344379130445597e-01}, 84*9f046c2dSEric Chamberland {1.00000000000000000e+00, 2.59875254283874868e-01, 7.57379172987105553e-01}, 85*9f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 1.00000000000000000e+00}, 86*9f046c2dSEric Chamberland {7.76318984007844159e-01, 1.00000000000000000e+00, 1.00000000000000000e+00}, 87*9f046c2dSEric Chamberland {7.40766735324061787e-01, 7.39643260505815636e-01, 7.33697838509081768e-01}, 88*9f046c2dSEric Chamberland {1.00000000000000000e+00, 7.59875254283874924e-01, 7.66479695645241543e-01}, 89*9f046c2dSEric Chamberland {7.68408704792055142e-01, 1.00000000000000000e+00, 7.58019637386860512e-01}, 90*9f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 1.00000000000000000e+00}, 91*9f046c2dSEric Chamberland {0.00000000000000000e+00, 7.81085527042108207e-01, 1.00000000000000000e+00}, 92*9f046c2dSEric Chamberland {0.00000000000000000e+00, 7.57448248192627016e-01, 7.42567232639678054e-01}, 93*9f046c2dSEric Chamberland {2.68408704792055197e-01, 1.00000000000000000e+00, 7.40660905690776916e-01}, 94*9f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 1.00000000000000000e+00}, 95*9f046c2dSEric Chamberland {7.24250931694056410e-01, 0.00000000000000000e+00, 0.00000000000000000e+00}, 96*9f046c2dSEric Chamberland {7.24250931694056410e-01, 2.75433417601945563e-01, 0.00000000000000000e+00}, 97*9f046c2dSEric Chamberland {4.44911591845989052e-01, 2.87419338850266937e-01, 0.00000000000000000e+00}, 98*9f046c2dSEric Chamberland {4.64891454921984804e-01, 0.00000000000000000e+00, 2.50940152310505593e-01}, 99*9f046c2dSEric Chamberland {4.74065698566411453e-01, 2.69894910706158231e-01, 2.33193089031005774e-01}, 100*9f046c2dSEric Chamberland {4.48501863388112709e-01, 0.00000000000000000e+00, 0.00000000000000000e+00}, 101*9f046c2dSEric Chamberland {7.20660660151932753e-01, 7.87419338850266937e-01, 0.00000000000000000e+00}, 102*9f046c2dSEric Chamberland {7.20660660151932753e-01, 5.62852756452212555e-01, 0.00000000000000000e+00}, 103*9f046c2dSEric Chamberland {2.20660660151932697e-01, 5.46655049273362614e-01, 0.00000000000000000e+00}, 104*9f046c2dSEric Chamberland {4.83239942210838158e-01, 5.39789821412316462e-01, 2.15446025751505982e-01}, 105*9f046c2dSEric Chamberland {7.41619971105419107e-01, 5.48153388719186951e-01, 2.48227882887665785e-01}, 106*9f046c2dSEric Chamberland {2.41619971105419079e-01, 5.16418250589278927e-01, 2.41674151578185781e-01}, 107*9f046c2dSEric Chamberland {4.41321320303865394e-01, 5.74838677700533873e-01, 0.00000000000000000e+00}, 108*9f046c2dSEric Chamberland {1.00000000000000000e+00, 7.75433417601945507e-01, 0.00000000000000000e+00}, 109*9f046c2dSEric Chamberland {1.00000000000000000e+00, 5.56516956026057219e-01, 2.81009740023825560e-01}, 110*9f046c2dSEric Chamberland {7.32445727460992457e-01, 2.78258478013028610e-01, 2.65974946167165549e-01}, 111*9f046c2dSEric Chamberland {1.00000000000000000e+00, 5.50866835203891125e-01, 0.00000000000000000e+00}, 112*9f046c2dSEric Chamberland {0.00000000000000000e+00, 2.59235710423095733e-01, 0.00000000000000000e+00}, 113*9f046c2dSEric Chamberland {0.00000000000000000e+00, 4.93046679766241558e-01, 2.67902277404865552e-01}, 114*9f046c2dSEric Chamberland {2.55042653233710115e-01, 7.46523339883120807e-01, 2.65995950455964469e-01}, 115*9f046c2dSEric Chamberland {0.00000000000000000e+00, 5.18471420846191466e-01, 0.00000000000000000e+00}, 116*9f046c2dSEric Chamberland {2.49586880891154023e-01, 1.00000000000000000e+00, 0.00000000000000000e+00}, 117*9f046c2dSEric Chamberland {2.49586880891154023e-01, 7.59235710423095789e-01, 0.00000000000000000e+00}, 118*9f046c2dSEric Chamberland {4.70247541043086748e-01, 7.87419338850266937e-01, 0.00000000000000000e+00}, 119*9f046c2dSEric Chamberland {5.10085306467420230e-01, 1.00000000000000000e+00, 2.64089623507063387e-01}, 120*9f046c2dSEric Chamberland {4.96662624339129222e-01, 7.69894910706158231e-01, 2.39767824629284698e-01}, 121*9f046c2dSEric Chamberland {4.99173761782308045e-01, 1.00000000000000000e+00, 0.00000000000000000e+00}, 122*9f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 7.70243537720639027e-01}, 123*9f046c2dSEric Chamberland {2.40640523227928449e-01, 0.00000000000000000e+00, 5.21183690031144620e-01}, 124*9f046c2dSEric Chamberland {2.62579282058905461e-01, 2.52370482562049525e-01, 4.85689563472144981e-01}, 125*9f046c2dSEric Chamberland {0.00000000000000000e+00, 2.33810969343145825e-01, 5.38145815125504523e-01}, 126*9f046c2dSEric Chamberland {0.00000000000000000e+00, 0.00000000000000000e+00, 5.40487075441278053e-01}, 127*9f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 7.33748605950385602e-01}, 128*9f046c2dSEric Chamberland {7.40640523227928504e-01, 0.00000000000000000e+00, 4.84688758260891195e-01}, 129*9f046c2dSEric Chamberland {1.00000000000000000e+00, 2.81083538424111656e-01, 5.14758345974211107e-01}, 130*9f046c2dSEric Chamberland {1.00000000000000000e+00, 0.00000000000000000e+00, 4.67497211900771203e-01}, 131*9f046c2dSEric Chamberland {7.38934376044866781e-01, 0.00000000000000000e+00, 1.00000000000000000e+00}, 132*9f046c2dSEric Chamberland {4.79574899272795285e-01, 0.00000000000000000e+00, 7.50940152310505593e-01}, 133*9f046c2dSEric Chamberland {4.77868752089733617e-01, 0.00000000000000000e+00, 1.00000000000000000e+00}, 134*9f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 7.51949651266657582e-01}, 135*9f046c2dSEric Chamberland {7.62579282058905461e-01, 7.52370482562049525e-01, 4.67395677018163536e-01}, 136*9f046c2dSEric Chamberland {1.00000000000000000e+00, 7.81083538424111712e-01, 5.32959391290483087e-01}, 137*9f046c2dSEric Chamberland {7.60498425576266124e-01, 1.00000000000000000e+00, 5.16039274773721024e-01}, 138*9f046c2dSEric Chamberland {1.00000000000000000e+00, 1.00000000000000000e+00, 5.03899302533315163e-01}, 139*9f046c2dSEric Chamberland {7.18954188589218113e-01, 7.26916038449581636e-01, 1.00000000000000000e+00}, 140*9f046c2dSEric Chamberland {4.81533470648123629e-01, 4.79286521011631217e-01, 7.15446025751505954e-01}, 141*9f046c2dSEric Chamberland {4.57888564634085005e-01, 2.26916038449581692e-01, 1.00000000000000000e+00}, 142*9f046c2dSEric Chamberland {4.95273172597062383e-01, 7.26916038449581636e-01, 1.00000000000000000e+00}, 143*9f046c2dSEric Chamberland {4.37908377178436337e-01, 4.53832076899163384e-01, 1.00000000000000000e+00}, 144*9f046c2dSEric Chamberland {1.00000000000000000e+00, 7.38666970143638135e-01, 1.00000000000000000e+00}, 145*9f046c2dSEric Chamberland {1.00000000000000000e+00, 5.19750508567749736e-01, 7.81009740023825616e-01}, 146*9f046c2dSEric Chamberland {7.38934376044866781e-01, 2.38666970143638080e-01, 1.00000000000000000e+00}, 147*9f046c2dSEric Chamberland {7.18954188589218113e-01, 4.65583008593219771e-01, 1.00000000000000000e+00}, 148*9f046c2dSEric Chamberland {1.00000000000000000e+00, 4.77333940287276159e-01, 1.00000000000000000e+00}, 149*9f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 7.17232187874490501e-01}, 150*9f046c2dSEric Chamberland {0.00000000000000000e+00, 7.33810969343145825e-01, 4.85134465279355998e-01}, 151*9f046c2dSEric Chamberland {2.60498425576266179e-01, 1.00000000000000000e+00, 4.81321811381553832e-01}, 152*9f046c2dSEric Chamberland {0.00000000000000000e+00, 1.00000000000000000e+00, 4.34464375748980947e-01}, 153*9f046c2dSEric Chamberland {0.00000000000000000e+00, 2.81085527042108152e-01, 1.00000000000000000e+00}, 154*9f046c2dSEric Chamberland {0.00000000000000000e+00, 5.14896496385254032e-01, 7.67902277404865607e-01}, 155*9f046c2dSEric Chamberland {2.76318984007844215e-01, 7.81085527042108207e-01, 1.00000000000000000e+00}, 156*9f046c2dSEric Chamberland {2.18954188589218168e-01, 5.08001565491689844e-01, 1.00000000000000000e+00}, 157*9f046c2dSEric Chamberland {0.00000000000000000e+00, 5.62171054084216304e-01, 1.00000000000000000e+00}, 158*9f046c2dSEric Chamberland {2.76318984007844215e-01, 1.00000000000000000e+00, 1.00000000000000000e+00}, 159*9f046c2dSEric Chamberland {5.36817409584110394e-01, 1.00000000000000000e+00, 7.64089623507063331e-01}, 160*9f046c2dSEric Chamberland {5.52637968015688430e-01, 1.00000000000000000e+00, 1.00000000000000000e+00}, 161*9f046c2dSEric Chamberland {5.03219805286833965e-01, 2.52370482562049525e-01, 4.66386178062011547e-01}, 162*9f046c2dSEric Chamberland {4.80554184960459430e-01, 2.39643260505815608e-01, 7.33193089031005774e-01}, 163*9f046c2dSEric Chamberland {4.81281046455856898e-01, 0.00000000000000000e+00, 5.01880304621011186e-01}, 164*9f046c2dSEric Chamberland {7.62579282058905461e-01, 5.33454020986161126e-01, 4.96455765775331570e-01}, 165*9f046c2dSEric Chamberland {2.62579282058905461e-01, 4.86181451905195350e-01, 4.83348303156371562e-01}, 166*9f046c2dSEric Chamberland {7.40766735324061787e-01, 4.99518514789690449e-01, 7.48227882887665841e-01}, 167*9f046c2dSEric Chamberland {2.40766735324061815e-01, 4.97091508698442541e-01, 7.41674151578185725e-01}, 168*9f046c2dSEric Chamberland {5.25158564117810922e-01, 5.04740965124099050e-01, 4.30892051503011964e-01}, 169*9f046c2dSEric Chamberland {7.40640523227928504e-01, 2.81083538424111656e-01, 5.31949892334331098e-01}, 170*9f046c2dSEric Chamberland {7.39787449636397643e-01, 2.59875254283874868e-01, 7.65974946167165549e-01}, 171*9f046c2dSEric Chamberland {1.00000000000000000e+00, 5.62167076848223313e-01, 5.62019480047651121e-01}, 172*9f046c2dSEric Chamberland {2.60498425576266179e-01, 7.33810969343145825e-01, 5.31991900911928939e-01}, 173*9f046c2dSEric Chamberland {2.68408704792055197e-01, 7.57448248192627016e-01, 7.65995950455964469e-01}, 174*9f046c2dSEric Chamberland {0.00000000000000000e+00, 4.67621938686291649e-01, 5.35804554809731104e-01}, 175*9f046c2dSEric Chamberland {5.23077707635171585e-01, 7.52370482562049525e-01, 4.79535649258569396e-01}, 176*9f046c2dSEric Chamberland {5.09175440116116929e-01, 7.39643260505815636e-01, 7.39767824629284698e-01}, 177*9f046c2dSEric Chamberland {5.20996851152532359e-01, 1.00000000000000000e+00, 5.28179247014126774e-01} 178*9f046c2dSEric Chamberland }; 179*9f046c2dSEric Chamberland 180*9f046c2dSEric Chamberland const PetscInt sConnectivityPrismsMesh[128][6] = { 181*9f046c2dSEric Chamberland /* rank 0 */ 182*9f046c2dSEric Chamberland {11, 7, 42, 8, 10, 9 }, 183*9f046c2dSEric Chamberland {47, 42, 43, 45, 9, 57 }, 184*9f046c2dSEric Chamberland {8, 10, 9, 77, 76, 75 }, 185*9f046c2dSEric Chamberland {45, 9, 57, 110, 75, 116}, 186*9f046c2dSEric Chamberland {17, 48, 55, 13, 14, 15 }, 187*9f046c2dSEric Chamberland {58, 55, 49, 56, 15, 52 }, 188*9f046c2dSEric Chamberland {13, 14, 15, 85, 82, 83 }, 189*9f046c2dSEric Chamberland {56, 15, 52, 118, 83, 111}, 190*9f046c2dSEric Chamberland {6, 0, 1, 2, 3, 4 }, 191*9f046c2dSEric Chamberland {54, 1, 44, 51, 4, 46 }, 192*9f046c2dSEric Chamberland {2, 3, 4, 73, 70, 71 }, 193*9f046c2dSEric Chamberland {51, 4, 46, 115, 71, 108}, 194*9f046c2dSEric Chamberland {58, 49, 43, 56, 52, 57 }, 195*9f046c2dSEric Chamberland {47, 43, 44, 45, 57, 46 }, 196*9f046c2dSEric Chamberland {56, 52, 57, 118, 111, 116}, 197*9f046c2dSEric Chamberland {45, 57, 46, 110, 116, 108}, 198*9f046c2dSEric Chamberland {77, 76, 75, 74, 31, 30 }, 199*9f046c2dSEric Chamberland {110, 75, 116, 79, 30, 117}, 200*9f046c2dSEric Chamberland {74, 31, 30, 32, 29, 78 }, 201*9f046c2dSEric Chamberland {79, 30, 117, 80, 78, 93 }, 202*9f046c2dSEric Chamberland {85, 82, 83, 81, 34, 35 }, 203*9f046c2dSEric Chamberland {118, 83, 111, 92, 35, 113}, 204*9f046c2dSEric Chamberland {81, 34, 35, 37, 86, 91 }, 205*9f046c2dSEric Chamberland {92, 35, 113, 95, 91, 94 }, 206*9f046c2dSEric Chamberland {73, 70, 71, 69, 25, 26 }, 207*9f046c2dSEric Chamberland {115, 71, 108, 87, 26, 109}, 208*9f046c2dSEric Chamberland {69, 25, 26, 28, 23, 24 }, 209*9f046c2dSEric Chamberland {87, 26, 109, 90, 24, 88 }, 210*9f046c2dSEric Chamberland {118, 111, 116, 92, 113, 117}, 211*9f046c2dSEric Chamberland {110, 116, 108, 79, 117, 109}, 212*9f046c2dSEric Chamberland {92, 113, 117, 95, 94, 93 }, 213*9f046c2dSEric Chamberland {79, 117, 109, 80, 93, 88 }, 214*9f046c2dSEric Chamberland {22, 18, 63, 19, 20, 21 }, 215*9f046c2dSEric Chamberland {68, 63, 64, 66, 21, 61 }, 216*9f046c2dSEric Chamberland {19, 20, 21, 99, 97, 98 }, 217*9f046c2dSEric Chamberland {66, 21, 61, 124, 98, 119}, 218*9f046c2dSEric Chamberland {6, 1, 59, 2, 4, 5 }, 219*9f046c2dSEric Chamberland {62, 59, 50, 60, 5, 53 }, 220*9f046c2dSEric Chamberland {2, 4, 5, 73, 71, 72 }, 221*9f046c2dSEric Chamberland {60, 5, 53, 121, 72, 112}, 222*9f046c2dSEric Chamberland {17, 12, 48, 13, 16, 14 }, 223*9f046c2dSEric Chamberland {54, 48, 65, 51, 14, 67 }, 224*9f046c2dSEric Chamberland {13, 16, 14, 85, 84, 82 }, 225*9f046c2dSEric Chamberland {51, 14, 67, 115, 82, 122}, 226*9f046c2dSEric Chamberland {62, 50, 64, 60, 53, 61 }, 227*9f046c2dSEric Chamberland {68, 64, 65, 66, 61, 67 }, 228*9f046c2dSEric Chamberland {60, 53, 61, 121, 112, 119}, 229*9f046c2dSEric Chamberland {66, 61, 67, 124, 119, 122}, 230*9f046c2dSEric Chamberland {99, 97, 98, 96, 39, 40 }, 231*9f046c2dSEric Chamberland {124, 98, 119, 106, 40, 120}, 232*9f046c2dSEric Chamberland {96, 39, 40, 41, 38, 105}, 233*9f046c2dSEric Chamberland {106, 40, 120, 107, 105, 102}, 234*9f046c2dSEric Chamberland {73, 71, 72, 69, 26, 27 }, 235*9f046c2dSEric Chamberland {121, 72, 112, 101, 27, 114}, 236*9f046c2dSEric Chamberland {69, 26, 27, 28, 24, 100}, 237*9f046c2dSEric Chamberland {101, 27, 114, 104, 100, 103}, 238*9f046c2dSEric Chamberland {85, 84, 82, 81, 36, 34 }, 239*9f046c2dSEric Chamberland {115, 82, 122, 87, 34, 123}, 240*9f046c2dSEric Chamberland {81, 36, 34, 37, 33, 86 }, 241*9f046c2dSEric Chamberland {87, 34, 123, 90, 86, 89 }, 242*9f046c2dSEric Chamberland {121, 112, 119, 101, 114, 120}, 243*9f046c2dSEric Chamberland {124, 119, 122, 106, 120, 123}, 244*9f046c2dSEric Chamberland {101, 114, 120, 104, 103, 102}, 245*9f046c2dSEric Chamberland {106, 120, 123, 107, 102, 89 }, 246*9f046c2dSEric Chamberland /* rank 1 */ 247*9f046c2dSEric Chamberland {58, 43, 7, 56, 57, 10 }, 248*9f046c2dSEric Chamberland {7, 43, 42, 10, 57, 9 }, 249*9f046c2dSEric Chamberland {56, 57, 10, 118, 116, 76 }, 250*9f046c2dSEric Chamberland {10, 57, 9, 76, 116, 75 }, 251*9f046c2dSEric Chamberland {54, 49, 48, 51, 52, 14 }, 252*9f046c2dSEric Chamberland {48, 49, 55, 14, 52, 15 }, 253*9f046c2dSEric Chamberland {51, 52, 14, 115, 111, 82 }, 254*9f046c2dSEric Chamberland {14, 52, 15, 82, 111, 83 }, 255*9f046c2dSEric Chamberland {47, 44, 0, 45, 46, 3 }, 256*9f046c2dSEric Chamberland {0, 44, 1, 3, 46, 4 }, 257*9f046c2dSEric Chamberland {45, 46, 3, 110, 108, 70 }, 258*9f046c2dSEric Chamberland {3, 46, 4, 70, 108, 71 }, 259*9f046c2dSEric Chamberland {54, 44, 49, 51, 46, 52 }, 260*9f046c2dSEric Chamberland {49, 44, 43, 52, 46, 57 }, 261*9f046c2dSEric Chamberland {51, 46, 52, 115, 108, 111}, 262*9f046c2dSEric Chamberland {52, 46, 57, 111, 108, 116}, 263*9f046c2dSEric Chamberland {118, 116, 76, 92, 117, 31 }, 264*9f046c2dSEric Chamberland {76, 116, 75, 31, 117, 30 }, 265*9f046c2dSEric Chamberland {92, 117, 31, 95, 93, 29 }, 266*9f046c2dSEric Chamberland {31, 117, 30, 29, 93, 78 }, 267*9f046c2dSEric Chamberland {115, 111, 82, 87, 113, 34 }, 268*9f046c2dSEric Chamberland {82, 111, 83, 34, 113, 35 }, 269*9f046c2dSEric Chamberland {87, 113, 34, 90, 94, 86 }, 270*9f046c2dSEric Chamberland {34, 113, 35, 86, 94, 91 }, 271*9f046c2dSEric Chamberland {110, 108, 70, 79, 109, 25 }, 272*9f046c2dSEric Chamberland {70, 108, 71, 25, 109, 26 }, 273*9f046c2dSEric Chamberland {79, 109, 25, 80, 88, 23 }, 274*9f046c2dSEric Chamberland {25, 109, 26, 23, 88, 24 }, 275*9f046c2dSEric Chamberland {115, 108, 111, 87, 109, 113}, 276*9f046c2dSEric Chamberland {111, 108, 116, 113, 109, 117}, 277*9f046c2dSEric Chamberland {87, 109, 113, 90, 88, 94 }, 278*9f046c2dSEric Chamberland {113, 109, 117, 94, 88, 93 }, 279*9f046c2dSEric Chamberland {62, 64, 18, 60, 61, 20 }, 280*9f046c2dSEric Chamberland {18, 64, 63, 20, 61, 21 }, 281*9f046c2dSEric Chamberland {60, 61, 20, 121, 119, 97 }, 282*9f046c2dSEric Chamberland {20, 61, 21, 97, 119, 98 }, 283*9f046c2dSEric Chamberland {54, 50, 1, 51, 53, 4 }, 284*9f046c2dSEric Chamberland {1, 50, 59, 4, 53, 5 }, 285*9f046c2dSEric Chamberland {51, 53, 4, 115, 112, 71 }, 286*9f046c2dSEric Chamberland {4, 53, 5, 71, 112, 72 }, 287*9f046c2dSEric Chamberland {68, 65, 12, 66, 67, 16 }, 288*9f046c2dSEric Chamberland {12, 65, 48, 16, 67, 14 }, 289*9f046c2dSEric Chamberland {66, 67, 16, 124, 122, 84 }, 290*9f046c2dSEric Chamberland {16, 67, 14, 84, 122, 82 }, 291*9f046c2dSEric Chamberland {54, 65, 50, 51, 67, 53 }, 292*9f046c2dSEric Chamberland {50, 65, 64, 53, 67, 61 }, 293*9f046c2dSEric Chamberland {51, 67, 53, 115, 122, 112}, 294*9f046c2dSEric Chamberland {53, 67, 61, 112, 122, 119}, 295*9f046c2dSEric Chamberland {121, 119, 97, 101, 120, 39 }, 296*9f046c2dSEric Chamberland {97, 119, 98, 39, 120, 40 }, 297*9f046c2dSEric Chamberland {101, 120, 39, 104, 102, 38 }, 298*9f046c2dSEric Chamberland {39, 120, 40, 38, 102, 105}, 299*9f046c2dSEric Chamberland {115, 112, 71, 87, 114, 26 }, 300*9f046c2dSEric Chamberland {71, 112, 72, 26, 114, 27 }, 301*9f046c2dSEric Chamberland {87, 114, 26, 90, 103, 24 }, 302*9f046c2dSEric Chamberland {26, 114, 27, 24, 103, 100}, 303*9f046c2dSEric Chamberland {124, 122, 84, 106, 123, 36 }, 304*9f046c2dSEric Chamberland {84, 122, 82, 36, 123, 34 }, 305*9f046c2dSEric Chamberland {106, 123, 36, 107, 89, 33 }, 306*9f046c2dSEric Chamberland {36, 123, 34, 33, 89, 86 }, 307*9f046c2dSEric Chamberland {115, 122, 112, 87, 123, 114}, 308*9f046c2dSEric Chamberland {112, 122, 119, 114, 123, 120}, 309*9f046c2dSEric Chamberland {87, 123, 114, 90, 89, 103}, 310*9f046c2dSEric Chamberland {114, 123, 120, 103, 89, 102} 311*9f046c2dSEric Chamberland }; 312*9f046c2dSEric Chamberland 313*9f046c2dSEric Chamberland /* Partitions of prisms mesh : */ 314*9f046c2dSEric Chamberland const PetscInt sInitialPartitionPrismsMesh[2][64] = { 315*9f046c2dSEric 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 }, 316*9f046c2dSEric 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, 317*9f046c2dSEric 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} 318*9f046c2dSEric Chamberland }; 319*9f046c2dSEric Chamberland 320*9f046c2dSEric Chamberland const PetscInt sNLoclCellsPrismsMesh = 64; 321*9f046c2dSEric Chamberland const PetscInt sNGlobVertsPrismsMesh = 125; 322*9f046c2dSEric Chamberland 323d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 324d71ae5a4SJacob Faibussowitsch { 325*9f046c2dSEric Chamberland PetscInt Nc = 0; 326*9f046c2dSEric Chamberland const PetscInt *InitPartForRank[2]; 327231de6a2SMatthew G. Knepley DM dm, idm, ddm; 328231de6a2SMatthew G. Knepley PetscSF sfVert, sfMig, sfPart; 329231de6a2SMatthew G. Knepley PetscPartitioner part; 330231de6a2SMatthew G. Knepley PetscSection s; 331231de6a2SMatthew G. Knepley PetscInt *cells, c; 332231de6a2SMatthew G. Knepley PetscMPIInt size, rank; 333*9f046c2dSEric Chamberland PetscBool box = PETSC_FALSE, field = PETSC_FALSE, quadsmesh = PETSC_FALSE, prismsmesh = PETSC_FALSE; 334231de6a2SMatthew G. Knepley 335327415f7SBarry Smith PetscFunctionBeginUser; 336231de6a2SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 337231de6a2SMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 338231de6a2SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 339231de6a2SMatthew G. Knepley PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a 2 processors example only"); 340231de6a2SMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-field", &field, NULL)); 341*9f046c2dSEric Chamberland PetscCall(PetscOptionsGetBool(NULL, NULL, "-box", &box, NULL)); 342*9f046c2dSEric Chamberland PetscCall(PetscOptionsGetBool(NULL, NULL, "-quadsmesh", &quadsmesh, NULL)); 343*9f046c2dSEric Chamberland PetscCall(PetscOptionsGetBool(NULL, NULL, "-prismsmesh", &prismsmesh, NULL)); 344*9f046c2dSEric Chamberland PetscCheck(1 == (box ? 1 : 0) + (quadsmesh ? 1 : 0) + (prismsmesh ? 1 : 0), PETSC_COMM_WORLD, PETSC_ERR_SUP, "Specify one and only one of -box, -quadsmesh or -prismsmesh"); 345231de6a2SMatthew G. Knepley 346231de6a2SMatthew G. Knepley PetscCall(DMPlexCreate(PETSC_COMM_WORLD, &dm)); 347231de6a2SMatthew G. Knepley if (box) { 348231de6a2SMatthew G. Knepley PetscCall(DMSetType(dm, DMPLEX)); 349231de6a2SMatthew G. Knepley PetscCall(DMSetFromOptions(dm)); 350231de6a2SMatthew G. Knepley } else { 351*9f046c2dSEric Chamberland if (quadsmesh) { 352*9f046c2dSEric Chamberland Nc = sNLoclCells2x5Mesh; //Same on each rank for this example... 353*9f046c2dSEric Chamberland PetscInt Nv = sNGlobVerts2x5Mesh; 354*9f046c2dSEric Chamberland InitPartForRank[0] = &sInitialPartition2x5Mesh[0][0]; 355*9f046c2dSEric Chamberland InitPartForRank[1] = &sInitialPartition2x5Mesh[1][0]; 356*9f046c2dSEric Chamberland const PetscInt(*Conn)[4] = sConnectivity2x5Mesh; 357*9f046c2dSEric Chamberland 358*9f046c2dSEric Chamberland const PetscInt Ncor = 4; 359*9f046c2dSEric Chamberland const PetscInt dim = 2; 360*9f046c2dSEric Chamberland 361231de6a2SMatthew G. Knepley PetscCall(PetscMalloc1(Nc * Ncor, &cells)); 362231de6a2SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 363231de6a2SMatthew G. Knepley PetscInt cell = (InitPartForRank[rank])[c], cor; 364231de6a2SMatthew G. Knepley 365ad540459SPierre Jolivet for (cor = 0; cor < Ncor; ++cor) cells[c * Ncor + cor] = Conn[cell][cor]; 366231de6a2SMatthew G. Knepley } 367231de6a2SMatthew G. Knepley PetscCall(DMSetDimension(dm, dim)); 368231de6a2SMatthew G. Knepley PetscCall(DMPlexBuildFromCellListParallel(dm, Nc, PETSC_DECIDE, Nv, Ncor, cells, &sfVert, NULL)); 369*9f046c2dSEric Chamberland } else if (prismsmesh) { 370*9f046c2dSEric Chamberland Nc = sNLoclCellsPrismsMesh; //Same on each rank for this example... 371*9f046c2dSEric Chamberland PetscInt Nv = sNGlobVertsPrismsMesh; 372*9f046c2dSEric Chamberland InitPartForRank[0] = &sInitialPartitionPrismsMesh[0][0]; 373*9f046c2dSEric Chamberland InitPartForRank[1] = &sInitialPartitionPrismsMesh[1][0]; 374*9f046c2dSEric Chamberland const PetscInt(*Conn)[6] = sConnectivityPrismsMesh; 375*9f046c2dSEric Chamberland 376*9f046c2dSEric Chamberland const PetscInt Ncor = 6; 377*9f046c2dSEric Chamberland const PetscInt dim = 3; 378*9f046c2dSEric Chamberland 379*9f046c2dSEric Chamberland PetscCall(PetscMalloc1(Nc * Ncor, &cells)); 380*9f046c2dSEric Chamberland for (c = 0; c < Nc; ++c) { 381*9f046c2dSEric Chamberland PetscInt cell = (InitPartForRank[rank])[c], cor; 382*9f046c2dSEric Chamberland 383*9f046c2dSEric Chamberland for (cor = 0; cor < Ncor; ++cor) cells[c * Ncor + cor] = Conn[cell][cor]; 384*9f046c2dSEric Chamberland } 385*9f046c2dSEric Chamberland PetscCall(DMSetDimension(dm, dim)); 386*9f046c2dSEric Chamberland PetscCall(DMPlexBuildFromCellListParallel(dm, Nc, PETSC_DECIDE, Nv, Ncor, cells, &sfVert, NULL)); 387*9f046c2dSEric Chamberland } 388231de6a2SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfVert)); 389231de6a2SMatthew G. Knepley PetscCall(PetscFree(cells)); 390231de6a2SMatthew G. Knepley PetscCall(DMPlexInterpolate(dm, &idm)); 391231de6a2SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 392231de6a2SMatthew G. Knepley dm = idm; 393231de6a2SMatthew G. Knepley } 394231de6a2SMatthew G. Knepley PetscCall(DMSetUseNatural(dm, PETSC_TRUE)); 395231de6a2SMatthew G. Knepley PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 396231de6a2SMatthew G. Knepley 397231de6a2SMatthew G. Knepley if (field) { 398231de6a2SMatthew G. Knepley const PetscInt Nf = 1; 399231de6a2SMatthew G. Knepley const PetscInt numBC = 0; 400*9f046c2dSEric Chamberland const PetscInt numComp[1] = {1}; 401*9f046c2dSEric Chamberland PetscInt numDof[4] = {0, 0, 0, 0}; 402*9f046c2dSEric Chamberland PetscInt dim; 403*9f046c2dSEric Chamberland 404*9f046c2dSEric Chamberland PetscCall(DMGetDimension(dm, &dim)); 405*9f046c2dSEric Chamberland numDof[dim] = 1; 406231de6a2SMatthew G. Knepley 407231de6a2SMatthew G. Knepley PetscCall(DMSetNumFields(dm, Nf)); 408231de6a2SMatthew G. Knepley PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &s)); 409231de6a2SMatthew G. Knepley PetscCall(DMSetLocalSection(dm, s)); 410231de6a2SMatthew G. Knepley PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD)); 411231de6a2SMatthew G. Knepley PetscCall(PetscSectionDestroy(&s)); 412231de6a2SMatthew G. Knepley } 413231de6a2SMatthew G. Knepley 414231de6a2SMatthew G. Knepley PetscCall(DMPlexGetPartitioner(dm, &part)); 415231de6a2SMatthew G. Knepley PetscCall(PetscPartitionerSetFromOptions(part)); 416231de6a2SMatthew G. Knepley 417231de6a2SMatthew G. Knepley PetscCall(DMPlexDistribute(dm, 0, &sfMig, &ddm)); 418231de6a2SMatthew G. Knepley PetscCall(PetscSFView(sfMig, PETSC_VIEWER_STDOUT_WORLD)); 419231de6a2SMatthew G. Knepley PetscCall(PetscSFCreateInverseSF(sfMig, &sfPart)); 420231de6a2SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sfPart, "Inverse Migration SF")); 421231de6a2SMatthew G. Knepley PetscCall(PetscSFView(sfPart, PETSC_VIEWER_STDOUT_WORLD)); 422231de6a2SMatthew G. Knepley 423231de6a2SMatthew G. Knepley Vec lGlobalVec, lNatVec; 424231de6a2SMatthew G. Knepley PetscScalar *lNatVecArray; 425231de6a2SMatthew G. Knepley 426231de6a2SMatthew G. Knepley { 427231de6a2SMatthew G. Knepley PetscSection s; 428231de6a2SMatthew G. Knepley 429231de6a2SMatthew G. Knepley PetscCall(DMGetGlobalSection(dm, &s)); 430231de6a2SMatthew G. Knepley PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD)); 431231de6a2SMatthew G. Knepley } 432231de6a2SMatthew G. Knepley PetscCall(DMGetGlobalVector(dm, &lNatVec)); 433231de6a2SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)lNatVec, "Natural Vector (initial partition)")); 434231de6a2SMatthew G. Knepley 435231de6a2SMatthew G. Knepley //Copying the initial partition into the "natural" vector: 436*9f046c2dSEric Chamberland PetscCall(VecZeroEntries(lNatVec)); 437231de6a2SMatthew G. Knepley PetscCall(VecGetArray(lNatVec, &lNatVecArray)); 438231de6a2SMatthew G. Knepley for (c = 0; c < Nc; ++c) lNatVecArray[c] = (InitPartForRank[rank])[c]; 439231de6a2SMatthew G. Knepley PetscCall(VecRestoreArray(lNatVec, &lNatVecArray)); 440231de6a2SMatthew G. Knepley 441231de6a2SMatthew G. Knepley PetscCall(DMGetGlobalVector(ddm, &lGlobalVec)); 442231de6a2SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)lGlobalVec, "Global Vector (reordered element numbers in the petsc distributed order)")); 443231de6a2SMatthew G. Knepley PetscCall(VecZeroEntries(lGlobalVec)); 444231de6a2SMatthew G. Knepley 445231de6a2SMatthew G. Knepley // The call to DMPlexNaturalToGlobalBegin/End does not produce our expected result... 446231de6a2SMatthew G. Knepley // In lGlobalVec, we expect to have: 447231de6a2SMatthew G. Knepley /* 448231de6a2SMatthew G. Knepley * Process [0] 449231de6a2SMatthew G. Knepley * 2. 450231de6a2SMatthew G. Knepley * 4. 451231de6a2SMatthew G. Knepley * 8. 452231de6a2SMatthew G. Knepley * 3. 453231de6a2SMatthew G. Knepley * 9. 454231de6a2SMatthew G. Knepley * Process [1] 455231de6a2SMatthew G. Knepley * 1. 456231de6a2SMatthew G. Knepley * 5. 457231de6a2SMatthew G. Knepley * 7. 458231de6a2SMatthew G. Knepley * 0. 459231de6a2SMatthew G. Knepley * 6. 460231de6a2SMatthew G. Knepley * 461231de6a2SMatthew G. Knepley * but we obtained: 462231de6a2SMatthew G. Knepley * 463231de6a2SMatthew G. Knepley * Process [0] 464231de6a2SMatthew G. Knepley * 2. 465231de6a2SMatthew G. Knepley * 4. 466231de6a2SMatthew G. Knepley * 8. 467231de6a2SMatthew G. Knepley * 0. 468231de6a2SMatthew G. Knepley * 0. 469231de6a2SMatthew G. Knepley * Process [1] 470231de6a2SMatthew G. Knepley * 0. 471231de6a2SMatthew G. Knepley * 0. 472231de6a2SMatthew G. Knepley * 0. 473231de6a2SMatthew G. Knepley * 0. 474231de6a2SMatthew G. Knepley * 0. 475231de6a2SMatthew G. Knepley */ 476231de6a2SMatthew G. Knepley 477231de6a2SMatthew G. Knepley { 478231de6a2SMatthew G. Knepley PetscSF nsf; 479231de6a2SMatthew G. Knepley 480231de6a2SMatthew G. Knepley PetscCall(DMPlexGetGlobalToNaturalSF(ddm, &nsf)); 481231de6a2SMatthew G. Knepley PetscCall(PetscSFView(nsf, NULL)); 482231de6a2SMatthew G. Knepley } 483231de6a2SMatthew G. Knepley PetscCall(DMPlexNaturalToGlobalBegin(ddm, lNatVec, lGlobalVec)); 484231de6a2SMatthew G. Knepley PetscCall(DMPlexNaturalToGlobalEnd(ddm, lNatVec, lGlobalVec)); 485231de6a2SMatthew G. Knepley 486231de6a2SMatthew G. Knepley PetscCall(VecView(lNatVec, PETSC_VIEWER_STDOUT_WORLD)); 487231de6a2SMatthew G. Knepley PetscCall(VecView(lGlobalVec, PETSC_VIEWER_STDOUT_WORLD)); 488231de6a2SMatthew G. Knepley 489231de6a2SMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dm, &lNatVec)); 490231de6a2SMatthew G. Knepley PetscCall(DMRestoreGlobalVector(ddm, &lGlobalVec)); 491231de6a2SMatthew G. Knepley 492*9f046c2dSEric Chamberland const PetscBool lUseCone = PETSC_FALSE; 493*9f046c2dSEric Chamberland const PetscBool lUseClosure = PETSC_TRUE; 494*9f046c2dSEric Chamberland PetscCall(DMSetBasicAdjacency(ddm, lUseCone, lUseClosure)); 495*9f046c2dSEric Chamberland const PetscInt lNbCellsInOverlap = 1; 496*9f046c2dSEric Chamberland PetscSF lSFMigrationOvl; 497*9f046c2dSEric Chamberland DM ddm_with_overlap; 498*9f046c2dSEric Chamberland 499*9f046c2dSEric Chamberland PetscCall(DMPlexDistributeOverlap(ddm, lNbCellsInOverlap, &lSFMigrationOvl, &ddm_with_overlap)); 500*9f046c2dSEric Chamberland 501*9f046c2dSEric Chamberland IS lISCellWithOvl = 0; 502*9f046c2dSEric Chamberland /* This is the buggy call with prisms since commit 5ae96e2b862 */ 503*9f046c2dSEric Chamberland PetscCall(DMPlexGetCellNumbering(ddm_with_overlap, &lISCellWithOvl)); 504*9f046c2dSEric Chamberland /* Here, we can see the elements in the overlap within the IS: they are the ones with negative indices */ 505*9f046c2dSEric Chamberland PetscCall(ISView(lISCellWithOvl, PETSC_VIEWER_STDOUT_WORLD)); 506*9f046c2dSEric Chamberland 507*9f046c2dSEric Chamberland PetscCall(PetscSFDestroy(&lSFMigrationOvl)); 508*9f046c2dSEric Chamberland PetscCall(DMDestroy(&ddm_with_overlap)); 509231de6a2SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfMig)); 510231de6a2SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfPart)); 511231de6a2SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 512231de6a2SMatthew G. Knepley PetscCall(DMDestroy(&ddm)); 513231de6a2SMatthew G. Knepley PetscCall(PetscFinalize()); 514231de6a2SMatthew G. Knepley return 0; 515231de6a2SMatthew G. Knepley } 516231de6a2SMatthew G. Knepley 517231de6a2SMatthew G. Knepley /*TEST 518231de6a2SMatthew G. Knepley 519231de6a2SMatthew G. Knepley testset: 520231de6a2SMatthew G. Knepley args: -field -petscpartitioner_type simple 521231de6a2SMatthew G. Knepley nsize: 2 522231de6a2SMatthew G. Knepley 523231de6a2SMatthew G. Knepley test: 524231de6a2SMatthew G. Knepley suffix: 0 525*9f046c2dSEric Chamberland args: -quadsmesh 526*9f046c2dSEric Chamberland output_file: output/ex47_0.out 527231de6a2SMatthew G. Knepley 528231de6a2SMatthew G. Knepley test: 529231de6a2SMatthew G. Knepley suffix: 1 530231de6a2SMatthew G. Knepley args: -box -dm_plex_simplex 0 -dm_plex_box_faces 2,5 -dm_distribute 531*9f046c2dSEric Chamberland output_file: output/ex47_1.out 532*9f046c2dSEric Chamberland 533*9f046c2dSEric Chamberland test: 534*9f046c2dSEric Chamberland suffix: 2 535*9f046c2dSEric Chamberland args: -prismsmesh 536*9f046c2dSEric Chamberland output_file: output/ex47_2.out 537231de6a2SMatthew G. Knepley 538231de6a2SMatthew G. Knepley TEST*/ 539