xref: /petsc/src/dm/impls/stag/tests/ex2.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1 static char help[] = "Directly check a DMStag local-to-global scatter";
2 
3 #include <petscdm.h>
4 #include <petscdmstag.h>
5 
6 static PetscErrorCode Test_3d_4x4x4_3x3x3(DM dmstag);
7 
main(int argc,char ** argv)8 int main(int argc, char **argv)
9 {
10   DM       dmstag;
11   PetscInt dim, dof[4], i, elx, ely, elz;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
15   dim = 3;
16   for (i = 0; i < 4; ++i) dof[i] = 1;
17   elx = ely = elz = 4;
18   switch (dim) {
19   case 3:
20     PetscCall(DMStagCreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, elx, ely, elz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof[0], dof[1], dof[2], dof[3], DMSTAG_STENCIL_BOX, 1, NULL, NULL, NULL, &dmstag));
21     break;
22   default:
23     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "No support for dimension %" PetscInt_FMT, dim);
24   }
25   PetscCall(DMSetFromOptions(dmstag));
26   PetscCall(DMSetUp(dmstag));
27   PetscCall(Test_3d_4x4x4_3x3x3(dmstag));
28   PetscCall(DMDestroy(&dmstag));
29   PetscCall(PetscFinalize());
30   return 0;
31 }
32 
Test_3d_4x4x4_3x3x3(DM dmstag)33 static PetscErrorCode Test_3d_4x4x4_3x3x3(DM dmstag)
34 {
35   Vec          vecLocal, vecGlobal;
36   PetscInt     i, low, high, n;
37   PetscScalar *arr;
38 
39   PetscFunctionBeginUser;
40   /* Check that grid and rank grid is as expected for this test */
41   {
42     PetscInt nRanks[3], n[3], dim;
43     PetscCall(DMGetDimension(dmstag, &dim));
44     PetscCheck(dim == 3, PetscObjectComm((PetscObject)dmstag), PETSC_ERR_SUP, "This is a 3d test");
45     PetscCall(DMStagGetNumRanks(dmstag, &nRanks[0], &nRanks[1], &nRanks[2]));
46     for (i = 0; i < 3; ++i) PetscCheck(nRanks[i] == 3, PetscObjectComm((PetscObject)dmstag), PETSC_ERR_SUP, "This test requires a 3x3x3 rank grid (run on 27 ranks)");
47     PetscCall(DMStagGetGlobalSizes(dmstag, &n[0], &n[1], &n[2]));
48     for (i = 0; i < 3; ++i) PetscCheck(n[i] == 4, PetscObjectComm((PetscObject)dmstag), PETSC_ERR_SUP, "This test requires a 4x4x4 element grid");
49   }
50 
51   /* Populate global vector by converting the global index number to a scalar value. */
52   PetscCall(DMCreateGlobalVector(dmstag, &vecGlobal));
53   PetscCall(VecGetOwnershipRange(vecGlobal, &low, &high));
54   n = high - low;
55   PetscCall(VecGetArray(vecGlobal, &arr));
56   for (i = 0; i < n; ++i) arr[i] = (PetscScalar)(i + low);
57   PetscCall(VecRestoreArray(vecGlobal, &arr));
58 
59   /* Populate a local vector initially with -1, then glocal->local scatter */
60   PetscCall(DMCreateLocalVector(dmstag, &vecLocal));
61   PetscCall(VecSet(vecLocal, -1.0));
62   PetscCall(DMGlobalToLocalBegin(dmstag, vecGlobal, INSERT_VALUES, vecLocal));
63   PetscCall(DMGlobalToLocalEnd(dmstag, vecGlobal, INSERT_VALUES, vecLocal));
64 
65   /* Check that entries are as expected */
66   {
67     PetscInt           entriesGhost, nerr;
68     const PetscInt     maxErrPerRank = 3;
69     PetscScalar       *arrLocalExpected;
70     const PetscScalar *arrLocal;
71     PetscMPIInt        rank;
72 
73     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dmstag), &rank));
74     PetscCall(VecGetSize(vecLocal, &entriesGhost)); /* entriesGhost happens to always be 216 here */
75     PetscCall(PetscMalloc1(entriesGhost, &arrLocalExpected));
76 
77     /* Hand-computed expected entries (27 blocks of 8 in all cases) */
78     if (rank == 0) {
79       const PetscScalar arrLocalExpectedHere[] = {
80         0,   1,   2,   3,   4,   5,   6,   7,   8,   9,   10,  11,  12,  13,  14,  15,  64,  65,  66,  67,  68,  69,  70,  71,  16,  17,  18,  19,  20,  21,  22,  23,  24,  25,  26,  27,  28,  29,  30,  31,  72,  73,  74,  75,
81         76,  77,  78,  79,  144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 176, 177, 178, 179, 180, 181, 182, 183, 32,  33,  34,  35,  36,  37,  38,  39,  40,  41,  42,  43,  44,  45,  46,  47,
82         80,  81,  82,  83,  84,  85,  86,  87,  48,  49,  50,  51,  52,  53,  54,  55,  56,  57,  58,  59,  60,  61,  62,  63,  88,  89,  90,  91,  92,  93,  94,  95,  160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171,
83         172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339, 356, 357, 358, 359, 360, 361, 362, 363, 340, 341, 342, 343, 344, 345, 346, 347,
84         348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419,
85       };
86       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
87     } else if (rank == 1) {
88       const PetscScalar arrLocalExpectedHere[] = {
89         8,   9,   10,  11,  12,  13,  14,  15,  64,  65,  66,  67,  68,  69,  70,  71,  96,  97,  98,  99,  100, 101, 102, 103, 24,  25,  26,  27,  28,  29,  30,  31,  72,  73,  74,  75,  76,  77,  78,  79,  108, 109, 110, 111,
90         112, 113, 114, 115, 152, 153, 154, 155, 156, 157, 158, 159, 176, 177, 178, 179, 180, 181, 182, 183, 192, 193, 194, 195, 196, 197, 198, 199, 40,  41,  42,  43,  44,  45,  46,  47,  80,  81,  82,  83,  84,  85,  86,  87,
91         120, 121, 122, 123, 124, 125, 126, 127, 56,  57,  58,  59,  60,  61,  62,  63,  88,  89,  90,  91,  92,  93,  94,  95,  132, 133, 134, 135, 136, 137, 138, 139, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187,
92         188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 332, 333, 334, 335, 336, 337, 338, 339, 356, 357, 358, 359, 360, 361, 362, 363, 372, 373, 374, 375, 376, 377, 378, 379, 348, 349, 350, 351, 352, 353, 354, 355,
93         364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427,
94       };
95       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
96     } else if (rank == 2) {
97       const PetscScalar arrLocalExpectedHere[] = {
98         64,  65,  66,  67,  68,  69,  70,  71,  96,  97,  98,  99,  100, 101, 102, 103, 104, -1,  105, -1,  106, -1,  107, -1,  72,  73,  74,  75,  76,  77,  78,  79,  108, 109, 110, 111, 112, 113, 114, 115, 116, -1,  117, -1,
99         118, -1,  119, -1,  176, 177, 178, 179, 180, 181, 182, 183, 192, 193, 194, 195, 196, 197, 198, 199, 200, -1,  201, -1,  202, -1,  203, -1,  80,  81,  82,  83,  84,  85,  86,  87,  120, 121, 122, 123, 124, 125, 126, 127,
100         128, -1,  129, -1,  130, -1,  131, -1,  88,  89,  90,  91,  92,  93,  94,  95,  132, 133, 134, 135, 136, 137, 138, 139, 140, -1,  141, -1,  142, -1,  143, -1,  184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207,
101         208, 209, 210, 211, 212, -1,  213, -1,  214, -1,  215, -1,  356, 357, 358, 359, 360, 361, 362, 363, 372, 373, 374, 375, 376, 377, 378, 379, 380, -1,  381, -1,  382, -1,  383, -1,  364, 365, 366, 367, 368, 369, 370, 371,
102         384, 385, 386, 387, 388, 389, 390, 391, 392, -1,  393, -1,  394, -1,  395, -1,  412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 428, -1,  429, -1,  430, -1,  431, -1,
103       };
104       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
105     } else if (rank == 3) {
106       const PetscScalar arrLocalExpectedHere[] = {
107         16,  17,  18,  19,  20,  21,  22,  23,  24,  25,  26,  27,  28,  29,  30,  31,  72,  73,  74,  75,  76,  77,  78,  79,  144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 176, 177, 178, 179,
108         180, 181, 182, 183, 216, 217, 218, 219, 220, 221, 222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 264, 265, 266, 267, 268, 269, 270, 271, 48,  49,  50,  51,  52,  53,  54,  55,  56,  57,  58,  59,  60,  61,  62,  63,
109         88,  89,  90,  91,  92,  93,  94,  95,  160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251,
110         252, 253, 254, 255, 276, 277, 278, 279, 280, 281, 282, 283, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 396, 397, 398, 399, 400, 401, 402, 403,
111         404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 432, 433, 434, 435, 436, 437, 438, 439, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463,
112       };
113       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
114     } else if (rank == 4) {
115       const PetscScalar arrLocalExpectedHere[] = {
116         24,  25,  26,  27,  28,  29,  30,  31,  72,  73,  74,  75,  76,  77,  78,  79,  108, 109, 110, 111, 112, 113, 114, 115, 152, 153, 154, 155, 156, 157, 158, 159, 176, 177, 178, 179, 180, 181, 182, 183, 192, 193, 194, 195,
117         196, 197, 198, 199, 224, 225, 226, 227, 228, 229, 230, 231, 264, 265, 266, 267, 268, 269, 270, 271, 288, 289, 290, 291, 292, 293, 294, 295, 56,  57,  58,  59,  60,  61,  62,  63,  88,  89,  90,  91,  92,  93,  94,  95,
118         132, 133, 134, 135, 136, 137, 138, 139, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 248, 249, 250, 251, 252, 253, 254, 255, 276, 277, 278, 279,
119         280, 281, 282, 283, 306, 307, 308, 309, 310, 311, 312, 313, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 404, 405, 406, 407, 408, 409, 410, 411,
120         412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475,
121       };
122       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
123     } else if (rank == 5) {
124       const PetscScalar arrLocalExpectedHere[] = {
125         72,  73,  74,  75,  76,  77,  78,  79,  108, 109, 110, 111, 112, 113, 114, 115, 116, -1,  117, -1,  118, -1,  119, -1,  176, 177, 178, 179, 180, 181, 182, 183, 192, 193, 194, 195, 196, 197, 198, 199, 200, -1,  201, -1,
126         202, -1,  203, -1,  264, 265, 266, 267, 268, 269, 270, 271, 288, 289, 290, 291, 292, 293, 294, 295, 296, -1,  297, -1,  298, -1,  299, -1,  88,  89,  90,  91,  92,  93,  94,  95,  132, 133, 134, 135, 136, 137, 138, 139,
127         140, -1,  141, -1,  142, -1,  143, -1,  184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 212, -1,  213, -1,  214, -1,  215, -1,  276, 277, 278, 279, 280, 281, 282, 283, 306, 307, 308, 309,
128         310, 311, 312, 313, 314, -1,  315, -1,  316, -1,  317, -1,  364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 392, -1,  393, -1,  394, -1,  395, -1,  412, 413, 414, 415, 416, 417, 418, 419,
129         420, 421, 422, 423, 424, 425, 426, 427, 428, -1,  429, -1,  430, -1,  431, -1,  456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 476, -1,  477, -1,  478, -1,  479, -1,
130       };
131       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
132     } else if (rank == 6) {
133       const PetscScalar arrLocalExpectedHere[] = {
134         144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 176, 177, 178, 179, 180, 181, 182, 183, 216, 217, 218, 219, 220, 221, 222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 264, 265, 266, 267,
135         268, 269, 270, 271, 232, 233, -1,  -1,  234, 235, -1,  -1,  236, 237, -1,  -1,  238, 239, -1,  -1,  272, 273, -1,  -1,  274, 275, -1,  -1,  160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175,
136         184, 185, 186, 187, 188, 189, 190, 191, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 276, 277, 278, 279, 280, 281, 282, 283, 256, 257, -1,  -1,  258, 259, -1,  -1,  260, 261, -1,  -1,
137         262, 263, -1,  -1,  284, 285, -1,  -1,  286, 287, -1,  -1,  396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 432, 433, 434, 435, 436, 437, 438, 439,
138         440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 448, 449, -1,  -1,  450, 451, -1,  -1,  452, 453, -1,  -1,  454, 455, -1,  -1,  464, 465, -1,  -1,  466, 467, -1,  -1,
139       };
140       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
141     } else if (rank == 7) {
142       const PetscScalar arrLocalExpectedHere[] = {
143         152, 153, 154, 155, 156, 157, 158, 159, 176, 177, 178, 179, 180, 181, 182, 183, 192, 193, 194, 195, 196, 197, 198, 199, 224, 225, 226, 227, 228, 229, 230, 231, 264, 265, 266, 267, 268, 269, 270, 271, 288, 289, 290, 291,
144         292, 293, 294, 295, 236, 237, -1,  -1,  238, 239, -1,  -1,  272, 273, -1,  -1,  274, 275, -1,  -1,  300, 301, -1,  -1,  302, 303, -1,  -1,  168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191,
145         204, 205, 206, 207, 208, 209, 210, 211, 248, 249, 250, 251, 252, 253, 254, 255, 276, 277, 278, 279, 280, 281, 282, 283, 306, 307, 308, 309, 310, 311, 312, 313, 260, 261, -1,  -1,  262, 263, -1,  -1,  284, 285, -1,  -1,
146         286, 287, -1,  -1,  318, 319, -1,  -1,  320, 321, -1,  -1,  404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 440, 441, 442, 443, 444, 445, 446, 447,
147         456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 452, 453, -1,  -1,  454, 455, -1,  -1,  464, 465, -1,  -1,  466, 467, -1,  -1,  480, 481, -1,  -1,  482, 483, -1,  -1,
148       };
149       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
150     } else if (rank == 8) {
151       const PetscScalar arrLocalExpectedHere[] = {
152         176, 177, 178, 179, 180, 181, 182, 183, 192, 193, 194, 195, 196, 197, 198, 199, 200, -1,  201, -1,  202, -1,  203, -1,  264, 265, 266, 267, 268, 269, 270, 271, 288, 289, 290, 291, 292, 293, 294, 295, 296, -1,  297, -1,
153         298, -1,  299, -1,  272, 273, -1,  -1,  274, 275, -1,  -1,  300, 301, -1,  -1,  302, 303, -1,  -1,  304, -1,  -1,  -1,  305, -1,  -1,  -1,  184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211,
154         212, -1,  213, -1,  214, -1,  215, -1,  276, 277, 278, 279, 280, 281, 282, 283, 306, 307, 308, 309, 310, 311, 312, 313, 314, -1,  315, -1,  316, -1,  317, -1,  284, 285, -1,  -1,  286, 287, -1,  -1,  318, 319, -1,  -1,
155         320, 321, -1,  -1,  322, -1,  -1,  -1,  323, -1,  -1,  -1,  412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 428, -1,  429, -1,  430, -1,  431, -1,  456, 457, 458, 459, 460, 461, 462, 463,
156         468, 469, 470, 471, 472, 473, 474, 475, 476, -1,  477, -1,  478, -1,  479, -1,  464, 465, -1,  -1,  466, 467, -1,  -1,  480, 481, -1,  -1,  482, 483, -1,  -1,  484, -1,  -1,  -1,  485, -1,  -1,  -1,
157       };
158       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
159     } else if (rank == 9) {
160       const PetscScalar arrLocalExpectedHere[] = {
161         32,  33,  34,  35,  36,  37,  38,  39,  40,  41,  42,  43,  44,  45,  46,  47,  80,  81,  82,  83,  84,  85,  86,  87,  48,  49,  50,  51,  52,  53,  54,  55,  56,  57,  58,  59,  60,  61,  62,  63,  88,  89,  90,  91,
162         92,  93,  94,  95,  160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339,
163         356, 357, 358, 359, 360, 361, 362, 363, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407,
164         408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 486, 487, 488, 489, 490, 491, 492, 493, 494, 495, 496, 497, 498, 499, 500, 501, 534, 535, 536, 537, 538, 539, 540, 541, 502, 503, 504, 505, 506, 507, 508, 509,
165         510, 511, 512, 513, 514, 515, 516, 517, 542, 543, 544, 545, 546, 547, 548, 549, 594, 595, 596, 597, 598, 599, 600, 601, 602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625,
166       };
167       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
168     } else if (rank == 10) {
169       const PetscScalar arrLocalExpectedHere[] = {
170         40,  41,  42,  43,  44,  45,  46,  47,  80,  81,  82,  83,  84,  85,  86,  87,  120, 121, 122, 123, 124, 125, 126, 127, 56,  57,  58,  59,  60,  61,  62,  63,  88,  89,  90,  91,  92,  93,  94,  95,  132, 133, 134, 135,
171         136, 137, 138, 139, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 332, 333, 334, 335, 336, 337, 338, 339, 356, 357, 358, 359, 360, 361, 362, 363,
172         372, 373, 374, 375, 376, 377, 378, 379, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415,
173         416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 494, 495, 496, 497, 498, 499, 500, 501, 534, 535, 536, 537, 538, 539, 540, 541, 558, 559, 560, 561, 562, 563, 564, 565, 510, 511, 512, 513, 514, 515, 516, 517,
174         542, 543, 544, 545, 546, 547, 548, 549, 570, 571, 572, 573, 574, 575, 576, 577, 602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637,
175       };
176       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
177     } else if (rank == 11) {
178       const PetscScalar arrLocalExpectedHere[] = {
179         80,  81,  82,  83,  84,  85,  86,  87,  120, 121, 122, 123, 124, 125, 126, 127, 128, -1,  129, -1,  130, -1,  131, -1,  88,  89,  90,  91,  92,  93,  94,  95,  132, 133, 134, 135, 136, 137, 138, 139, 140, -1,  141, -1,
180         142, -1,  143, -1,  184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 212, -1,  213, -1,  214, -1,  215, -1,  356, 357, 358, 359, 360, 361, 362, 363, 372, 373, 374, 375, 376, 377, 378, 379,
181         380, -1,  381, -1,  382, -1,  383, -1,  364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 392, -1,  393, -1,  394, -1,  395, -1,  412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423,
182         424, 425, 426, 427, 428, -1,  429, -1,  430, -1,  431, -1,  534, 535, 536, 537, 538, 539, 540, 541, 558, 559, 560, 561, 562, 563, 564, 565, 566, -1,  567, -1,  568, -1,  569, -1,  542, 543, 544, 545, 546, 547, 548, 549,
183         570, 571, 572, 573, 574, 575, 576, 577, 578, -1,  579, -1,  580, -1,  581, -1,  618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637, 638, -1,  639, -1,  640, -1,  641, -1,
184       };
185       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
186     } else if (rank == 12) {
187       const PetscScalar arrLocalExpectedHere[] = {
188         48,  49,  50,  51,  52,  53,  54,  55,  56,  57,  58,  59,  60,  61,  62,  63,  88,  89,  90,  91,  92,  93,  94,  95,  160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187,
189         188, 189, 190, 191, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 276, 277, 278, 279, 280, 281, 282, 283, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355,
190         364, 365, 366, 367, 368, 369, 370, 371, 396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 432, 433, 434, 435, 436, 437, 438, 439, 440, 441, 442, 443,
191         444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 502, 503, 504, 505, 506, 507, 508, 509, 510, 511, 512, 513, 514, 515, 516, 517, 542, 543, 544, 545, 546, 547, 548, 549, 594, 595, 596, 597, 598, 599, 600, 601,
192         602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625, 648, 649, 650, 651, 652, 653, 654, 655, 656, 657, 658, 659, 660, 661, 662, 663, 684, 685, 686, 687, 688, 689, 690, 691,
193       };
194       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
195     } else if (rank == 13) {
196       const PetscScalar arrLocalExpectedHere[] = {
197         56,  57,  58,  59,  60,  61,  62,  63,  88,  89,  90,  91,  92,  93,  94,  95,  132, 133, 134, 135, 136, 137, 138, 139, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207,
198         208, 209, 210, 211, 248, 249, 250, 251, 252, 253, 254, 255, 276, 277, 278, 279, 280, 281, 282, 283, 306, 307, 308, 309, 310, 311, 312, 313, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371,
199         384, 385, 386, 387, 388, 389, 390, 391, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459,
200         460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 510, 511, 512, 513, 514, 515, 516, 517, 542, 543, 544, 545, 546, 547, 548, 549, 570, 571, 572, 573, 574, 575, 576, 577, 602, 603, 604, 605, 606, 607, 608, 609,
201         618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637, 656, 657, 658, 659, 660, 661, 662, 663, 684, 685, 686, 687, 688, 689, 690, 691, 702, 703, 704, 705, 706, 707, 708, 709,
202       };
203       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
204     } else if (rank == 14) {
205       const PetscScalar arrLocalExpectedHere[] = {
206         88,  89,  90,  91,  92,  93,  94,  95,  132, 133, 134, 135, 136, 137, 138, 139, 140, -1,  141, -1,  142, -1,  143, -1,  184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 212, -1,  213, -1,
207         214, -1,  215, -1,  276, 277, 278, 279, 280, 281, 282, 283, 306, 307, 308, 309, 310, 311, 312, 313, 314, -1,  315, -1,  316, -1,  317, -1,  364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391,
208         392, -1,  393, -1,  394, -1,  395, -1,  412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 428, -1,  429, -1,  430, -1,  431, -1,  456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471,
209         472, 473, 474, 475, 476, -1,  477, -1,  478, -1,  479, -1,  542, 543, 544, 545, 546, 547, 548, 549, 570, 571, 572, 573, 574, 575, 576, 577, 578, -1,  579, -1,  580, -1,  581, -1,  618, 619, 620, 621, 622, 623, 624, 625,
210         630, 631, 632, 633, 634, 635, 636, 637, 638, -1,  639, -1,  640, -1,  641, -1,  684, 685, 686, 687, 688, 689, 690, 691, 702, 703, 704, 705, 706, 707, 708, 709, 710, -1,  711, -1,  712, -1,  713, -1,
211       };
212       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
213     } else if (rank == 15) {
214       const PetscScalar arrLocalExpectedHere[] = {
215         160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 276, 277, 278, 279,
216         280, 281, 282, 283, 256, 257, -1,  -1,  258, 259, -1,  -1,  260, 261, -1,  -1,  262, 263, -1,  -1,  284, 285, -1,  -1,  286, 287, -1,  -1,  396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407, 408, 409, 410, 411,
217         412, 413, 414, 415, 416, 417, 418, 419, 432, 433, 434, 435, 436, 437, 438, 439, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 448, 449, -1,  -1,  450, 451, -1,  -1,  452, 453, -1,  -1,
218         454, 455, -1,  -1,  464, 465, -1,  -1,  466, 467, -1,  -1,  594, 595, 596, 597, 598, 599, 600, 601, 602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625, 648, 649, 650, 651, 652, 653, 654, 655,
219         656, 657, 658, 659, 660, 661, 662, 663, 684, 685, 686, 687, 688, 689, 690, 691, 664, 665, -1,  -1,  666, 667, -1,  -1,  668, 669, -1,  -1,  670, 671, -1,  -1,  692, 693, -1,  -1,  694, 695, -1,  -1,
220       };
221       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
222     } else if (rank == 16) {
223       const PetscScalar arrLocalExpectedHere[] = {
224         168, 169, 170, 171, 172, 173, 174, 175, 184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 248, 249, 250, 251, 252, 253, 254, 255, 276, 277, 278, 279, 280, 281, 282, 283, 306, 307, 308, 309,
225         310, 311, 312, 313, 260, 261, -1,  -1,  262, 263, -1,  -1,  284, 285, -1,  -1,  286, 287, -1,  -1,  318, 319, -1,  -1,  320, 321, -1,  -1,  404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419,
226         420, 421, 422, 423, 424, 425, 426, 427, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 452, 453, -1,  -1,  454, 455, -1,  -1,  464, 465, -1,  -1,
227         466, 467, -1,  -1,  480, 481, -1,  -1,  482, 483, -1,  -1,  602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637, 656, 657, 658, 659, 660, 661, 662, 663,
228         684, 685, 686, 687, 688, 689, 690, 691, 702, 703, 704, 705, 706, 707, 708, 709, 668, 669, -1,  -1,  670, 671, -1,  -1,  692, 693, -1,  -1,  694, 695, -1,  -1,  714, 715, -1,  -1,  716, 717, -1,  -1,
229       };
230       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
231     } else if (rank == 17) {
232       const PetscScalar arrLocalExpectedHere[] = {
233         184, 185, 186, 187, 188, 189, 190, 191, 204, 205, 206, 207, 208, 209, 210, 211, 212, -1,  213, -1,  214, -1,  215, -1,  276, 277, 278, 279, 280, 281, 282, 283, 306, 307, 308, 309, 310, 311, 312, 313, 314, -1,  315, -1,
234         316, -1,  317, -1,  284, 285, -1,  -1,  286, 287, -1,  -1,  318, 319, -1,  -1,  320, 321, -1,  -1,  322, -1,  -1,  -1,  323, -1,  -1,  -1,  412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427,
235         428, -1,  429, -1,  430, -1,  431, -1,  456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 476, -1,  477, -1,  478, -1,  479, -1,  464, 465, -1,  -1,  466, 467, -1,  -1,  480, 481, -1,  -1,
236         482, 483, -1,  -1,  484, -1,  -1,  -1,  485, -1,  -1,  -1,  618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637, 638, -1,  639, -1,  640, -1,  641, -1,  684, 685, 686, 687, 688, 689, 690, 691,
237         702, 703, 704, 705, 706, 707, 708, 709, 710, -1,  711, -1,  712, -1,  713, -1,  692, 693, -1,  -1,  694, 695, -1,  -1,  714, 715, -1,  -1,  716, 717, -1,  -1,  718, -1,  -1,  -1,  719, -1,  -1,  -1,
238       };
239       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
240     } else if (rank == 18) {
241       const PetscScalar arrLocalExpectedHere[] = {
242         324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339, 356, 357, 358, 359, 360, 361, 362, 363, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367,
243         368, 369, 370, 371, 396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 486, 487, 488, 489, 490, 491, 492, 493, 494, 495, 496, 497, 498, 499, 500, 501,
244         534, 535, 536, 537, 538, 539, 540, 541, 502, 503, 504, 505, 506, 507, 508, 509, 510, 511, 512, 513, 514, 515, 516, 517, 542, 543, 544, 545, 546, 547, 548, 549, 594, 595, 596, 597, 598, 599, 600, 601, 602, 603, 604, 605,
245         606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625, 518, 519, 520, 521, -1,  -1,  -1,  -1,  522, 523, 524, 525, -1,  -1,  -1,  -1,  550, 551, 552, 553, -1,  -1,  -1,  -1,  526, 527, 528, 529, -1,  -1,  -1,  -1,
246         530, 531, 532, 533, -1,  -1,  -1,  -1,  554, 555, 556, 557, -1,  -1,  -1,  -1,  610, 611, 612, 613, -1,  -1,  -1,  -1,  614, 615, 616, 617, -1,  -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,
247       };
248       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
249     } else if (rank == 19) {
250       const PetscScalar arrLocalExpectedHere[] = {
251         332, 333, 334, 335, 336, 337, 338, 339, 356, 357, 358, 359, 360, 361, 362, 363, 372, 373, 374, 375, 376, 377, 378, 379, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387,
252         388, 389, 390, 391, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 494, 495, 496, 497, 498, 499, 500, 501, 534, 535, 536, 537, 538, 539, 540, 541,
253         558, 559, 560, 561, 562, 563, 564, 565, 510, 511, 512, 513, 514, 515, 516, 517, 542, 543, 544, 545, 546, 547, 548, 549, 570, 571, 572, 573, 574, 575, 576, 577, 602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621,
254         622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637, 522, 523, 524, 525, -1,  -1,  -1,  -1,  550, 551, 552, 553, -1,  -1,  -1,  -1,  582, 583, 584, 585, -1,  -1,  -1,  -1,  530, 531, 532, 533, -1,  -1,  -1,  -1,
255         554, 555, 556, 557, -1,  -1,  -1,  -1,  588, 589, 590, 591, -1,  -1,  -1,  -1,  614, 615, 616, 617, -1,  -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,  642, 643, 644, 645, -1,  -1,  -1,  -1,
256       };
257       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
258     } else if (rank == 20) {
259       const PetscScalar arrLocalExpectedHere[] = {
260         356, 357, 358, 359, 360, 361, 362, 363, 372, 373, 374, 375, 376, 377, 378, 379, 380, -1,  381, -1,  382, -1,  383, -1,  364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 392, -1,  393, -1,
261         394, -1,  395, -1,  412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 428, -1,  429, -1,  430, -1,  431, -1,  534, 535, 536, 537, 538, 539, 540, 541, 558, 559, 560, 561, 562, 563, 564, 565,
262         566, -1,  567, -1,  568, -1,  569, -1,  542, 543, 544, 545, 546, 547, 548, 549, 570, 571, 572, 573, 574, 575, 576, 577, 578, -1,  579, -1,  580, -1,  581, -1,  618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633,
263         634, 635, 636, 637, 638, -1,  639, -1,  640, -1,  641, -1,  550, 551, 552, 553, -1,  -1,  -1,  -1,  582, 583, 584, 585, -1,  -1,  -1,  -1,  586, -1,  587, -1,  -1,  -1,  -1,  -1,  554, 555, 556, 557, -1,  -1,  -1,  -1,
264         588, 589, 590, 591, -1,  -1,  -1,  -1,  592, -1,  593, -1,  -1,  -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,  642, 643, 644, 645, -1,  -1,  -1,  -1,  646, -1,  647, -1,  -1,  -1,  -1,  -1,
265       };
266       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
267     } else if (rank == 21) {
268       const PetscScalar arrLocalExpectedHere[] = {
269         340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415,
270         416, 417, 418, 419, 432, 433, 434, 435, 436, 437, 438, 439, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 502, 503, 504, 505, 506, 507, 508, 509, 510, 511, 512, 513, 514, 515, 516, 517,
271         542, 543, 544, 545, 546, 547, 548, 549, 594, 595, 596, 597, 598, 599, 600, 601, 602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625, 648, 649, 650, 651, 652, 653, 654, 655, 656, 657, 658, 659,
272         660, 661, 662, 663, 684, 685, 686, 687, 688, 689, 690, 691, 526, 527, 528, 529, -1,  -1,  -1,  -1,  530, 531, 532, 533, -1,  -1,  -1,  -1,  554, 555, 556, 557, -1,  -1,  -1,  -1,  610, 611, 612, 613, -1,  -1,  -1,  -1,
273         614, 615, 616, 617, -1,  -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,  672, 673, 674, 675, -1,  -1,  -1,  -1,  676, 677, 678, 679, -1,  -1,  -1,  -1,  696, 697, 698, 699, -1,  -1,  -1,  -1,
274       };
275       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
276     } else if (rank == 22) {
277       const PetscScalar arrLocalExpectedHere[] = {
278         348, 349, 350, 351, 352, 353, 354, 355, 364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423,
279         424, 425, 426, 427, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 510, 511, 512, 513, 514, 515, 516, 517, 542, 543, 544, 545, 546, 547, 548, 549,
280         570, 571, 572, 573, 574, 575, 576, 577, 602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637, 656, 657, 658, 659, 660, 661, 662, 663, 684, 685, 686, 687,
281         688, 689, 690, 691, 702, 703, 704, 705, 706, 707, 708, 709, 530, 531, 532, 533, -1,  -1,  -1,  -1,  554, 555, 556, 557, -1,  -1,  -1,  -1,  588, 589, 590, 591, -1,  -1,  -1,  -1,  614, 615, 616, 617, -1,  -1,  -1,  -1,
282         626, 627, 628, 629, -1,  -1,  -1,  -1,  642, 643, 644, 645, -1,  -1,  -1,  -1,  676, 677, 678, 679, -1,  -1,  -1,  -1,  696, 697, 698, 699, -1,  -1,  -1,  -1,  720, 721, 722, 723, -1,  -1,  -1,  -1,
283       };
284       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
285     } else if (rank == 23) {
286       const PetscScalar arrLocalExpectedHere[] = {
287         364, 365, 366, 367, 368, 369, 370, 371, 384, 385, 386, 387, 388, 389, 390, 391, 392, -1,  393, -1,  394, -1,  395, -1,  412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 428, -1,  429, -1,
288         430, -1,  431, -1,  456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 476, -1,  477, -1,  478, -1,  479, -1,  542, 543, 544, 545, 546, 547, 548, 549, 570, 571, 572, 573, 574, 575, 576, 577,
289         578, -1,  579, -1,  580, -1,  581, -1,  618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637, 638, -1,  639, -1,  640, -1,  641, -1,  684, 685, 686, 687, 688, 689, 690, 691, 702, 703, 704, 705,
290         706, 707, 708, 709, 710, -1,  711, -1,  712, -1,  713, -1,  554, 555, 556, 557, -1,  -1,  -1,  -1,  588, 589, 590, 591, -1,  -1,  -1,  -1,  592, -1,  593, -1,  -1,  -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,
291         642, 643, 644, 645, -1,  -1,  -1,  -1,  646, -1,  647, -1,  -1,  -1,  -1,  -1,  696, 697, 698, 699, -1,  -1,  -1,  -1,  720, 721, 722, 723, -1,  -1,  -1,  -1,  724, -1,  725, -1,  -1,  -1,  -1,  -1,
292       };
293       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
294     } else if (rank == 24) {
295       const PetscScalar arrLocalExpectedHere[] = {
296         396, 397, 398, 399, 400, 401, 402, 403, 404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 432, 433, 434, 435, 436, 437, 438, 439, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459,
297         460, 461, 462, 463, 448, 449, -1,  -1,  450, 451, -1,  -1,  452, 453, -1,  -1,  454, 455, -1,  -1,  464, 465, -1,  -1,  466, 467, -1,  -1,  594, 595, 596, 597, 598, 599, 600, 601, 602, 603, 604, 605, 606, 607, 608, 609,
298         618, 619, 620, 621, 622, 623, 624, 625, 648, 649, 650, 651, 652, 653, 654, 655, 656, 657, 658, 659, 660, 661, 662, 663, 684, 685, 686, 687, 688, 689, 690, 691, 664, 665, -1,  -1,  666, 667, -1,  -1,  668, 669, -1,  -1,
299         670, 671, -1,  -1,  692, 693, -1,  -1,  694, 695, -1,  -1,  610, 611, 612, 613, -1,  -1,  -1,  -1,  614, 615, 616, 617, -1,  -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,  672, 673, 674, 675, -1,  -1,  -1,  -1,
300         676, 677, 678, 679, -1,  -1,  -1,  -1,  696, 697, 698, 699, -1,  -1,  -1,  -1,  680, 681, -1,  -1,  -1,  -1,  -1,  -1,  682, 683, -1,  -1,  -1,  -1,  -1,  -1,  700, 701, -1,  -1,  -1,  -1,  -1,  -1,
301       };
302       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
303     } else if (rank == 25) {
304       const PetscScalar arrLocalExpectedHere[] = {
305         404, 405, 406, 407, 408, 409, 410, 411, 412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 440, 441, 442, 443, 444, 445, 446, 447, 456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471,
306         472, 473, 474, 475, 452, 453, -1,  -1,  454, 455, -1,  -1,  464, 465, -1,  -1,  466, 467, -1,  -1,  480, 481, -1,  -1,  482, 483, -1,  -1,  602, 603, 604, 605, 606, 607, 608, 609, 618, 619, 620, 621, 622, 623, 624, 625,
307         630, 631, 632, 633, 634, 635, 636, 637, 656, 657, 658, 659, 660, 661, 662, 663, 684, 685, 686, 687, 688, 689, 690, 691, 702, 703, 704, 705, 706, 707, 708, 709, 668, 669, -1,  -1,  670, 671, -1,  -1,  692, 693, -1,  -1,
308         694, 695, -1,  -1,  714, 715, -1,  -1,  716, 717, -1,  -1,  614, 615, 616, 617, -1,  -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,  642, 643, 644, 645, -1,  -1,  -1,  -1,  676, 677, 678, 679, -1,  -1,  -1,  -1,
309         696, 697, 698, 699, -1,  -1,  -1,  -1,  720, 721, 722, 723, -1,  -1,  -1,  -1,  682, 683, -1,  -1,  -1,  -1,  -1,  -1,  700, 701, -1,  -1,  -1,  -1,  -1,  -1,  726, 727, -1,  -1,  -1,  -1,  -1,  -1,
310       };
311       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
312     } else if (rank == 26) {
313       const PetscScalar arrLocalExpectedHere[] = {
314         412, 413, 414, 415, 416, 417, 418, 419, 420, 421, 422, 423, 424, 425, 426, 427, 428, -1,  429, -1,  430, -1,  431, -1,  456, 457, 458, 459, 460, 461, 462, 463, 468, 469, 470, 471, 472, 473, 474, 475, 476, -1,  477, -1,
315         478, -1,  479, -1,  464, 465, -1,  -1,  466, 467, -1,  -1,  480, 481, -1,  -1,  482, 483, -1,  -1,  484, -1,  -1,  -1,  485, -1,  -1,  -1,  618, 619, 620, 621, 622, 623, 624, 625, 630, 631, 632, 633, 634, 635, 636, 637,
316         638, -1,  639, -1,  640, -1,  641, -1,  684, 685, 686, 687, 688, 689, 690, 691, 702, 703, 704, 705, 706, 707, 708, 709, 710, -1,  711, -1,  712, -1,  713, -1,  692, 693, -1,  -1,  694, 695, -1,  -1,  714, 715, -1,  -1,
317         716, 717, -1,  -1,  718, -1,  -1,  -1,  719, -1,  -1,  -1,  626, 627, 628, 629, -1,  -1,  -1,  -1,  642, 643, 644, 645, -1,  -1,  -1,  -1,  646, -1,  647, -1,  -1,  -1,  -1,  -1,  696, 697, 698, 699, -1,  -1,  -1,  -1,
318         720, 721, 722, 723, -1,  -1,  -1,  -1,  724, -1,  725, -1,  -1,  -1,  -1,  -1,  700, 701, -1,  -1,  -1,  -1,  -1,  -1,  726, 727, -1,  -1,  -1,  -1,  -1,  -1,  728, -1,  -1,  -1,  -1,  -1,  -1,  -1,
319       };
320       for (i = 0; i < entriesGhost; ++i) arrLocalExpected[i] = arrLocalExpectedHere[i];
321     }
322 
323     PetscCall(VecGetArrayRead(vecLocal, &arrLocal));
324     for (i = 0, nerr = 0; i < entriesGhost; ++i) {
325       if (arrLocal[i] != arrLocalExpected[i]) {
326         ++nerr;
327         if (nerr <= maxErrPerRank) {
328           PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Entry %" PetscInt_FMT " has value %g instead of the expected %g\n", rank, i, (double)PetscRealPart(arrLocal[i]), (double)PetscRealPart(arrLocalExpected[i])));
329           if (nerr == maxErrPerRank + 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Skipping additional errors on this rank\n", rank));
330         }
331       }
332     }
333     if (nerr > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] %" PetscInt_FMT " incorrect values on this rank\n", rank, nerr));
334     PetscCall(VecRestoreArrayRead(vecLocal, &arrLocal));
335     PetscCall(PetscFree(arrLocalExpected));
336   }
337 
338   PetscCall(VecDestroy(&vecLocal));
339   PetscCall(VecDestroy(&vecGlobal));
340   PetscFunctionReturn(PETSC_SUCCESS);
341 }
342 
343 /*TEST
344 
345    test:
346       nsize: 27
347       output_file: output/empty.out
348 
349 TEST*/
350