Actual source code: ex56.c
petsc-3.10.5 2019-03-28
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 coeficient in embedded circle\n\n";
8: #include <petscksp.h>
10: static PetscBool log_stages = PETSC_TRUE;
11: static PetscErrorCode MaybeLogStagePush(PetscLogStage stage) { return log_stages ? PetscLogStagePush(stage) : 0; }
12: static PetscErrorCode MaybeLogStagePop() { return log_stages ? PetscLogStagePop() : 0; }
13: PetscErrorCode elem_3d_elast_v_25(PetscScalar *);
15: int main(int argc,char **args)
16: {
17: Mat Amat;
19: PetscInt m,nn,M,Istart,Iend,i,j,k,ii,jj,kk,ic,ne=4,id;
20: PetscReal x,y,z,h,*coords,soft_alpha=1.e-3;
21: PetscBool two_solves=PETSC_FALSE,test_nonzero_cols=PETSC_FALSE,use_nearnullspace=PETSC_FALSE,test_late_bs=PETSC_FALSE;
22: Vec xx,bb;
23: KSP ksp;
24: MPI_Comm comm;
25: PetscMPIInt npe,mype;
26: PetscScalar DD[24][24],DD2[24][24];
27: PetscLogStage stage[6];
28: PetscScalar DD1[24][24];
30: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
31: comm = PETSC_COMM_WORLD;
32: MPI_Comm_rank(comm, &mype);
33: MPI_Comm_size(comm, &npe);
35: PetscOptionsBegin(comm,NULL,"3D bilinear Q1 elasticity options","");
36: {
37: char nestring[256];
38: PetscSNPrintf(nestring,sizeof nestring,"number of elements in each direction, ne+1 must be a multiple of %D (sizes^{1/3})",(PetscInt)(PetscPowReal((PetscReal)npe,1./3.) + .5));
39: PetscOptionsInt("-ne",nestring,"",ne,&ne,NULL);
40: PetscOptionsBool("-log_stages","Log stages of solve separately","",log_stages,&log_stages,NULL);
41: PetscOptionsReal("-alpha","material coefficient inside circle","",soft_alpha,&soft_alpha,NULL);
42: PetscOptionsBool("-two_solves","solve additional variant of the problem","",two_solves,&two_solves,NULL);
43: PetscOptionsBool("-test_nonzero_cols","nonzero test","",test_nonzero_cols,&test_nonzero_cols,NULL);
44: PetscOptionsBool("-use_mat_nearnullspace","MatNearNullSpace API test","",use_nearnullspace,&use_nearnullspace,NULL);
45: PetscOptionsBool("-test_late_bs","","",test_late_bs,&test_late_bs,NULL);
46: }
47: PetscOptionsEnd();
49: if (log_stages) {
50: PetscLogStageRegister("Setup", &stage[0]);
51: PetscLogStageRegister("Solve", &stage[1]);
52: PetscLogStageRegister("2nd Setup", &stage[2]);
53: PetscLogStageRegister("2nd Solve", &stage[3]);
54: PetscLogStageRegister("3rd Setup", &stage[4]);
55: PetscLogStageRegister("3rd Solve", &stage[5]);
56: } else {
57: for (i=0; i<(PetscInt)(sizeof(stage)/sizeof(stage[0])); i++) stage[i] = -1;
58: }
60: h = 1./ne; nn = ne+1;
61: /* ne*ne; number of global elements */
62: M = 3*nn*nn*nn; /* global number of equations */
63: if (npe==2) {
64: if (mype==1) m=0;
65: else m = nn*nn*nn;
66: npe = 1;
67: } else {
68: m = nn*nn*nn/npe;
69: if (mype==npe-1) m = nn*nn*nn - (npe-1)*m;
70: }
71: m *= 3; /* number of equations local*/
72: /* Setup solver */
73: KSPCreate(PETSC_COMM_WORLD,&ksp);
74: KSPSetComputeSingularValues(ksp, PETSC_TRUE);
75: KSPSetFromOptions(ksp);
76: {
77: /* configuration */
78: const PetscInt NP = (PetscInt)(PetscPowReal((PetscReal)npe,1./3.) + .5);
79: const PetscInt ipx = mype%NP, ipy = (mype%(NP*NP))/NP, ipz = mype/(NP*NP);
80: const PetscInt Ni0 = ipx*(nn/NP), Nj0 = ipy*(nn/NP), Nk0 = ipz*(nn/NP);
81: const PetscInt Ni1 = Ni0 + (m>0 ? (nn/NP) : 0), Nj1 = Nj0 + (nn/NP), Nk1 = Nk0 + (nn/NP);
82: const PetscInt NN = nn/NP, id0 = ipz*nn*nn*NN + ipy*nn*NN*NN + ipx*NN*NN*NN;
83: PetscInt *d_nnz, *o_nnz,osz[4]={0,9,15,19},nbc;
84: PetscScalar vv[24], v2[24];
85: if (npe!=NP*NP*NP) SETERRQ1(comm,PETSC_ERR_ARG_WRONG, "npe=%d: npe^{1/3} must be integer",npe);
86: if (nn!=NP*(nn/NP)) SETERRQ1(comm,PETSC_ERR_ARG_WRONG, "-ne %d: (ne+1)%(npe^{1/3}) must equal zero",ne);
88: /* count nnz */
89: PetscMalloc1(m+1, &d_nnz);
90: PetscMalloc1(m+1, &o_nnz);
91: for (i=Ni0,ic=0; i<Ni1; i++) {
92: for (j=Nj0; j<Nj1; j++) {
93: for (k=Nk0; k<Nk1; k++) {
94: nbc = 0;
95: if (i==Ni0 || i==Ni1-1) nbc++;
96: if (j==Nj0 || j==Nj1-1) nbc++;
97: if (k==Nk0 || k==Nk1-1) nbc++;
98: for (jj=0; jj<3; jj++,ic++) {
99: d_nnz[ic] = 3*(27-osz[nbc]);
100: o_nnz[ic] = 3*osz[nbc];
101: }
102: }
103: }
104: }
105: if (ic != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ic %D does not equal m %D",ic,m);
107: /* create stiffness matrix */
108: MatCreate(comm,&Amat);
109: MatSetSizes(Amat,m,m,M,M);
110: if (!test_late_bs) {
111: MatSetBlockSize(Amat,3);
112: }
113: MatSetType(Amat,MATAIJ);
114: MatSetFromOptions(Amat);
115: MatSeqAIJSetPreallocation(Amat,0,d_nnz);
116: MatMPIAIJSetPreallocation(Amat,0,d_nnz,0,o_nnz);
118: PetscFree(d_nnz);
119: PetscFree(o_nnz);
120: MatCreateVecs(Amat,&bb,&xx);
122: MatGetOwnershipRange(Amat,&Istart,&Iend);
124: if (m != Iend - Istart) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"m %D does not equal Iend %D - Istart %D",m,Iend,Istart);
125: /* generate element matrices */
126: {
127: PetscBool hasData = PETSC_TRUE;
128: if (!hasData) {
129: PetscPrintf(PETSC_COMM_WORLD,"\t No data is provided\n");
130: for (i=0; i<24; i++) {
131: for (j=0; j<24; j++) {
132: if (i==j) DD1[i][j] = 1.0;
133: else DD1[i][j] = -.25;
134: }
135: }
136: } else {
137: /* Get array data */
138: elem_3d_elast_v_25((PetscScalar*)DD1);
139: }
141: /* BC version of element */
142: for (i=0; i<24; i++) {
143: for (j=0; j<24; j++) {
144: if (i<12 || (j < 12 && !test_nonzero_cols)) {
145: if (i==j) DD2[i][j] = 0.1*DD1[i][j];
146: else DD2[i][j] = 0.0;
147: } else DD2[i][j] = DD1[i][j];
148: }
149: }
150: /* element residual/load vector */
151: for (i=0; i<24; i++) {
152: if (i%3==0) vv[i] = h*h;
153: else if (i%3==1) vv[i] = 2.0*h*h;
154: else vv[i] = .0;
155: }
156: for (i=0; i<24; i++) {
157: if (i%3==0 && i>=12) v2[i] = h*h;
158: else if (i%3==1 && i>=12) v2[i] = 2.0*h*h;
159: else v2[i] = .0;
160: }
161: }
163: PetscMalloc1(m+1, &coords);
164: coords[m] = -99.0;
166: /* forms the element stiffness and coordinates */
167: for (i=Ni0,ic=0,ii=0; i<Ni1; i++,ii++) {
168: for (j=Nj0,jj=0; j<Nj1; j++,jj++) {
169: for (k=Nk0,kk=0; k<Nk1; k++,kk++,ic++) {
170: /* coords */
171: x = coords[3*ic] = h*(PetscReal)i;
172: y = coords[3*ic+1] = h*(PetscReal)j;
173: z = coords[3*ic+2] = h*(PetscReal)k;
174: /* matrix */
175: id = id0 + ii + NN*jj + NN*NN*kk;
176: if (i<ne && j<ne && k<ne) {
177: /* radius */
178: 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));
179: PetscReal alpha = 1.0;
180: PetscInt jx,ix,idx[8],idx3[24];
181: idx[0] = id;
182: idx[1] = id+1;
183: idx[2] = id+NN+1;
184: idx[3] = id+NN;
185: idx[4] = id + NN*NN;
186: idx[5] = id+1 + NN*NN;
187: idx[6] = id+NN+1 + NN*NN;
188: idx[7] = id+NN + NN*NN;
190: /* correct indices */
191: if (i==Ni1-1 && Ni1!=nn) {
192: idx[1] += NN*(NN*NN-1);
193: idx[2] += NN*(NN*NN-1);
194: idx[5] += NN*(NN*NN-1);
195: idx[6] += NN*(NN*NN-1);
196: }
197: if (j==Nj1-1 && Nj1!=nn) {
198: idx[2] += NN*NN*(nn-1);
199: idx[3] += NN*NN*(nn-1);
200: idx[6] += NN*NN*(nn-1);
201: idx[7] += NN*NN*(nn-1);
202: }
203: if (k==Nk1-1 && Nk1!=nn) {
204: idx[4] += NN*(nn*nn-NN*NN);
205: idx[5] += NN*(nn*nn-NN*NN);
206: idx[6] += NN*(nn*nn-NN*NN);
207: idx[7] += NN*(nn*nn-NN*NN);
208: }
210: if (radius < 0.25) alpha = soft_alpha;
212: for (ix=0; ix<24; ix++) {
213: for (jx=0;jx<24;jx++) DD[ix][jx] = alpha*DD1[ix][jx];
214: }
215: if (k>0) {
216: if (!test_late_bs) {
217: MatSetValuesBlocked(Amat,8,idx,8,idx,(const PetscScalar*)DD,ADD_VALUES);
218: VecSetValuesBlocked(bb,8,idx,(const PetscScalar*)vv,ADD_VALUES);
219: } else {
220: for (ix=0; ix<8; ix++) { idx3[3*ix] = 3*idx[ix]; idx3[3*ix+1] = 3*idx[ix]+1; idx3[3*ix+2] = 3*idx[ix]+2; }
221: MatSetValues(Amat,24,idx3,24,idx3,(const PetscScalar*)DD,ADD_VALUES);
222: VecSetValues(bb,24,idx3,(const PetscScalar*)vv,ADD_VALUES);
223: }
224: } else {
225: /* a BC */
226: for (ix=0;ix<24;ix++) {
227: for (jx=0;jx<24;jx++) DD[ix][jx] = alpha*DD2[ix][jx];
228: }
229: if (!test_late_bs) {
230: MatSetValuesBlocked(Amat,8,idx,8,idx,(const PetscScalar*)DD,ADD_VALUES);
231: VecSetValuesBlocked(bb,8,idx,(const PetscScalar*)v2,ADD_VALUES);
232: } else {
233: for (ix=0; ix<8; ix++) { idx3[3*ix] = 3*idx[ix]; idx3[3*ix+1] = 3*idx[ix]+1; idx3[3*ix+2] = 3*idx[ix]+2; }
234: MatSetValues(Amat,24,idx3,24,idx3,(const PetscScalar*)DD,ADD_VALUES);
235: VecSetValues(bb,24,idx3,(const PetscScalar*)v2,ADD_VALUES);
236: }
237: }
238: }
239: }
240: }
242: }
243: MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
244: MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
245: VecAssemblyBegin(bb);
246: VecAssemblyEnd(bb);
247: }
248: MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
249: MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
250: VecAssemblyBegin(bb);
251: VecAssemblyEnd(bb);
252: if (test_late_bs) {
253: VecSetBlockSize(xx,3);
254: VecSetBlockSize(bb,3);
255: MatSetBlockSize(Amat,3);
256: }
258: if (!PETSC_TRUE) {
259: PetscViewer viewer;
260: PetscViewerASCIIOpen(comm, "Amat.m", &viewer);
261: PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
262: MatView(Amat,viewer);
263: PetscViewerPopFormat(viewer);
264: PetscViewerDestroy(&viewer);
265: }
267: /* finish KSP/PC setup */
268: KSPSetOperators(ksp, Amat, Amat);
269: if (use_nearnullspace) {
270: MatNullSpace matnull;
271: Vec vec_coords;
272: PetscScalar *c;
274: VecCreate(MPI_COMM_WORLD,&vec_coords);
275: VecSetBlockSize(vec_coords,3);
276: VecSetSizes(vec_coords,m,PETSC_DECIDE);
277: VecSetUp(vec_coords);
278: VecGetArray(vec_coords,&c);
279: for (i=0; i<m; i++) c[i] = coords[i]; /* Copy since Scalar type might be Complex */
280: VecRestoreArray(vec_coords,&c);
281: MatNullSpaceCreateRigidBody(vec_coords,&matnull);
282: MatSetNearNullSpace(Amat,matnull);
283: MatNullSpaceDestroy(&matnull);
284: VecDestroy(&vec_coords);
285: } else {
286: PC pc;
287: KSPGetPC(ksp,&pc);
288: PCSetCoordinates(pc, 3, m/3, coords);
289: }
291: MaybeLogStagePush(stage[0]);
293: /* PC setup basically */
294: KSPSetUp(ksp);
296: MaybeLogStagePop();
297: MaybeLogStagePush(stage[1]);
299: /* test BCs */
300: if (test_nonzero_cols) {
301: VecZeroEntries(xx);
302: if (mype==0) VecSetValue(xx,0,1.0,INSERT_VALUES);
303: VecAssemblyBegin(xx);
304: VecAssemblyEnd(xx);
305: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
306: }
308: /* 1st solve */
309: KSPSolve(ksp, bb, xx);
311: MaybeLogStagePop();
313: /* 2nd solve */
314: if (two_solves) {
315: PetscReal emax, emin;
316: PetscReal norm,norm2;
317: Vec res;
319: MaybeLogStagePush(stage[2]);
320: /* PC setup basically */
321: MatScale(Amat, 100000.0);
322: KSPSetOperators(ksp, Amat, Amat);
323: KSPSetUp(ksp);
325: MaybeLogStagePop();
326: MaybeLogStagePush(stage[3]);
327: KSPSolve(ksp, bb, xx);
328: KSPComputeExtremeSingularValues(ksp, &emax, &emin);
330: MaybeLogStagePop();
331: MaybeLogStagePush(stage[4]);
333: MaybeLogStagePop();
334: MaybeLogStagePush(stage[5]);
336: /* 3rd solve */
337: KSPSolve(ksp, bb, xx);
339: MaybeLogStagePop();
342: VecNorm(bb, NORM_2, &norm2);
344: VecDuplicate(xx, &res);
345: MatMult(Amat, xx, res);
346: VecAXPY(bb, -1.0, res);
347: VecDestroy(&res);
348: VecNorm(bb, NORM_2, &norm);
349: 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);
350: }
352: /* Free work space */
353: KSPDestroy(&ksp);
354: VecDestroy(&xx);
355: VecDestroy(&bb);
356: MatDestroy(&Amat);
357: PetscFree(coords);
359: PetscFinalize();
360: return ierr;
361: }
363: /* Data was previously provided in the file data/elem_3d_elast_v_25.tx */
364: PetscErrorCode elem_3d_elast_v_25(PetscScalar *dd)
365: {
367: PetscScalar DD[] = {
368: 0.18981481481481474 ,
369: 5.27777777777777568E-002,
370: 5.27777777777777568E-002,
371: -5.64814814814814659E-002,
372: -1.38888888888889072E-002,
373: -1.38888888888889089E-002,
374: -8.24074074074073876E-002,
375: -5.27777777777777429E-002,
376: 1.38888888888888725E-002,
377: 4.90740740740740339E-002,
378: 1.38888888888889124E-002,
379: 4.72222222222222071E-002,
380: 4.90740740740740339E-002,
381: 4.72222222222221932E-002,
382: 1.38888888888888968E-002,
383: -8.24074074074073876E-002,
384: 1.38888888888888673E-002,
385: -5.27777777777777429E-002,
386: -7.87037037037036785E-002,
387: -4.72222222222221932E-002,
388: -4.72222222222222071E-002,
389: 1.20370370370370180E-002,
390: -1.38888888888888742E-002,
391: -1.38888888888888829E-002,
392: 5.27777777777777568E-002,
393: 0.18981481481481474 ,
394: 5.27777777777777568E-002,
395: 1.38888888888889124E-002,
396: 4.90740740740740269E-002,
397: 4.72222222222221932E-002,
398: -5.27777777777777637E-002,
399: -8.24074074074073876E-002,
400: 1.38888888888888725E-002,
401: -1.38888888888889037E-002,
402: -5.64814814814814728E-002,
403: -1.38888888888888985E-002,
404: 4.72222222222221932E-002,
405: 4.90740740740740478E-002,
406: 1.38888888888888968E-002,
407: -1.38888888888888673E-002,
408: 1.20370370370370058E-002,
409: -1.38888888888888742E-002,
410: -4.72222222222221932E-002,
411: -7.87037037037036785E-002,
412: -4.72222222222222002E-002,
413: 1.38888888888888742E-002,
414: -8.24074074074073598E-002,
415: -5.27777777777777568E-002,
416: 5.27777777777777568E-002,
417: 5.27777777777777568E-002,
418: 0.18981481481481474 ,
419: 1.38888888888889055E-002,
420: 4.72222222222222002E-002,
421: 4.90740740740740269E-002,
422: -1.38888888888888829E-002,
423: -1.38888888888888829E-002,
424: 1.20370370370370180E-002,
425: 4.72222222222222002E-002,
426: 1.38888888888888985E-002,
427: 4.90740740740740339E-002,
428: -1.38888888888888985E-002,
429: -1.38888888888888968E-002,
430: -5.64814814814814520E-002,
431: -5.27777777777777568E-002,
432: 1.38888888888888777E-002,
433: -8.24074074074073876E-002,
434: -4.72222222222222002E-002,
435: -4.72222222222221932E-002,
436: -7.87037037037036646E-002,
437: 1.38888888888888794E-002,
438: -5.27777777777777568E-002,
439: -8.24074074074073598E-002,
440: -5.64814814814814659E-002,
441: 1.38888888888889124E-002,
442: 1.38888888888889055E-002,
443: 0.18981481481481474 ,
444: -5.27777777777777568E-002,
445: -5.27777777777777499E-002,
446: 4.90740740740740269E-002,
447: -1.38888888888889072E-002,
448: -4.72222222222221932E-002,
449: -8.24074074074073876E-002,
450: 5.27777777777777568E-002,
451: -1.38888888888888812E-002,
452: -8.24074074074073876E-002,
453: -1.38888888888888742E-002,
454: 5.27777777777777499E-002,
455: 4.90740740740740269E-002,
456: -4.72222222222221863E-002,
457: -1.38888888888889089E-002,
458: 1.20370370370370162E-002,
459: 1.38888888888888673E-002,
460: 1.38888888888888742E-002,
461: -7.87037037037036785E-002,
462: 4.72222222222222002E-002,
463: 4.72222222222222071E-002,
464: -1.38888888888889072E-002,
465: 4.90740740740740269E-002,
466: 4.72222222222222002E-002,
467: -5.27777777777777568E-002,
468: 0.18981481481481480 ,
469: 5.27777777777777568E-002,
470: 1.38888888888889020E-002,
471: -5.64814814814814728E-002,
472: -1.38888888888888951E-002,
473: 5.27777777777777637E-002,
474: -8.24074074074073876E-002,
475: 1.38888888888888881E-002,
476: 1.38888888888888742E-002,
477: 1.20370370370370232E-002,
478: -1.38888888888888812E-002,
479: -4.72222222222221863E-002,
480: 4.90740740740740339E-002,
481: 1.38888888888888933E-002,
482: -1.38888888888888812E-002,
483: -8.24074074074073876E-002,
484: -5.27777777777777568E-002,
485: 4.72222222222222071E-002,
486: -7.87037037037036924E-002,
487: -4.72222222222222140E-002,
488: -1.38888888888889089E-002,
489: 4.72222222222221932E-002,
490: 4.90740740740740269E-002,
491: -5.27777777777777499E-002,
492: 5.27777777777777568E-002,
493: 0.18981481481481477 ,
494: -4.72222222222222071E-002,
495: 1.38888888888888968E-002,
496: 4.90740740740740131E-002,
497: 1.38888888888888812E-002,
498: -1.38888888888888708E-002,
499: 1.20370370370370267E-002,
500: 5.27777777777777568E-002,
501: 1.38888888888888812E-002,
502: -8.24074074074073876E-002,
503: 1.38888888888889124E-002,
504: -1.38888888888889055E-002,
505: -5.64814814814814589E-002,
506: -1.38888888888888812E-002,
507: -5.27777777777777568E-002,
508: -8.24074074074073737E-002,
509: 4.72222222222222002E-002,
510: -4.72222222222222002E-002,
511: -7.87037037037036924E-002,
512: -8.24074074074073876E-002,
513: -5.27777777777777637E-002,
514: -1.38888888888888829E-002,
515: 4.90740740740740269E-002,
516: 1.38888888888889020E-002,
517: -4.72222222222222071E-002,
518: 0.18981481481481480 ,
519: 5.27777777777777637E-002,
520: -5.27777777777777637E-002,
521: -5.64814814814814728E-002,
522: -1.38888888888889037E-002,
523: 1.38888888888888951E-002,
524: -7.87037037037036785E-002,
525: -4.72222222222222002E-002,
526: 4.72222222222221932E-002,
527: 1.20370370370370128E-002,
528: -1.38888888888888725E-002,
529: 1.38888888888888812E-002,
530: 4.90740740740740408E-002,
531: 4.72222222222222002E-002,
532: -1.38888888888888951E-002,
533: -8.24074074074073876E-002,
534: 1.38888888888888812E-002,
535: 5.27777777777777637E-002,
536: -5.27777777777777429E-002,
537: -8.24074074074073876E-002,
538: -1.38888888888888829E-002,
539: -1.38888888888889072E-002,
540: -5.64814814814814728E-002,
541: 1.38888888888888968E-002,
542: 5.27777777777777637E-002,
543: 0.18981481481481480 ,
544: -5.27777777777777568E-002,
545: 1.38888888888888916E-002,
546: 4.90740740740740339E-002,
547: -4.72222222222222210E-002,
548: -4.72222222222221932E-002,
549: -7.87037037037036924E-002,
550: 4.72222222222222002E-002,
551: 1.38888888888888742E-002,
552: -8.24074074074073876E-002,
553: 5.27777777777777429E-002,
554: 4.72222222222222002E-002,
555: 4.90740740740740269E-002,
556: -1.38888888888888951E-002,
557: -1.38888888888888846E-002,
558: 1.20370370370370267E-002,
559: 1.38888888888888916E-002,
560: 1.38888888888888725E-002,
561: 1.38888888888888725E-002,
562: 1.20370370370370180E-002,
563: -4.72222222222221932E-002,
564: -1.38888888888888951E-002,
565: 4.90740740740740131E-002,
566: -5.27777777777777637E-002,
567: -5.27777777777777568E-002,
568: 0.18981481481481480 ,
569: -1.38888888888888968E-002,
570: -4.72222222222221932E-002,
571: 4.90740740740740339E-002,
572: 4.72222222222221932E-002,
573: 4.72222222222222071E-002,
574: -7.87037037037036646E-002,
575: -1.38888888888888742E-002,
576: 5.27777777777777499E-002,
577: -8.24074074074073737E-002,
578: 1.38888888888888933E-002,
579: 1.38888888888889020E-002,
580: -5.64814814814814589E-002,
581: 5.27777777777777568E-002,
582: -1.38888888888888794E-002,
583: -8.24074074074073876E-002,
584: 4.90740740740740339E-002,
585: -1.38888888888889037E-002,
586: 4.72222222222222002E-002,
587: -8.24074074074073876E-002,
588: 5.27777777777777637E-002,
589: 1.38888888888888812E-002,
590: -5.64814814814814728E-002,
591: 1.38888888888888916E-002,
592: -1.38888888888888968E-002,
593: 0.18981481481481480 ,
594: -5.27777777777777499E-002,
595: 5.27777777777777707E-002,
596: 1.20370370370370180E-002,
597: 1.38888888888888812E-002,
598: -1.38888888888888812E-002,
599: -7.87037037037036785E-002,
600: 4.72222222222222002E-002,
601: -4.72222222222222071E-002,
602: -8.24074074074073876E-002,
603: -1.38888888888888742E-002,
604: -5.27777777777777568E-002,
605: 4.90740740740740616E-002,
606: -4.72222222222222002E-002,
607: 1.38888888888888846E-002,
608: 1.38888888888889124E-002,
609: -5.64814814814814728E-002,
610: 1.38888888888888985E-002,
611: 5.27777777777777568E-002,
612: -8.24074074074073876E-002,
613: -1.38888888888888708E-002,
614: -1.38888888888889037E-002,
615: 4.90740740740740339E-002,
616: -4.72222222222221932E-002,
617: -5.27777777777777499E-002,
618: 0.18981481481481480 ,
619: -5.27777777777777568E-002,
620: -1.38888888888888673E-002,
621: -8.24074074074073598E-002,
622: 5.27777777777777429E-002,
623: 4.72222222222222002E-002,
624: -7.87037037037036785E-002,
625: 4.72222222222222002E-002,
626: 1.38888888888888708E-002,
627: 1.20370370370370128E-002,
628: 1.38888888888888760E-002,
629: -4.72222222222222002E-002,
630: 4.90740740740740478E-002,
631: -1.38888888888888951E-002,
632: 4.72222222222222071E-002,
633: -1.38888888888888985E-002,
634: 4.90740740740740339E-002,
635: -1.38888888888888812E-002,
636: 1.38888888888888881E-002,
637: 1.20370370370370267E-002,
638: 1.38888888888888951E-002,
639: -4.72222222222222210E-002,
640: 4.90740740740740339E-002,
641: 5.27777777777777707E-002,
642: -5.27777777777777568E-002,
643: 0.18981481481481477 ,
644: 1.38888888888888829E-002,
645: 5.27777777777777707E-002,
646: -8.24074074074073598E-002,
647: -4.72222222222222140E-002,
648: 4.72222222222222140E-002,
649: -7.87037037037036646E-002,
650: -5.27777777777777707E-002,
651: -1.38888888888888829E-002,
652: -8.24074074074073876E-002,
653: -1.38888888888888881E-002,
654: 1.38888888888888881E-002,
655: -5.64814814814814589E-002,
656: 4.90740740740740339E-002,
657: 4.72222222222221932E-002,
658: -1.38888888888888985E-002,
659: -8.24074074074073876E-002,
660: 1.38888888888888742E-002,
661: 5.27777777777777568E-002,
662: -7.87037037037036785E-002,
663: -4.72222222222221932E-002,
664: 4.72222222222221932E-002,
665: 1.20370370370370180E-002,
666: -1.38888888888888673E-002,
667: 1.38888888888888829E-002,
668: 0.18981481481481469 ,
669: 5.27777777777777429E-002,
670: -5.27777777777777429E-002,
671: -5.64814814814814659E-002,
672: -1.38888888888889055E-002,
673: 1.38888888888889055E-002,
674: -8.24074074074074153E-002,
675: -5.27777777777777429E-002,
676: -1.38888888888888760E-002,
677: 4.90740740740740408E-002,
678: 1.38888888888888968E-002,
679: -4.72222222222222071E-002,
680: 4.72222222222221932E-002,
681: 4.90740740740740478E-002,
682: -1.38888888888888968E-002,
683: -1.38888888888888742E-002,
684: 1.20370370370370232E-002,
685: 1.38888888888888812E-002,
686: -4.72222222222222002E-002,
687: -7.87037037037036924E-002,
688: 4.72222222222222071E-002,
689: 1.38888888888888812E-002,
690: -8.24074074074073598E-002,
691: 5.27777777777777707E-002,
692: 5.27777777777777429E-002,
693: 0.18981481481481477 ,
694: -5.27777777777777499E-002,
695: 1.38888888888889107E-002,
696: 4.90740740740740478E-002,
697: -4.72222222222221932E-002,
698: -5.27777777777777568E-002,
699: -8.24074074074074153E-002,
700: -1.38888888888888812E-002,
701: -1.38888888888888846E-002,
702: -5.64814814814814659E-002,
703: 1.38888888888888812E-002,
704: 1.38888888888888968E-002,
705: 1.38888888888888968E-002,
706: -5.64814814814814520E-002,
707: 5.27777777777777499E-002,
708: -1.38888888888888812E-002,
709: -8.24074074074073876E-002,
710: 4.72222222222221932E-002,
711: 4.72222222222222002E-002,
712: -7.87037037037036646E-002,
713: -1.38888888888888812E-002,
714: 5.27777777777777429E-002,
715: -8.24074074074073598E-002,
716: -5.27777777777777429E-002,
717: -5.27777777777777499E-002,
718: 0.18981481481481474 ,
719: -1.38888888888888985E-002,
720: -4.72222222222221863E-002,
721: 4.90740740740740339E-002,
722: 1.38888888888888829E-002,
723: 1.38888888888888777E-002,
724: 1.20370370370370249E-002,
725: -4.72222222222222002E-002,
726: -1.38888888888888933E-002,
727: 4.90740740740740339E-002,
728: -8.24074074074073876E-002,
729: -1.38888888888888673E-002,
730: -5.27777777777777568E-002,
731: 4.90740740740740269E-002,
732: -4.72222222222221863E-002,
733: 1.38888888888889124E-002,
734: 1.20370370370370128E-002,
735: 1.38888888888888742E-002,
736: -1.38888888888888742E-002,
737: -7.87037037037036785E-002,
738: 4.72222222222222002E-002,
739: -4.72222222222222140E-002,
740: -5.64814814814814659E-002,
741: 1.38888888888889107E-002,
742: -1.38888888888888985E-002,
743: 0.18981481481481474 ,
744: -5.27777777777777499E-002,
745: 5.27777777777777499E-002,
746: 4.90740740740740339E-002,
747: -1.38888888888889055E-002,
748: 4.72222222222221932E-002,
749: -8.24074074074074153E-002,
750: 5.27777777777777499E-002,
751: 1.38888888888888829E-002,
752: 1.38888888888888673E-002,
753: 1.20370370370370058E-002,
754: 1.38888888888888777E-002,
755: -4.72222222222221863E-002,
756: 4.90740740740740339E-002,
757: -1.38888888888889055E-002,
758: -1.38888888888888725E-002,
759: -8.24074074074073876E-002,
760: 5.27777777777777499E-002,
761: 4.72222222222222002E-002,
762: -7.87037037037036785E-002,
763: 4.72222222222222140E-002,
764: -1.38888888888889055E-002,
765: 4.90740740740740478E-002,
766: -4.72222222222221863E-002,
767: -5.27777777777777499E-002,
768: 0.18981481481481469 ,
769: -5.27777777777777499E-002,
770: 1.38888888888889072E-002,
771: -5.64814814814814659E-002,
772: 1.38888888888889003E-002,
773: 5.27777777777777429E-002,
774: -8.24074074074074153E-002,
775: -1.38888888888888812E-002,
776: -5.27777777777777429E-002,
777: -1.38888888888888742E-002,
778: -8.24074074074073876E-002,
779: -1.38888888888889089E-002,
780: 1.38888888888888933E-002,
781: -5.64814814814814589E-002,
782: 1.38888888888888812E-002,
783: 5.27777777777777429E-002,
784: -8.24074074074073737E-002,
785: -4.72222222222222071E-002,
786: 4.72222222222222002E-002,
787: -7.87037037037036646E-002,
788: 1.38888888888889055E-002,
789: -4.72222222222221932E-002,
790: 4.90740740740740339E-002,
791: 5.27777777777777499E-002,
792: -5.27777777777777499E-002,
793: 0.18981481481481474 ,
794: 4.72222222222222002E-002,
795: -1.38888888888888985E-002,
796: 4.90740740740740339E-002,
797: -1.38888888888888846E-002,
798: 1.38888888888888812E-002,
799: 1.20370370370370284E-002,
800: -7.87037037037036785E-002,
801: -4.72222222222221932E-002,
802: -4.72222222222222002E-002,
803: 1.20370370370370162E-002,
804: -1.38888888888888812E-002,
805: -1.38888888888888812E-002,
806: 4.90740740740740408E-002,
807: 4.72222222222222002E-002,
808: 1.38888888888888933E-002,
809: -8.24074074074073876E-002,
810: 1.38888888888888708E-002,
811: -5.27777777777777707E-002,
812: -8.24074074074074153E-002,
813: -5.27777777777777568E-002,
814: 1.38888888888888829E-002,
815: 4.90740740740740339E-002,
816: 1.38888888888889072E-002,
817: 4.72222222222222002E-002,
818: 0.18981481481481477 ,
819: 5.27777777777777429E-002,
820: 5.27777777777777568E-002,
821: -5.64814814814814659E-002,
822: -1.38888888888888846E-002,
823: -1.38888888888888881E-002,
824: -4.72222222222221932E-002,
825: -7.87037037037036785E-002,
826: -4.72222222222221932E-002,
827: 1.38888888888888673E-002,
828: -8.24074074074073876E-002,
829: -5.27777777777777568E-002,
830: 4.72222222222222002E-002,
831: 4.90740740740740269E-002,
832: 1.38888888888889020E-002,
833: -1.38888888888888742E-002,
834: 1.20370370370370128E-002,
835: -1.38888888888888829E-002,
836: -5.27777777777777429E-002,
837: -8.24074074074074153E-002,
838: 1.38888888888888777E-002,
839: -1.38888888888889055E-002,
840: -5.64814814814814659E-002,
841: -1.38888888888888985E-002,
842: 5.27777777777777429E-002,
843: 0.18981481481481469 ,
844: 5.27777777777777429E-002,
845: 1.38888888888888933E-002,
846: 4.90740740740740339E-002,
847: 4.72222222222222071E-002,
848: -4.72222222222222071E-002,
849: -4.72222222222222002E-002,
850: -7.87037037037036646E-002,
851: 1.38888888888888742E-002,
852: -5.27777777777777568E-002,
853: -8.24074074074073737E-002,
854: -1.38888888888888951E-002,
855: -1.38888888888888951E-002,
856: -5.64814814814814589E-002,
857: -5.27777777777777568E-002,
858: 1.38888888888888760E-002,
859: -8.24074074074073876E-002,
860: -1.38888888888888760E-002,
861: -1.38888888888888812E-002,
862: 1.20370370370370249E-002,
863: 4.72222222222221932E-002,
864: 1.38888888888889003E-002,
865: 4.90740740740740339E-002,
866: 5.27777777777777568E-002,
867: 5.27777777777777429E-002,
868: 0.18981481481481474 ,
869: 1.38888888888888933E-002,
870: 4.72222222222222071E-002,
871: 4.90740740740740339E-002,
872: 1.20370370370370180E-002,
873: 1.38888888888888742E-002,
874: 1.38888888888888794E-002,
875: -7.87037037037036785E-002,
876: 4.72222222222222071E-002,
877: 4.72222222222222002E-002,
878: -8.24074074074073876E-002,
879: -1.38888888888888846E-002,
880: 5.27777777777777568E-002,
881: 4.90740740740740616E-002,
882: -4.72222222222222002E-002,
883: -1.38888888888888881E-002,
884: 4.90740740740740408E-002,
885: -1.38888888888888846E-002,
886: -4.72222222222222002E-002,
887: -8.24074074074074153E-002,
888: 5.27777777777777429E-002,
889: -1.38888888888888846E-002,
890: -5.64814814814814659E-002,
891: 1.38888888888888933E-002,
892: 1.38888888888888933E-002,
893: 0.18981481481481477 ,
894: -5.27777777777777568E-002,
895: -5.27777777777777637E-002,
896: -1.38888888888888742E-002,
897: -8.24074074074073598E-002,
898: -5.27777777777777568E-002,
899: 4.72222222222222002E-002,
900: -7.87037037037036924E-002,
901: -4.72222222222222002E-002,
902: 1.38888888888888812E-002,
903: 1.20370370370370267E-002,
904: -1.38888888888888794E-002,
905: -4.72222222222222002E-002,
906: 4.90740740740740478E-002,
907: 1.38888888888888881E-002,
908: 1.38888888888888968E-002,
909: -5.64814814814814659E-002,
910: -1.38888888888888933E-002,
911: 5.27777777777777499E-002,
912: -8.24074074074074153E-002,
913: 1.38888888888888812E-002,
914: -1.38888888888888846E-002,
915: 4.90740740740740339E-002,
916: 4.72222222222222071E-002,
917: -5.27777777777777568E-002,
918: 0.18981481481481477 ,
919: 5.27777777777777637E-002,
920: -1.38888888888888829E-002,
921: -5.27777777777777568E-002,
922: -8.24074074074073598E-002,
923: 4.72222222222222071E-002,
924: -4.72222222222222140E-002,
925: -7.87037037037036924E-002,
926: 5.27777777777777637E-002,
927: 1.38888888888888916E-002,
928: -8.24074074074073876E-002,
929: 1.38888888888888846E-002,
930: -1.38888888888888951E-002,
931: -5.64814814814814589E-002,
932: -4.72222222222222071E-002,
933: 1.38888888888888812E-002,
934: 4.90740740740740339E-002,
935: 1.38888888888888829E-002,
936: -1.38888888888888812E-002,
937: 1.20370370370370284E-002,
938: -1.38888888888888881E-002,
939: 4.72222222222222071E-002,
940: 4.90740740740740339E-002,
941: -5.27777777777777637E-002,
942: 5.27777777777777637E-002,
943: 0.18981481481481477 ,
944: };
947: PetscMemcpy(dd,DD,sizeof(PetscScalar)*576);
948: return(0);
949: }
952: /*TEST
954: test:
955: nsize: 8
956: 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 -ksp_view -use_mat_nearnullspace -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -pc_gamg_square_graph 1 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -gamg_est_ksp_type cg -gamg_est_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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short
957: filter: grep -v variant
959: test:
960: suffix: 2
961: nsize: 8
962: args: -ne 31 -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_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -pc_gamg_square_graph 1 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -gamg_est_ksp_type cg -gamg_est_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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg
963: filter: grep -v variant
965: test:
966: suffix: latebs
967: filter: grep -v variant
968: nsize: 8
969: 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 -ksp_view -use_mat_nearnullspace -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -pc_gamg_square_graph 1 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -gamg_est_ksp_type cg -gamg_est_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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short -ksp_view
971: test:
972: suffix: latebs-2
973: filter: grep -v variant
974: nsize: 8
975: 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 -ksp_view -use_mat_nearnullspace -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -pc_gamg_square_graph 1 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -gamg_est_ksp_type cg -gamg_est_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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short -ksp_view
977: test:
978: suffix: ml
979: nsize: 8
980: 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 -mg_levels_esteig_ksp_type cg
981: requires: ml
983: test:
984: suffix: nns
985: args: -ne 9 -alpha 1.e-3 -ksp_converged_reason -ksp_type cg -ksp_max_it 50 -pc_type gamg -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 -mg_levels_esteig_ksp_type cg
987: test:
988: suffix: nns_telescope
989: nsize: 2
990: args: -use_mat_nearnullspace -ksp_monitor_short -pc_type telescope -pc_telescope_reduction_factor 2 -telescope_pc_type gamg
992: test:
993: suffix: seqaijmkl
994: nsize: 8
995: requires: mkl_sparse
996: 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_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -pc_gamg_square_graph 1 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_pc_type jacobi -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.05 -gamg_est_ksp_type cg -gamg_est_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
998: TEST*/