Actual source code: ex56.c
1: static char help[] = "3D, tri-linear quadrilateral (Q1), displacement finite element formulation\n\
2: of linear elasticity. E=1.0, nu=0.25.\n\
3: Unit square domain with Dirichelet boundary condition on the y=0 side only.\n\
4: Load of 1.0 in x + 2y direction on all nodes (not a true uniform load).\n\
5: -ne <size> : number of (square) quadrilateral elements in each dimension\n\
6: -alpha <v> : scaling of material coefficient in embedded circle\n\n";
8: #include <petscksp.h>
10: static PetscBool log_stages = PETSC_TRUE;
12: static PetscErrorCode MaybeLogStagePush(PetscLogStage stage)
13: {
14: return log_stages ? PetscLogStagePush(stage) : PETSC_SUCCESS;
15: }
17: static PetscErrorCode MaybeLogStagePop()
18: {
19: return log_stages ? PetscLogStagePop() : PETSC_SUCCESS;
20: }
22: PetscErrorCode elem_3d_elast_v_25(PetscScalar *);
24: int main(int argc, char **args)
25: {
26: Mat Amat;
27: PetscInt m, nn, M, Istart, Iend, i, j, k, ii, jj, kk, ic, ne = 4, id;
28: PetscReal x, y, z, h, *coords, soft_alpha = 1.e-3;
29: PetscBool two_solves = PETSC_FALSE, test_nonzero_cols = PETSC_FALSE, use_nearnullspace = PETSC_FALSE, test_late_bs = PETSC_FALSE;
30: Vec xx, bb;
31: KSP ksp;
32: MPI_Comm comm;
33: PetscMPIInt npe, mype;
34: PetscScalar DD[24][24], DD2[24][24];
35: PetscLogStage stage[6];
36: PetscScalar DD1[24][24];
38: PetscFunctionBeginUser;
39: PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
40: comm = PETSC_COMM_WORLD;
41: PetscCallMPI(MPI_Comm_rank(comm, &mype));
42: PetscCallMPI(MPI_Comm_size(comm, &npe));
44: PetscOptionsBegin(comm, NULL, "3D bilinear Q1 elasticity options", "");
45: {
46: char nestring[256];
47: PetscCall(PetscSNPrintf(nestring, sizeof nestring, "number of elements in each direction, ne+1 must be a multiple of %" PetscInt_FMT " (sizes^{1/3})", (PetscInt)(PetscPowReal((PetscReal)npe, 1. / 3.) + .5)));
48: PetscCall(PetscOptionsInt("-ne", nestring, "", ne, &ne, NULL));
49: PetscCall(PetscOptionsBool("-log_stages", "Log stages of solve separately", "", log_stages, &log_stages, NULL));
50: PetscCall(PetscOptionsReal("-alpha", "material coefficient inside circle", "", soft_alpha, &soft_alpha, NULL));
51: PetscCall(PetscOptionsBool("-two_solves", "solve additional variant of the problem", "", two_solves, &two_solves, NULL));
52: PetscCall(PetscOptionsBool("-test_nonzero_cols", "nonzero test", "", test_nonzero_cols, &test_nonzero_cols, NULL));
53: PetscCall(PetscOptionsBool("-use_mat_nearnullspace", "MatNearNullSpace API test", "", use_nearnullspace, &use_nearnullspace, NULL));
54: PetscCall(PetscOptionsBool("-test_late_bs", "", "", test_late_bs, &test_late_bs, NULL));
55: }
56: PetscOptionsEnd();
58: if (log_stages) {
59: PetscCall(PetscLogStageRegister("Setup", &stage[0]));
60: PetscCall(PetscLogStageRegister("Solve", &stage[1]));
61: PetscCall(PetscLogStageRegister("2nd Setup", &stage[2]));
62: PetscCall(PetscLogStageRegister("2nd Solve", &stage[3]));
63: PetscCall(PetscLogStageRegister("3rd Setup", &stage[4]));
64: PetscCall(PetscLogStageRegister("3rd Solve", &stage[5]));
65: } else {
66: for (i = 0; i < (PetscInt)PETSC_STATIC_ARRAY_LENGTH(stage); i++) stage[i] = -1;
67: }
69: h = 1. / ne;
70: nn = ne + 1;
71: /* ne*ne; number of global elements */
72: M = 3 * nn * nn * nn; /* global number of equations */
73: if (npe == 2) {
74: if (mype == 1) m = 0;
75: else m = nn * nn * nn;
76: npe = 1;
77: } else {
78: m = nn * nn * nn / npe;
79: if (mype == npe - 1) m = nn * nn * nn - (npe - 1) * m;
80: }
81: m *= 3; /* number of equations local*/
82: /* Setup solver */
83: PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
84: PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
85: PetscCall(KSPSetFromOptions(ksp));
86: {
87: /* configuration */
88: const PetscInt NP = (PetscInt)(PetscPowReal((PetscReal)npe, 1. / 3.) + .5);
89: const PetscInt ipx = mype % NP, ipy = (mype % (NP * NP)) / NP, ipz = mype / (NP * NP);
90: const PetscInt Ni0 = ipx * (nn / NP), Nj0 = ipy * (nn / NP), Nk0 = ipz * (nn / NP);
91: const PetscInt Ni1 = Ni0 + (m > 0 ? (nn / NP) : 0), Nj1 = Nj0 + (nn / NP), Nk1 = Nk0 + (nn / NP);
92: const PetscInt NN = nn / NP, id0 = ipz * nn * nn * NN + ipy * nn * NN * NN + ipx * NN * NN * NN;
93: PetscInt *d_nnz, *o_nnz, osz[4] = {0, 9, 15, 19}, nbc;
94: PetscScalar vv[24], v2[24];
95: PetscCheck(npe == NP * NP * NP, comm, PETSC_ERR_ARG_WRONG, "npe=%d: npe^{1/3} must be integer", npe);
96: PetscCheck(nn == NP * (nn / NP), comm, PETSC_ERR_ARG_WRONG, "-ne %" PetscInt_FMT ": (ne+1)%%(npe^{1/3}) must equal zero", ne);
98: /* count nnz */
99: PetscCall(PetscMalloc1(m + 1, &d_nnz));
100: PetscCall(PetscMalloc1(m + 1, &o_nnz));
101: for (i = Ni0, ic = 0; i < Ni1; i++) {
102: for (j = Nj0; j < Nj1; j++) {
103: for (k = Nk0; k < Nk1; k++) {
104: nbc = 0;
105: if (i == Ni0 || i == Ni1 - 1) nbc++;
106: if (j == Nj0 || j == Nj1 - 1) nbc++;
107: if (k == Nk0 || k == Nk1 - 1) nbc++;
108: for (jj = 0; jj < 3; jj++, ic++) {
109: d_nnz[ic] = 3 * (27 - osz[nbc]);
110: o_nnz[ic] = 3 * osz[nbc];
111: }
112: }
113: }
114: }
115: PetscCheck(ic == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ic %" PetscInt_FMT " does not equal m %" PetscInt_FMT, ic, m);
117: /* create stiffness matrix */
118: PetscCall(MatCreate(comm, &Amat));
119: PetscCall(MatSetSizes(Amat, m, m, M, M));
120: if (!test_late_bs) PetscCall(MatSetBlockSize(Amat, 3));
121: PetscCall(MatSetType(Amat, MATAIJ));
122: PetscCall(MatSetOption(Amat, MAT_SPD, PETSC_TRUE));
123: PetscCall(MatSetOption(Amat, MAT_SPD_ETERNAL, PETSC_TRUE));
124: PetscCall(MatSetFromOptions(Amat));
125: PetscCall(MatSeqAIJSetPreallocation(Amat, 0, d_nnz));
126: PetscCall(MatMPIAIJSetPreallocation(Amat, 0, d_nnz, 0, o_nnz));
128: PetscCall(PetscFree(d_nnz));
129: PetscCall(PetscFree(o_nnz));
130: PetscCall(MatCreateVecs(Amat, &bb, &xx));
132: PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
134: PetscCheck(m == Iend - Istart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "m %" PetscInt_FMT " does not equal Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT, m, Iend, Istart);
135: /* generate element matrices */
136: {
137: PetscBool hasData = PETSC_TRUE;
138: if (!hasData) {
139: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t No data is provided\n"));
140: for (i = 0; i < 24; i++) {
141: for (j = 0; j < 24; j++) {
142: if (i == j) DD1[i][j] = 1.0;
143: else DD1[i][j] = -.25;
144: }
145: }
146: } else {
147: /* Get array data */
148: PetscCall(elem_3d_elast_v_25((PetscScalar *)DD1));
149: }
151: /* BC version of element */
152: for (i = 0; i < 24; i++) {
153: for (j = 0; j < 24; j++) {
154: if (i < 12 || (j < 12 && !test_nonzero_cols)) {
155: if (i == j) DD2[i][j] = 0.1 * DD1[i][j];
156: else DD2[i][j] = 0.0;
157: } else DD2[i][j] = DD1[i][j];
158: }
159: }
160: /* element residual/load vector */
161: for (i = 0; i < 24; i++) {
162: if (i % 3 == 0) vv[i] = h * h;
163: else if (i % 3 == 1) vv[i] = 2.0 * h * h;
164: else vv[i] = .0;
165: }
166: for (i = 0; i < 24; i++) {
167: if (i % 3 == 0 && i >= 12) v2[i] = h * h;
168: else if (i % 3 == 1 && i >= 12) v2[i] = 2.0 * h * h;
169: else v2[i] = .0;
170: }
171: }
173: PetscCall(PetscMalloc1(m + 1, &coords));
174: coords[m] = -99.0;
176: /* forms the element stiffness and coordinates */
177: for (i = Ni0, ic = 0, ii = 0; i < Ni1; i++, ii++) {
178: for (j = Nj0, jj = 0; j < Nj1; j++, jj++) {
179: for (k = Nk0, kk = 0; k < Nk1; k++, kk++, ic++) {
180: /* coords */
181: x = coords[3 * ic] = h * (PetscReal)i;
182: y = coords[3 * ic + 1] = h * (PetscReal)j;
183: z = coords[3 * ic + 2] = h * (PetscReal)k;
184: /* matrix */
185: id = id0 + ii + NN * jj + NN * NN * kk;
186: if (i < ne && j < ne && k < ne) {
187: /* radius */
188: PetscReal radius = PetscSqrtReal((x - .5 + h / 2) * (x - .5 + h / 2) + (y - .5 + h / 2) * (y - .5 + h / 2) + (z - .5 + h / 2) * (z - .5 + h / 2));
189: PetscReal alpha = 1.0;
190: PetscInt jx, ix, idx[8], idx3[24];
191: idx[0] = id;
192: idx[1] = id + 1;
193: idx[2] = id + NN + 1;
194: idx[3] = id + NN;
195: idx[4] = id + NN * NN;
196: idx[5] = id + 1 + NN * NN;
197: idx[6] = id + NN + 1 + NN * NN;
198: idx[7] = id + NN + NN * NN;
200: /* correct indices */
201: if (i == Ni1 - 1 && Ni1 != nn) {
202: idx[1] += NN * (NN * NN - 1);
203: idx[2] += NN * (NN * NN - 1);
204: idx[5] += NN * (NN * NN - 1);
205: idx[6] += NN * (NN * NN - 1);
206: }
207: if (j == Nj1 - 1 && Nj1 != nn) {
208: idx[2] += NN * NN * (nn - 1);
209: idx[3] += NN * NN * (nn - 1);
210: idx[6] += NN * NN * (nn - 1);
211: idx[7] += NN * NN * (nn - 1);
212: }
213: if (k == Nk1 - 1 && Nk1 != nn) {
214: idx[4] += NN * (nn * nn - NN * NN);
215: idx[5] += NN * (nn * nn - NN * NN);
216: idx[6] += NN * (nn * nn - NN * NN);
217: idx[7] += NN * (nn * nn - NN * NN);
218: }
220: if (radius < 0.25) alpha = soft_alpha;
222: for (ix = 0; ix < 24; ix++) {
223: for (jx = 0; jx < 24; jx++) DD[ix][jx] = alpha * DD1[ix][jx];
224: }
225: if (k > 0) {
226: if (!test_late_bs) {
227: PetscCall(MatSetValuesBlocked(Amat, 8, idx, 8, idx, (const PetscScalar *)DD, ADD_VALUES));
228: PetscCall(VecSetValuesBlocked(bb, 8, idx, (const PetscScalar *)vv, ADD_VALUES));
229: } else {
230: for (ix = 0; ix < 8; ix++) {
231: idx3[3 * ix] = 3 * idx[ix];
232: idx3[3 * ix + 1] = 3 * idx[ix] + 1;
233: idx3[3 * ix + 2] = 3 * idx[ix] + 2;
234: }
235: PetscCall(MatSetValues(Amat, 24, idx3, 24, idx3, (const PetscScalar *)DD, ADD_VALUES));
236: PetscCall(VecSetValues(bb, 24, idx3, (const PetscScalar *)vv, ADD_VALUES));
237: }
238: } else {
239: /* a BC */
240: for (ix = 0; ix < 24; ix++) {
241: for (jx = 0; jx < 24; jx++) DD[ix][jx] = alpha * DD2[ix][jx];
242: }
243: if (!test_late_bs) {
244: PetscCall(MatSetValuesBlocked(Amat, 8, idx, 8, idx, (const PetscScalar *)DD, ADD_VALUES));
245: PetscCall(VecSetValuesBlocked(bb, 8, idx, (const PetscScalar *)v2, ADD_VALUES));
246: } else {
247: for (ix = 0; ix < 8; ix++) {
248: idx3[3 * ix] = 3 * idx[ix];
249: idx3[3 * ix + 1] = 3 * idx[ix] + 1;
250: idx3[3 * ix + 2] = 3 * idx[ix] + 2;
251: }
252: PetscCall(MatSetValues(Amat, 24, idx3, 24, idx3, (const PetscScalar *)DD, ADD_VALUES));
253: PetscCall(VecSetValues(bb, 24, idx3, (const PetscScalar *)v2, ADD_VALUES));
254: }
255: }
256: }
257: }
258: }
259: }
260: PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY));
261: PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY));
262: PetscCall(VecAssemblyBegin(bb));
263: PetscCall(VecAssemblyEnd(bb));
264: }
265: PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY));
266: PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY));
267: PetscCall(VecAssemblyBegin(bb));
268: PetscCall(VecAssemblyEnd(bb));
269: if (test_late_bs) {
270: PetscCall(VecSetBlockSize(xx, 3));
271: PetscCall(VecSetBlockSize(bb, 3));
272: PetscCall(MatSetBlockSize(Amat, 3));
273: }
275: if (!PETSC_TRUE) {
276: PetscViewer viewer;
277: PetscCall(PetscViewerASCIIOpen(comm, "Amat.m", &viewer));
278: PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
279: PetscCall(MatView(Amat, viewer));
280: PetscCall(PetscViewerPopFormat(viewer));
281: PetscCall(PetscViewerDestroy(&viewer));
282: }
284: /* finish KSP/PC setup */
285: PetscCall(KSPSetOperators(ksp, Amat, Amat));
286: if (use_nearnullspace) {
287: MatNullSpace matnull;
288: Vec vec_coords;
289: PetscScalar *c;
291: PetscCall(VecCreate(MPI_COMM_WORLD, &vec_coords));
292: PetscCall(VecSetBlockSize(vec_coords, 3));
293: PetscCall(VecSetSizes(vec_coords, m, PETSC_DECIDE));
294: PetscCall(VecSetUp(vec_coords));
295: PetscCall(VecGetArray(vec_coords, &c));
296: for (i = 0; i < m; i++) c[i] = coords[i]; /* Copy since Scalar type might be Complex */
297: PetscCall(VecRestoreArray(vec_coords, &c));
298: PetscCall(MatNullSpaceCreateRigidBody(vec_coords, &matnull));
299: PetscCall(MatSetNearNullSpace(Amat, matnull));
300: PetscCall(MatNullSpaceDestroy(&matnull));
301: PetscCall(VecDestroy(&vec_coords));
302: } else {
303: PC pc;
304: PetscCall(KSPGetPC(ksp, &pc));
305: PetscCall(PCSetCoordinates(pc, 3, m / 3, coords));
306: PetscCall(PCGAMGSetUseSAEstEig(pc, PETSC_FALSE));
307: }
309: PetscCall(MaybeLogStagePush(stage[0]));
311: /* PC setup basically */
312: PetscCall(KSPSetUp(ksp));
314: PetscCall(MaybeLogStagePop());
315: PetscCall(MaybeLogStagePush(stage[1]));
317: /* test BCs */
318: if (test_nonzero_cols) {
319: PetscCall(VecZeroEntries(xx));
320: if (mype == 0) PetscCall(VecSetValue(xx, 0, 1.0, INSERT_VALUES));
321: PetscCall(VecAssemblyBegin(xx));
322: PetscCall(VecAssemblyEnd(xx));
323: PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
324: }
326: /* 1st solve */
327: PetscCall(KSPSolve(ksp, bb, xx));
329: PetscCall(MaybeLogStagePop());
331: /* 2nd solve */
332: if (two_solves) {
333: PetscReal emax, emin;
334: PetscReal norm, norm2;
335: Vec res;
337: PetscCall(MaybeLogStagePush(stage[2]));
338: /* PC setup basically */
339: PetscCall(MatScale(Amat, 100000.0));
340: PetscCall(KSPSetOperators(ksp, Amat, Amat));
341: PetscCall(KSPSetUp(ksp));
343: PetscCall(MaybeLogStagePop());
344: PetscCall(MaybeLogStagePush(stage[3]));
345: PetscCall(KSPSolve(ksp, bb, xx));
346: PetscCall(KSPComputeExtremeSingularValues(ksp, &emax, &emin));
348: PetscCall(MaybeLogStagePop());
349: PetscCall(MaybeLogStagePush(stage[4]));
351: PetscCall(MaybeLogStagePop());
352: PetscCall(MaybeLogStagePush(stage[5]));
354: /* 3rd solve */
355: PetscCall(KSPSolve(ksp, bb, xx));
357: PetscCall(MaybeLogStagePop());
359: PetscCall(VecNorm(bb, NORM_2, &norm2));
361: PetscCall(VecDuplicate(xx, &res));
362: PetscCall(MatMult(Amat, xx, res));
363: PetscCall(VecAXPY(bb, -1.0, res));
364: PetscCall(VecDestroy(&res));
365: PetscCall(VecNorm(bb, NORM_2, &norm));
366: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d]%s |b-Ax|/|b|=%e, |b|=%e, emax=%e\n", 0, PETSC_FUNCTION_NAME, (double)(norm / norm2), (double)norm2, (double)emax));
367: }
369: /* Free work space */
370: PetscCall(KSPDestroy(&ksp));
371: PetscCall(VecDestroy(&xx));
372: PetscCall(VecDestroy(&bb));
373: PetscCall(MatDestroy(&Amat));
374: PetscCall(PetscFree(coords));
376: PetscCall(PetscFinalize());
377: return 0;
378: }
380: /* Data was previously provided in the file data/elem_3d_elast_v_25.tx */
381: PetscErrorCode elem_3d_elast_v_25(PetscScalar *dd)
382: {
383: PetscScalar DD[] = {
384: 0.18981481481481474, 5.27777777777777568E-002, 5.27777777777777568E-002, -5.64814814814814659E-002, -1.38888888888889072E-002, -1.38888888888889089E-002, -8.24074074074073876E-002, -5.27777777777777429E-002, 1.38888888888888725E-002,
385: 4.90740740740740339E-002, 1.38888888888889124E-002, 4.72222222222222071E-002, 4.90740740740740339E-002, 4.72222222222221932E-002, 1.38888888888888968E-002, -8.24074074074073876E-002, 1.38888888888888673E-002, -5.27777777777777429E-002,
386: -7.87037037037036785E-002, -4.72222222222221932E-002, -4.72222222222222071E-002, 1.20370370370370180E-002, -1.38888888888888742E-002, -1.38888888888888829E-002, 5.27777777777777568E-002, 0.18981481481481474, 5.27777777777777568E-002,
387: 1.38888888888889124E-002, 4.90740740740740269E-002, 4.72222222222221932E-002, -5.27777777777777637E-002, -8.24074074074073876E-002, 1.38888888888888725E-002, -1.38888888888889037E-002, -5.64814814814814728E-002, -1.38888888888888985E-002,
388: 4.72222222222221932E-002, 4.90740740740740478E-002, 1.38888888888888968E-002, -1.38888888888888673E-002, 1.20370370370370058E-002, -1.38888888888888742E-002, -4.72222222222221932E-002, -7.87037037037036785E-002, -4.72222222222222002E-002,
389: 1.38888888888888742E-002, -8.24074074074073598E-002, -5.27777777777777568E-002, 5.27777777777777568E-002, 5.27777777777777568E-002, 0.18981481481481474, 1.38888888888889055E-002, 4.72222222222222002E-002, 4.90740740740740269E-002,
390: -1.38888888888888829E-002, -1.38888888888888829E-002, 1.20370370370370180E-002, 4.72222222222222002E-002, 1.38888888888888985E-002, 4.90740740740740339E-002, -1.38888888888888985E-002, -1.38888888888888968E-002, -5.64814814814814520E-002,
391: -5.27777777777777568E-002, 1.38888888888888777E-002, -8.24074074074073876E-002, -4.72222222222222002E-002, -4.72222222222221932E-002, -7.87037037037036646E-002, 1.38888888888888794E-002, -5.27777777777777568E-002, -8.24074074074073598E-002,
392: -5.64814814814814659E-002, 1.38888888888889124E-002, 1.38888888888889055E-002, 0.18981481481481474, -5.27777777777777568E-002, -5.27777777777777499E-002, 4.90740740740740269E-002, -1.38888888888889072E-002, -4.72222222222221932E-002,
393: -8.24074074074073876E-002, 5.27777777777777568E-002, -1.38888888888888812E-002, -8.24074074074073876E-002, -1.38888888888888742E-002, 5.27777777777777499E-002, 4.90740740740740269E-002, -4.72222222222221863E-002, -1.38888888888889089E-002,
394: 1.20370370370370162E-002, 1.38888888888888673E-002, 1.38888888888888742E-002, -7.87037037037036785E-002, 4.72222222222222002E-002, 4.72222222222222071E-002, -1.38888888888889072E-002, 4.90740740740740269E-002, 4.72222222222222002E-002,
395: -5.27777777777777568E-002, 0.18981481481481480, 5.27777777777777568E-002, 1.38888888888889020E-002, -5.64814814814814728E-002, -1.38888888888888951E-002, 5.27777777777777637E-002, -8.24074074074073876E-002, 1.38888888888888881E-002,
396: 1.38888888888888742E-002, 1.20370370370370232E-002, -1.38888888888888812E-002, -4.72222222222221863E-002, 4.90740740740740339E-002, 1.38888888888888933E-002, -1.38888888888888812E-002, -8.24074074074073876E-002, -5.27777777777777568E-002,
397: 4.72222222222222071E-002, -7.87037037037036924E-002, -4.72222222222222140E-002, -1.38888888888889089E-002, 4.72222222222221932E-002, 4.90740740740740269E-002, -5.27777777777777499E-002, 5.27777777777777568E-002, 0.18981481481481477,
398: -4.72222222222222071E-002, 1.38888888888888968E-002, 4.90740740740740131E-002, 1.38888888888888812E-002, -1.38888888888888708E-002, 1.20370370370370267E-002, 5.27777777777777568E-002, 1.38888888888888812E-002, -8.24074074074073876E-002,
399: 1.38888888888889124E-002, -1.38888888888889055E-002, -5.64814814814814589E-002, -1.38888888888888812E-002, -5.27777777777777568E-002, -8.24074074074073737E-002, 4.72222222222222002E-002, -4.72222222222222002E-002, -7.87037037037036924E-002,
400: -8.24074074074073876E-002, -5.27777777777777637E-002, -1.38888888888888829E-002, 4.90740740740740269E-002, 1.38888888888889020E-002, -4.72222222222222071E-002, 0.18981481481481480, 5.27777777777777637E-002, -5.27777777777777637E-002,
401: -5.64814814814814728E-002, -1.38888888888889037E-002, 1.38888888888888951E-002, -7.87037037037036785E-002, -4.72222222222222002E-002, 4.72222222222221932E-002, 1.20370370370370128E-002, -1.38888888888888725E-002, 1.38888888888888812E-002,
402: 4.90740740740740408E-002, 4.72222222222222002E-002, -1.38888888888888951E-002, -8.24074074074073876E-002, 1.38888888888888812E-002, 5.27777777777777637E-002, -5.27777777777777429E-002, -8.24074074074073876E-002, -1.38888888888888829E-002,
403: -1.38888888888889072E-002, -5.64814814814814728E-002, 1.38888888888888968E-002, 5.27777777777777637E-002, 0.18981481481481480, -5.27777777777777568E-002, 1.38888888888888916E-002, 4.90740740740740339E-002, -4.72222222222222210E-002,
404: -4.72222222222221932E-002, -7.87037037037036924E-002, 4.72222222222222002E-002, 1.38888888888888742E-002, -8.24074074074073876E-002, 5.27777777777777429E-002, 4.72222222222222002E-002, 4.90740740740740269E-002, -1.38888888888888951E-002,
405: -1.38888888888888846E-002, 1.20370370370370267E-002, 1.38888888888888916E-002, 1.38888888888888725E-002, 1.38888888888888725E-002, 1.20370370370370180E-002, -4.72222222222221932E-002, -1.38888888888888951E-002, 4.90740740740740131E-002,
406: -5.27777777777777637E-002, -5.27777777777777568E-002, 0.18981481481481480, -1.38888888888888968E-002, -4.72222222222221932E-002, 4.90740740740740339E-002, 4.72222222222221932E-002, 4.72222222222222071E-002, -7.87037037037036646E-002,
407: -1.38888888888888742E-002, 5.27777777777777499E-002, -8.24074074074073737E-002, 1.38888888888888933E-002, 1.38888888888889020E-002, -5.64814814814814589E-002, 5.27777777777777568E-002, -1.38888888888888794E-002, -8.24074074074073876E-002,
408: 4.90740740740740339E-002, -1.38888888888889037E-002, 4.72222222222222002E-002, -8.24074074074073876E-002, 5.27777777777777637E-002, 1.38888888888888812E-002, -5.64814814814814728E-002, 1.38888888888888916E-002, -1.38888888888888968E-002,
409: 0.18981481481481480, -5.27777777777777499E-002, 5.27777777777777707E-002, 1.20370370370370180E-002, 1.38888888888888812E-002, -1.38888888888888812E-002, -7.87037037037036785E-002, 4.72222222222222002E-002, -4.72222222222222071E-002,
410: -8.24074074074073876E-002, -1.38888888888888742E-002, -5.27777777777777568E-002, 4.90740740740740616E-002, -4.72222222222222002E-002, 1.38888888888888846E-002, 1.38888888888889124E-002, -5.64814814814814728E-002, 1.38888888888888985E-002,
411: 5.27777777777777568E-002, -8.24074074074073876E-002, -1.38888888888888708E-002, -1.38888888888889037E-002, 4.90740740740740339E-002, -4.72222222222221932E-002, -5.27777777777777499E-002, 0.18981481481481480, -5.27777777777777568E-002,
412: -1.38888888888888673E-002, -8.24074074074073598E-002, 5.27777777777777429E-002, 4.72222222222222002E-002, -7.87037037037036785E-002, 4.72222222222222002E-002, 1.38888888888888708E-002, 1.20370370370370128E-002, 1.38888888888888760E-002,
413: -4.72222222222222002E-002, 4.90740740740740478E-002, -1.38888888888888951E-002, 4.72222222222222071E-002, -1.38888888888888985E-002, 4.90740740740740339E-002, -1.38888888888888812E-002, 1.38888888888888881E-002, 1.20370370370370267E-002,
414: 1.38888888888888951E-002, -4.72222222222222210E-002, 4.90740740740740339E-002, 5.27777777777777707E-002, -5.27777777777777568E-002, 0.18981481481481477, 1.38888888888888829E-002, 5.27777777777777707E-002, -8.24074074074073598E-002,
415: -4.72222222222222140E-002, 4.72222222222222140E-002, -7.87037037037036646E-002, -5.27777777777777707E-002, -1.38888888888888829E-002, -8.24074074074073876E-002, -1.38888888888888881E-002, 1.38888888888888881E-002, -5.64814814814814589E-002,
416: 4.90740740740740339E-002, 4.72222222222221932E-002, -1.38888888888888985E-002, -8.24074074074073876E-002, 1.38888888888888742E-002, 5.27777777777777568E-002, -7.87037037037036785E-002, -4.72222222222221932E-002, 4.72222222222221932E-002,
417: 1.20370370370370180E-002, -1.38888888888888673E-002, 1.38888888888888829E-002, 0.18981481481481469, 5.27777777777777429E-002, -5.27777777777777429E-002, -5.64814814814814659E-002, -1.38888888888889055E-002, 1.38888888888889055E-002,
418: -8.24074074074074153E-002, -5.27777777777777429E-002, -1.38888888888888760E-002, 4.90740740740740408E-002, 1.38888888888888968E-002, -4.72222222222222071E-002, 4.72222222222221932E-002, 4.90740740740740478E-002, -1.38888888888888968E-002,
419: -1.38888888888888742E-002, 1.20370370370370232E-002, 1.38888888888888812E-002, -4.72222222222222002E-002, -7.87037037037036924E-002, 4.72222222222222071E-002, 1.38888888888888812E-002, -8.24074074074073598E-002, 5.27777777777777707E-002,
420: 5.27777777777777429E-002, 0.18981481481481477, -5.27777777777777499E-002, 1.38888888888889107E-002, 4.90740740740740478E-002, -4.72222222222221932E-002, -5.27777777777777568E-002, -8.24074074074074153E-002, -1.38888888888888812E-002,
421: -1.38888888888888846E-002, -5.64814814814814659E-002, 1.38888888888888812E-002, 1.38888888888888968E-002, 1.38888888888888968E-002, -5.64814814814814520E-002, 5.27777777777777499E-002, -1.38888888888888812E-002, -8.24074074074073876E-002,
422: 4.72222222222221932E-002, 4.72222222222222002E-002, -7.87037037037036646E-002, -1.38888888888888812E-002, 5.27777777777777429E-002, -8.24074074074073598E-002, -5.27777777777777429E-002, -5.27777777777777499E-002, 0.18981481481481474,
423: -1.38888888888888985E-002, -4.72222222222221863E-002, 4.90740740740740339E-002, 1.38888888888888829E-002, 1.38888888888888777E-002, 1.20370370370370249E-002, -4.72222222222222002E-002, -1.38888888888888933E-002, 4.90740740740740339E-002,
424: -8.24074074074073876E-002, -1.38888888888888673E-002, -5.27777777777777568E-002, 4.90740740740740269E-002, -4.72222222222221863E-002, 1.38888888888889124E-002, 1.20370370370370128E-002, 1.38888888888888742E-002, -1.38888888888888742E-002,
425: -7.87037037037036785E-002, 4.72222222222222002E-002, -4.72222222222222140E-002, -5.64814814814814659E-002, 1.38888888888889107E-002, -1.38888888888888985E-002, 0.18981481481481474, -5.27777777777777499E-002, 5.27777777777777499E-002,
426: 4.90740740740740339E-002, -1.38888888888889055E-002, 4.72222222222221932E-002, -8.24074074074074153E-002, 5.27777777777777499E-002, 1.38888888888888829E-002, 1.38888888888888673E-002, 1.20370370370370058E-002, 1.38888888888888777E-002,
427: -4.72222222222221863E-002, 4.90740740740740339E-002, -1.38888888888889055E-002, -1.38888888888888725E-002, -8.24074074074073876E-002, 5.27777777777777499E-002, 4.72222222222222002E-002, -7.87037037037036785E-002, 4.72222222222222140E-002,
428: -1.38888888888889055E-002, 4.90740740740740478E-002, -4.72222222222221863E-002, -5.27777777777777499E-002, 0.18981481481481469, -5.27777777777777499E-002, 1.38888888888889072E-002, -5.64814814814814659E-002, 1.38888888888889003E-002,
429: 5.27777777777777429E-002, -8.24074074074074153E-002, -1.38888888888888812E-002, -5.27777777777777429E-002, -1.38888888888888742E-002, -8.24074074074073876E-002, -1.38888888888889089E-002, 1.38888888888888933E-002, -5.64814814814814589E-002,
430: 1.38888888888888812E-002, 5.27777777777777429E-002, -8.24074074074073737E-002, -4.72222222222222071E-002, 4.72222222222222002E-002, -7.87037037037036646E-002, 1.38888888888889055E-002, -4.72222222222221932E-002, 4.90740740740740339E-002,
431: 5.27777777777777499E-002, -5.27777777777777499E-002, 0.18981481481481474, 4.72222222222222002E-002, -1.38888888888888985E-002, 4.90740740740740339E-002, -1.38888888888888846E-002, 1.38888888888888812E-002, 1.20370370370370284E-002,
432: -7.87037037037036785E-002, -4.72222222222221932E-002, -4.72222222222222002E-002, 1.20370370370370162E-002, -1.38888888888888812E-002, -1.38888888888888812E-002, 4.90740740740740408E-002, 4.72222222222222002E-002, 1.38888888888888933E-002,
433: -8.24074074074073876E-002, 1.38888888888888708E-002, -5.27777777777777707E-002, -8.24074074074074153E-002, -5.27777777777777568E-002, 1.38888888888888829E-002, 4.90740740740740339E-002, 1.38888888888889072E-002, 4.72222222222222002E-002,
434: 0.18981481481481477, 5.27777777777777429E-002, 5.27777777777777568E-002, -5.64814814814814659E-002, -1.38888888888888846E-002, -1.38888888888888881E-002, -4.72222222222221932E-002, -7.87037037037036785E-002, -4.72222222222221932E-002,
435: 1.38888888888888673E-002, -8.24074074074073876E-002, -5.27777777777777568E-002, 4.72222222222222002E-002, 4.90740740740740269E-002, 1.38888888888889020E-002, -1.38888888888888742E-002, 1.20370370370370128E-002, -1.38888888888888829E-002,
436: -5.27777777777777429E-002, -8.24074074074074153E-002, 1.38888888888888777E-002, -1.38888888888889055E-002, -5.64814814814814659E-002, -1.38888888888888985E-002, 5.27777777777777429E-002, 0.18981481481481469, 5.27777777777777429E-002,
437: 1.38888888888888933E-002, 4.90740740740740339E-002, 4.72222222222222071E-002, -4.72222222222222071E-002, -4.72222222222222002E-002, -7.87037037037036646E-002, 1.38888888888888742E-002, -5.27777777777777568E-002, -8.24074074074073737E-002,
438: -1.38888888888888951E-002, -1.38888888888888951E-002, -5.64814814814814589E-002, -5.27777777777777568E-002, 1.38888888888888760E-002, -8.24074074074073876E-002, -1.38888888888888760E-002, -1.38888888888888812E-002, 1.20370370370370249E-002,
439: 4.72222222222221932E-002, 1.38888888888889003E-002, 4.90740740740740339E-002, 5.27777777777777568E-002, 5.27777777777777429E-002, 0.18981481481481474, 1.38888888888888933E-002, 4.72222222222222071E-002, 4.90740740740740339E-002,
440: 1.20370370370370180E-002, 1.38888888888888742E-002, 1.38888888888888794E-002, -7.87037037037036785E-002, 4.72222222222222071E-002, 4.72222222222222002E-002, -8.24074074074073876E-002, -1.38888888888888846E-002, 5.27777777777777568E-002,
441: 4.90740740740740616E-002, -4.72222222222222002E-002, -1.38888888888888881E-002, 4.90740740740740408E-002, -1.38888888888888846E-002, -4.72222222222222002E-002, -8.24074074074074153E-002, 5.27777777777777429E-002, -1.38888888888888846E-002,
442: -5.64814814814814659E-002, 1.38888888888888933E-002, 1.38888888888888933E-002, 0.18981481481481477, -5.27777777777777568E-002, -5.27777777777777637E-002, -1.38888888888888742E-002, -8.24074074074073598E-002, -5.27777777777777568E-002,
443: 4.72222222222222002E-002, -7.87037037037036924E-002, -4.72222222222222002E-002, 1.38888888888888812E-002, 1.20370370370370267E-002, -1.38888888888888794E-002, -4.72222222222222002E-002, 4.90740740740740478E-002, 1.38888888888888881E-002,
444: 1.38888888888888968E-002, -5.64814814814814659E-002, -1.38888888888888933E-002, 5.27777777777777499E-002, -8.24074074074074153E-002, 1.38888888888888812E-002, -1.38888888888888846E-002, 4.90740740740740339E-002, 4.72222222222222071E-002,
445: -5.27777777777777568E-002, 0.18981481481481477, 5.27777777777777637E-002, -1.38888888888888829E-002, -5.27777777777777568E-002, -8.24074074074073598E-002, 4.72222222222222071E-002, -4.72222222222222140E-002, -7.87037037037036924E-002,
446: 5.27777777777777637E-002, 1.38888888888888916E-002, -8.24074074074073876E-002, 1.38888888888888846E-002, -1.38888888888888951E-002, -5.64814814814814589E-002, -4.72222222222222071E-002, 1.38888888888888812E-002, 4.90740740740740339E-002,
447: 1.38888888888888829E-002, -1.38888888888888812E-002, 1.20370370370370284E-002, -1.38888888888888881E-002, 4.72222222222222071E-002, 4.90740740740740339E-002, -5.27777777777777637E-002, 5.27777777777777637E-002, 0.18981481481481477,
448: };
450: PetscFunctionBeginUser;
451: PetscCall(PetscArraycpy(dd, DD, 576));
452: PetscFunctionReturn(PETSC_SUCCESS);
453: }
455: /*TEST
457: test:
458: nsize: 8
459: args: -ne 13 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -two_solves -ksp_converged_reason -use_mat_nearnullspace -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_threshold -0.01 -pc_gamg_coarse_eq_limit 200 -pc_gamg_process_eq_limit 30 -pc_gamg_repartition false -pc_mg_cycle_type v -pc_gamg_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short -pc_gamg_rank_reduction_factors 2,2,2
460: filter: grep -v variant
462: test:
463: suffix: 2
464: nsize: 1
465: args: -ne 11 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -two_solves -ksp_converged_reason -use_mat_nearnullspace -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_asm_use_agg true -mg_levels_sub_pc_type lu -mg_levels_pc_asm_overlap 0 -pc_gamg_threshold -0.01 -pc_gamg_coarse_eq_limit 200 -pc_gamg_process_eq_limit 30 -pc_gamg_repartition false -pc_mg_cycle_type v -pc_gamg_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg
466: filter: grep -v variant
468: test:
469: suffix: latebs
470: filter: grep -v variant
471: nsize: 8
472: args: -test_late_bs 0 -ne 9 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -two_solves -ksp_converged_reason -use_mat_nearnullspace false -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_threshold -0.01 -pc_gamg_coarse_eq_limit 200 -pc_gamg_process_eq_limit 30 -pc_gamg_repartition false -pc_mg_cycle_type v -pc_gamg_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short -ksp_view
474: test:
475: suffix: latebs-2
476: filter: grep -v variant
477: nsize: 8
478: args: -test_late_bs -ne 9 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -two_solves -ksp_converged_reason -use_mat_nearnullspace -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_threshold -0.01 -pc_gamg_coarse_eq_limit 200 -pc_gamg_process_eq_limit 30 -pc_gamg_repartition false -pc_mg_cycle_type v -pc_gamg_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short -ksp_view
480: test:
481: suffix: ml
482: nsize: 8
483: args: -ne 9 -alpha 1.e-3 -ksp_type cg -pc_type ml -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type sor -ksp_monitor_short
484: requires: ml
486: test:
487: suffix: nns
488: args: -ne 9 -alpha 1.e-3 -ksp_converged_reason -ksp_type cg -ksp_max_it 50 -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -mg_levels_ksp_type chebyshev -mg_levels_pc_type sor -pc_gamg_reuse_interpolation true -two_solves -use_mat_nearnullspace -pc_gamg_use_sa_esteig 0 -mg_levels_esteig_ksp_max_it 10
490: test:
491: suffix: nns_telescope
492: nsize: 2
493: args: -use_mat_nearnullspace -ksp_monitor_short -pc_type telescope -pc_telescope_reduction_factor 2 -telescope_pc_type gamg -telescope_pc_gamg_esteig_ksp_type cg -telescope_pc_gamg_esteig_ksp_max_it 10
495: test:
496: suffix: nns_gdsw
497: filter: grep -v "variant HERMITIAN"
498: nsize: 8
499: args: -use_mat_nearnullspace -ksp_monitor_short -pc_type mg -pc_mg_levels 2 -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_galerkin -mg_levels_pc_type bjacobi -ne 3 -ksp_view
501: test:
502: suffix: seqaijmkl
503: nsize: 8
504: requires: mkl_sparse
505: args: -ne 9 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -two_solves -ksp_converged_reason -use_mat_nearnullspace -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_pc_type jacobi -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_threshold 0.01 -pc_gamg_coarse_eq_limit 2000 -pc_gamg_process_eq_limit 200 -pc_gamg_repartition false -pc_mg_cycle_type v -ksp_monitor_short -mat_seqaij_type seqaijmkl
507: TEST*/