Actual source code: ex56.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  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:     MatSetOption(Amat,MAT_SPD,PETSC_TRUE);
115:     MatSetFromOptions(Amat);
116:     MatSeqAIJSetPreallocation(Amat,0,d_nnz);
117:     MatMPIAIJSetPreallocation(Amat,0,d_nnz,0,o_nnz);

119:     PetscFree(d_nnz);
120:     PetscFree(o_nnz);
121:     MatCreateVecs(Amat,&bb,&xx);

123:     MatGetOwnershipRange(Amat,&Istart,&Iend);

125:     if (m != Iend - Istart) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"m %D does not equal Iend %D - Istart %D",m,Iend,Istart);
126:     /* generate element matrices */
127:     {
128:       PetscBool hasData = PETSC_TRUE;
129:       if (!hasData) {
130:         PetscPrintf(PETSC_COMM_WORLD,"\t No data is provided\n");
131:         for (i=0; i<24; i++) {
132:           for (j=0; j<24; j++) {
133:             if (i==j) DD1[i][j] = 1.0;
134:             else DD1[i][j] = -.25;
135:           }
136:         }
137:       } else {
138:         /* Get array data */
139:         elem_3d_elast_v_25((PetscScalar*)DD1);
140:       }

142:       /* BC version of element */
143:       for (i=0; i<24; i++) {
144:         for (j=0; j<24; j++) {
145:           if (i<12 || (j < 12 && !test_nonzero_cols)) {
146:             if (i==j) DD2[i][j] = 0.1*DD1[i][j];
147:             else DD2[i][j] = 0.0;
148:           } else DD2[i][j] = DD1[i][j];
149:         }
150:       }
151:       /* element residual/load vector */
152:       for (i=0; i<24; i++) {
153:         if (i%3==0) vv[i] = h*h;
154:         else if (i%3==1) vv[i] = 2.0*h*h;
155:         else vv[i] = .0;
156:       }
157:       for (i=0; i<24; i++) {
158:         if (i%3==0 && i>=12) v2[i] = h*h;
159:         else if (i%3==1 && i>=12) v2[i] = 2.0*h*h;
160:         else v2[i] = .0;
161:       }
162:     }

164:     PetscMalloc1(m+1, &coords);
165:     coords[m] = -99.0;

167:     /* forms the element stiffness and coordinates */
168:     for (i=Ni0,ic=0,ii=0; i<Ni1; i++,ii++) {
169:       for (j=Nj0,jj=0; j<Nj1; j++,jj++) {
170:         for (k=Nk0,kk=0; k<Nk1; k++,kk++,ic++) {
171:           /* coords */
172:           x = coords[3*ic] = h*(PetscReal)i;
173:           y = coords[3*ic+1] = h*(PetscReal)j;
174:           z = coords[3*ic+2] = h*(PetscReal)k;
175:           /* matrix */
176:           id = id0 + ii + NN*jj + NN*NN*kk;
177:           if (i<ne && j<ne && k<ne) {
178:             /* radius */
179:             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));
180:             PetscReal alpha = 1.0;
181:             PetscInt  jx,ix,idx[8],idx3[24];
182:             idx[0] = id;
183:             idx[1] = id+1;
184:             idx[2] = id+NN+1;
185:             idx[3] = id+NN;
186:             idx[4] = id + NN*NN;
187:             idx[5] = id+1 + NN*NN;
188:             idx[6] = id+NN+1 + NN*NN;
189:             idx[7] = id+NN + NN*NN;

191:             /* correct indices */
192:             if (i==Ni1-1 && Ni1!=nn) {
193:               idx[1] += NN*(NN*NN-1);
194:               idx[2] += NN*(NN*NN-1);
195:               idx[5] += NN*(NN*NN-1);
196:               idx[6] += NN*(NN*NN-1);
197:             }
198:             if (j==Nj1-1 && Nj1!=nn) {
199:               idx[2] += NN*NN*(nn-1);
200:               idx[3] += NN*NN*(nn-1);
201:               idx[6] += NN*NN*(nn-1);
202:               idx[7] += NN*NN*(nn-1);
203:             }
204:             if (k==Nk1-1 && Nk1!=nn) {
205:               idx[4] += NN*(nn*nn-NN*NN);
206:               idx[5] += NN*(nn*nn-NN*NN);
207:               idx[6] += NN*(nn*nn-NN*NN);
208:               idx[7] += NN*(nn*nn-NN*NN);
209:             }

211:             if (radius < 0.25) alpha = soft_alpha;

213:             for (ix=0; ix<24; ix++) {
214:               for (jx=0;jx<24;jx++) DD[ix][jx] = alpha*DD1[ix][jx];
215:             }
216:             if (k>0) {
217:               if (!test_late_bs) {
218:                 MatSetValuesBlocked(Amat,8,idx,8,idx,(const PetscScalar*)DD,ADD_VALUES);
219:                 VecSetValuesBlocked(bb,8,idx,(const PetscScalar*)vv,ADD_VALUES);
220:               } else {
221:                 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; }
222:                 MatSetValues(Amat,24,idx3,24,idx3,(const PetscScalar*)DD,ADD_VALUES);
223:                 VecSetValues(bb,24,idx3,(const PetscScalar*)vv,ADD_VALUES);
224:               }
225:             } else {
226:               /* a BC */
227:               for (ix=0;ix<24;ix++) {
228:                 for (jx=0;jx<24;jx++) DD[ix][jx] = alpha*DD2[ix][jx];
229:               }
230:               if (!test_late_bs) {
231:                 MatSetValuesBlocked(Amat,8,idx,8,idx,(const PetscScalar*)DD,ADD_VALUES);
232:                 VecSetValuesBlocked(bb,8,idx,(const PetscScalar*)v2,ADD_VALUES);
233:               } else {
234:                 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; }
235:                 MatSetValues(Amat,24,idx3,24,idx3,(const PetscScalar*)DD,ADD_VALUES);
236:                 VecSetValues(bb,24,idx3,(const PetscScalar*)v2,ADD_VALUES);
237:               }
238:             }
239:           }
240:         }
241:       }

243:     }
244:     MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
245:     MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
246:     VecAssemblyBegin(bb);
247:     VecAssemblyEnd(bb);
248:   }
249:   MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
250:   MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
251:   VecAssemblyBegin(bb);
252:   VecAssemblyEnd(bb);
253:   if (test_late_bs) {
254:     VecSetBlockSize(xx,3);
255:     VecSetBlockSize(bb,3);
256:     MatSetBlockSize(Amat,3);
257:   }

259:   if (!PETSC_TRUE) {
260:     PetscViewer viewer;
261:     PetscViewerASCIIOpen(comm, "Amat.m", &viewer);
262:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
263:     MatView(Amat,viewer);
264:     PetscViewerPopFormat(viewer);
265:     PetscViewerDestroy(&viewer);
266:   }

268:   /* finish KSP/PC setup */
269:   KSPSetOperators(ksp, Amat, Amat);
270:   if (use_nearnullspace) {
271:     MatNullSpace matnull;
272:     Vec          vec_coords;
273:     PetscScalar  *c;

275:     VecCreate(MPI_COMM_WORLD,&vec_coords);
276:     VecSetBlockSize(vec_coords,3);
277:     VecSetSizes(vec_coords,m,PETSC_DECIDE);
278:     VecSetUp(vec_coords);
279:     VecGetArray(vec_coords,&c);
280:     for (i=0; i<m; i++) c[i] = coords[i]; /* Copy since Scalar type might be Complex */
281:     VecRestoreArray(vec_coords,&c);
282:     MatNullSpaceCreateRigidBody(vec_coords,&matnull);
283:     MatSetNearNullSpace(Amat,matnull);
284:     MatNullSpaceDestroy(&matnull);
285:     VecDestroy(&vec_coords);
286:   } else {
287:     PC             pc;
288:     KSPGetPC(ksp,&pc);
289:     PCSetCoordinates(pc, 3, m/3, coords);
290:   }

292:   MaybeLogStagePush(stage[0]);

294:   /* PC setup basically */
295:   KSPSetUp(ksp);

297:   MaybeLogStagePop();
298:   MaybeLogStagePush(stage[1]);

300:   /* test BCs */
301:   if (test_nonzero_cols) {
302:     VecZeroEntries(xx);
303:     if (mype==0) VecSetValue(xx,0,1.0,INSERT_VALUES);
304:     VecAssemblyBegin(xx);
305:     VecAssemblyEnd(xx);
306:     KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
307:   }

309:   /* 1st solve */
310:   KSPSolve(ksp, bb, xx);

312:   MaybeLogStagePop();

314:   /* 2nd solve */
315:   if (two_solves) {
316:     PetscReal emax, emin;
317:     PetscReal norm,norm2;
318:     Vec       res;

320:     MaybeLogStagePush(stage[2]);
321:     /* PC setup basically */
322:     MatScale(Amat, 100000.0);
323:     KSPSetOperators(ksp, Amat, Amat);
324:     KSPSetUp(ksp);

326:     MaybeLogStagePop();
327:     MaybeLogStagePush(stage[3]);
328:     KSPSolve(ksp, bb, xx);
329:     KSPComputeExtremeSingularValues(ksp, &emax, &emin);

331:     MaybeLogStagePop();
332:     MaybeLogStagePush(stage[4]);

334:     MaybeLogStagePop();
335:     MaybeLogStagePush(stage[5]);

337:     /* 3rd solve */
338:     KSPSolve(ksp, bb, xx);

340:     MaybeLogStagePop();


343:     VecNorm(bb, NORM_2, &norm2);

345:     VecDuplicate(xx, &res);
346:     MatMult(Amat, xx, res);
347:     VecAXPY(bb, -1.0, res);
348:     VecDestroy(&res);
349:     VecNorm(bb, NORM_2, &norm);
350:     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);
351:   }

353:   /* Free work space */
354:   KSPDestroy(&ksp);
355:   VecDestroy(&xx);
356:   VecDestroy(&bb);
357:   MatDestroy(&Amat);
358:   PetscFree(coords);

360:   PetscFinalize();
361:   return ierr;
362: }

364: /* Data was previously provided in the file data/elem_3d_elast_v_25.tx */
365: PetscErrorCode elem_3d_elast_v_25(PetscScalar *dd)
366: {
368:   PetscScalar    DD[] = {
369:   0.18981481481481474     ,
370:   5.27777777777777568E-002,
371:   5.27777777777777568E-002,
372:  -5.64814814814814659E-002,
373:  -1.38888888888889072E-002,
374:  -1.38888888888889089E-002,
375:  -8.24074074074073876E-002,
376:  -5.27777777777777429E-002,
377:   1.38888888888888725E-002,
378:   4.90740740740740339E-002,
379:   1.38888888888889124E-002,
380:   4.72222222222222071E-002,
381:   4.90740740740740339E-002,
382:   4.72222222222221932E-002,
383:   1.38888888888888968E-002,
384:  -8.24074074074073876E-002,
385:   1.38888888888888673E-002,
386:  -5.27777777777777429E-002,
387:  -7.87037037037036785E-002,
388:  -4.72222222222221932E-002,
389:  -4.72222222222222071E-002,
390:   1.20370370370370180E-002,
391:  -1.38888888888888742E-002,
392:  -1.38888888888888829E-002,
393:   5.27777777777777568E-002,
394:   0.18981481481481474     ,
395:   5.27777777777777568E-002,
396:   1.38888888888889124E-002,
397:   4.90740740740740269E-002,
398:   4.72222222222221932E-002,
399:  -5.27777777777777637E-002,
400:  -8.24074074074073876E-002,
401:   1.38888888888888725E-002,
402:  -1.38888888888889037E-002,
403:  -5.64814814814814728E-002,
404:  -1.38888888888888985E-002,
405:   4.72222222222221932E-002,
406:   4.90740740740740478E-002,
407:   1.38888888888888968E-002,
408:  -1.38888888888888673E-002,
409:   1.20370370370370058E-002,
410:  -1.38888888888888742E-002,
411:  -4.72222222222221932E-002,
412:  -7.87037037037036785E-002,
413:  -4.72222222222222002E-002,
414:   1.38888888888888742E-002,
415:  -8.24074074074073598E-002,
416:  -5.27777777777777568E-002,
417:   5.27777777777777568E-002,
418:   5.27777777777777568E-002,
419:   0.18981481481481474     ,
420:   1.38888888888889055E-002,
421:   4.72222222222222002E-002,
422:   4.90740740740740269E-002,
423:  -1.38888888888888829E-002,
424:  -1.38888888888888829E-002,
425:   1.20370370370370180E-002,
426:   4.72222222222222002E-002,
427:   1.38888888888888985E-002,
428:   4.90740740740740339E-002,
429:  -1.38888888888888985E-002,
430:  -1.38888888888888968E-002,
431:  -5.64814814814814520E-002,
432:  -5.27777777777777568E-002,
433:   1.38888888888888777E-002,
434:  -8.24074074074073876E-002,
435:  -4.72222222222222002E-002,
436:  -4.72222222222221932E-002,
437:  -7.87037037037036646E-002,
438:   1.38888888888888794E-002,
439:  -5.27777777777777568E-002,
440:  -8.24074074074073598E-002,
441:  -5.64814814814814659E-002,
442:   1.38888888888889124E-002,
443:   1.38888888888889055E-002,
444:   0.18981481481481474     ,
445:  -5.27777777777777568E-002,
446:  -5.27777777777777499E-002,
447:   4.90740740740740269E-002,
448:  -1.38888888888889072E-002,
449:  -4.72222222222221932E-002,
450:  -8.24074074074073876E-002,
451:   5.27777777777777568E-002,
452:  -1.38888888888888812E-002,
453:  -8.24074074074073876E-002,
454:  -1.38888888888888742E-002,
455:   5.27777777777777499E-002,
456:   4.90740740740740269E-002,
457:  -4.72222222222221863E-002,
458:  -1.38888888888889089E-002,
459:   1.20370370370370162E-002,
460:   1.38888888888888673E-002,
461:   1.38888888888888742E-002,
462:  -7.87037037037036785E-002,
463:   4.72222222222222002E-002,
464:   4.72222222222222071E-002,
465:  -1.38888888888889072E-002,
466:   4.90740740740740269E-002,
467:   4.72222222222222002E-002,
468:  -5.27777777777777568E-002,
469:   0.18981481481481480     ,
470:   5.27777777777777568E-002,
471:   1.38888888888889020E-002,
472:  -5.64814814814814728E-002,
473:  -1.38888888888888951E-002,
474:   5.27777777777777637E-002,
475:  -8.24074074074073876E-002,
476:   1.38888888888888881E-002,
477:   1.38888888888888742E-002,
478:   1.20370370370370232E-002,
479:  -1.38888888888888812E-002,
480:  -4.72222222222221863E-002,
481:   4.90740740740740339E-002,
482:   1.38888888888888933E-002,
483:  -1.38888888888888812E-002,
484:  -8.24074074074073876E-002,
485:  -5.27777777777777568E-002,
486:   4.72222222222222071E-002,
487:  -7.87037037037036924E-002,
488:  -4.72222222222222140E-002,
489:  -1.38888888888889089E-002,
490:   4.72222222222221932E-002,
491:   4.90740740740740269E-002,
492:  -5.27777777777777499E-002,
493:   5.27777777777777568E-002,
494:   0.18981481481481477     ,
495:  -4.72222222222222071E-002,
496:   1.38888888888888968E-002,
497:   4.90740740740740131E-002,
498:   1.38888888888888812E-002,
499:  -1.38888888888888708E-002,
500:   1.20370370370370267E-002,
501:   5.27777777777777568E-002,
502:   1.38888888888888812E-002,
503:  -8.24074074074073876E-002,
504:   1.38888888888889124E-002,
505:  -1.38888888888889055E-002,
506:  -5.64814814814814589E-002,
507:  -1.38888888888888812E-002,
508:  -5.27777777777777568E-002,
509:  -8.24074074074073737E-002,
510:   4.72222222222222002E-002,
511:  -4.72222222222222002E-002,
512:  -7.87037037037036924E-002,
513:  -8.24074074074073876E-002,
514:  -5.27777777777777637E-002,
515:  -1.38888888888888829E-002,
516:   4.90740740740740269E-002,
517:   1.38888888888889020E-002,
518:  -4.72222222222222071E-002,
519:   0.18981481481481480     ,
520:   5.27777777777777637E-002,
521:  -5.27777777777777637E-002,
522:  -5.64814814814814728E-002,
523:  -1.38888888888889037E-002,
524:   1.38888888888888951E-002,
525:  -7.87037037037036785E-002,
526:  -4.72222222222222002E-002,
527:   4.72222222222221932E-002,
528:   1.20370370370370128E-002,
529:  -1.38888888888888725E-002,
530:   1.38888888888888812E-002,
531:   4.90740740740740408E-002,
532:   4.72222222222222002E-002,
533:  -1.38888888888888951E-002,
534:  -8.24074074074073876E-002,
535:   1.38888888888888812E-002,
536:   5.27777777777777637E-002,
537:  -5.27777777777777429E-002,
538:  -8.24074074074073876E-002,
539:  -1.38888888888888829E-002,
540:  -1.38888888888889072E-002,
541:  -5.64814814814814728E-002,
542:   1.38888888888888968E-002,
543:   5.27777777777777637E-002,
544:   0.18981481481481480     ,
545:  -5.27777777777777568E-002,
546:   1.38888888888888916E-002,
547:   4.90740740740740339E-002,
548:  -4.72222222222222210E-002,
549:  -4.72222222222221932E-002,
550:  -7.87037037037036924E-002,
551:   4.72222222222222002E-002,
552:   1.38888888888888742E-002,
553:  -8.24074074074073876E-002,
554:   5.27777777777777429E-002,
555:   4.72222222222222002E-002,
556:   4.90740740740740269E-002,
557:  -1.38888888888888951E-002,
558:  -1.38888888888888846E-002,
559:   1.20370370370370267E-002,
560:   1.38888888888888916E-002,
561:   1.38888888888888725E-002,
562:   1.38888888888888725E-002,
563:   1.20370370370370180E-002,
564:  -4.72222222222221932E-002,
565:  -1.38888888888888951E-002,
566:   4.90740740740740131E-002,
567:  -5.27777777777777637E-002,
568:  -5.27777777777777568E-002,
569:   0.18981481481481480     ,
570:  -1.38888888888888968E-002,
571:  -4.72222222222221932E-002,
572:   4.90740740740740339E-002,
573:   4.72222222222221932E-002,
574:   4.72222222222222071E-002,
575:  -7.87037037037036646E-002,
576:  -1.38888888888888742E-002,
577:   5.27777777777777499E-002,
578:  -8.24074074074073737E-002,
579:   1.38888888888888933E-002,
580:   1.38888888888889020E-002,
581:  -5.64814814814814589E-002,
582:   5.27777777777777568E-002,
583:  -1.38888888888888794E-002,
584:  -8.24074074074073876E-002,
585:   4.90740740740740339E-002,
586:  -1.38888888888889037E-002,
587:   4.72222222222222002E-002,
588:  -8.24074074074073876E-002,
589:   5.27777777777777637E-002,
590:   1.38888888888888812E-002,
591:  -5.64814814814814728E-002,
592:   1.38888888888888916E-002,
593:  -1.38888888888888968E-002,
594:   0.18981481481481480     ,
595:  -5.27777777777777499E-002,
596:   5.27777777777777707E-002,
597:   1.20370370370370180E-002,
598:   1.38888888888888812E-002,
599:  -1.38888888888888812E-002,
600:  -7.87037037037036785E-002,
601:   4.72222222222222002E-002,
602:  -4.72222222222222071E-002,
603:  -8.24074074074073876E-002,
604:  -1.38888888888888742E-002,
605:  -5.27777777777777568E-002,
606:   4.90740740740740616E-002,
607:  -4.72222222222222002E-002,
608:   1.38888888888888846E-002,
609:   1.38888888888889124E-002,
610:  -5.64814814814814728E-002,
611:   1.38888888888888985E-002,
612:   5.27777777777777568E-002,
613:  -8.24074074074073876E-002,
614:  -1.38888888888888708E-002,
615:  -1.38888888888889037E-002,
616:   4.90740740740740339E-002,
617:  -4.72222222222221932E-002,
618:  -5.27777777777777499E-002,
619:   0.18981481481481480     ,
620:  -5.27777777777777568E-002,
621:  -1.38888888888888673E-002,
622:  -8.24074074074073598E-002,
623:   5.27777777777777429E-002,
624:   4.72222222222222002E-002,
625:  -7.87037037037036785E-002,
626:   4.72222222222222002E-002,
627:   1.38888888888888708E-002,
628:   1.20370370370370128E-002,
629:   1.38888888888888760E-002,
630:  -4.72222222222222002E-002,
631:   4.90740740740740478E-002,
632:  -1.38888888888888951E-002,
633:   4.72222222222222071E-002,
634:  -1.38888888888888985E-002,
635:   4.90740740740740339E-002,
636:  -1.38888888888888812E-002,
637:   1.38888888888888881E-002,
638:   1.20370370370370267E-002,
639:   1.38888888888888951E-002,
640:  -4.72222222222222210E-002,
641:   4.90740740740740339E-002,
642:   5.27777777777777707E-002,
643:  -5.27777777777777568E-002,
644:   0.18981481481481477     ,
645:   1.38888888888888829E-002,
646:   5.27777777777777707E-002,
647:  -8.24074074074073598E-002,
648:  -4.72222222222222140E-002,
649:   4.72222222222222140E-002,
650:  -7.87037037037036646E-002,
651:  -5.27777777777777707E-002,
652:  -1.38888888888888829E-002,
653:  -8.24074074074073876E-002,
654:  -1.38888888888888881E-002,
655:   1.38888888888888881E-002,
656:  -5.64814814814814589E-002,
657:   4.90740740740740339E-002,
658:   4.72222222222221932E-002,
659:  -1.38888888888888985E-002,
660:  -8.24074074074073876E-002,
661:   1.38888888888888742E-002,
662:   5.27777777777777568E-002,
663:  -7.87037037037036785E-002,
664:  -4.72222222222221932E-002,
665:   4.72222222222221932E-002,
666:   1.20370370370370180E-002,
667:  -1.38888888888888673E-002,
668:   1.38888888888888829E-002,
669:   0.18981481481481469     ,
670:   5.27777777777777429E-002,
671:  -5.27777777777777429E-002,
672:  -5.64814814814814659E-002,
673:  -1.38888888888889055E-002,
674:   1.38888888888889055E-002,
675:  -8.24074074074074153E-002,
676:  -5.27777777777777429E-002,
677:  -1.38888888888888760E-002,
678:   4.90740740740740408E-002,
679:   1.38888888888888968E-002,
680:  -4.72222222222222071E-002,
681:   4.72222222222221932E-002,
682:   4.90740740740740478E-002,
683:  -1.38888888888888968E-002,
684:  -1.38888888888888742E-002,
685:   1.20370370370370232E-002,
686:   1.38888888888888812E-002,
687:  -4.72222222222222002E-002,
688:  -7.87037037037036924E-002,
689:   4.72222222222222071E-002,
690:   1.38888888888888812E-002,
691:  -8.24074074074073598E-002,
692:   5.27777777777777707E-002,
693:   5.27777777777777429E-002,
694:   0.18981481481481477     ,
695:  -5.27777777777777499E-002,
696:   1.38888888888889107E-002,
697:   4.90740740740740478E-002,
698:  -4.72222222222221932E-002,
699:  -5.27777777777777568E-002,
700:  -8.24074074074074153E-002,
701:  -1.38888888888888812E-002,
702:  -1.38888888888888846E-002,
703:  -5.64814814814814659E-002,
704:   1.38888888888888812E-002,
705:   1.38888888888888968E-002,
706:   1.38888888888888968E-002,
707:  -5.64814814814814520E-002,
708:   5.27777777777777499E-002,
709:  -1.38888888888888812E-002,
710:  -8.24074074074073876E-002,
711:   4.72222222222221932E-002,
712:   4.72222222222222002E-002,
713:  -7.87037037037036646E-002,
714:  -1.38888888888888812E-002,
715:   5.27777777777777429E-002,
716:  -8.24074074074073598E-002,
717:  -5.27777777777777429E-002,
718:  -5.27777777777777499E-002,
719:   0.18981481481481474     ,
720:  -1.38888888888888985E-002,
721:  -4.72222222222221863E-002,
722:   4.90740740740740339E-002,
723:   1.38888888888888829E-002,
724:   1.38888888888888777E-002,
725:   1.20370370370370249E-002,
726:  -4.72222222222222002E-002,
727:  -1.38888888888888933E-002,
728:   4.90740740740740339E-002,
729:  -8.24074074074073876E-002,
730:  -1.38888888888888673E-002,
731:  -5.27777777777777568E-002,
732:   4.90740740740740269E-002,
733:  -4.72222222222221863E-002,
734:   1.38888888888889124E-002,
735:   1.20370370370370128E-002,
736:   1.38888888888888742E-002,
737:  -1.38888888888888742E-002,
738:  -7.87037037037036785E-002,
739:   4.72222222222222002E-002,
740:  -4.72222222222222140E-002,
741:  -5.64814814814814659E-002,
742:   1.38888888888889107E-002,
743:  -1.38888888888888985E-002,
744:   0.18981481481481474     ,
745:  -5.27777777777777499E-002,
746:   5.27777777777777499E-002,
747:   4.90740740740740339E-002,
748:  -1.38888888888889055E-002,
749:   4.72222222222221932E-002,
750:  -8.24074074074074153E-002,
751:   5.27777777777777499E-002,
752:   1.38888888888888829E-002,
753:   1.38888888888888673E-002,
754:   1.20370370370370058E-002,
755:   1.38888888888888777E-002,
756:  -4.72222222222221863E-002,
757:   4.90740740740740339E-002,
758:  -1.38888888888889055E-002,
759:  -1.38888888888888725E-002,
760:  -8.24074074074073876E-002,
761:   5.27777777777777499E-002,
762:   4.72222222222222002E-002,
763:  -7.87037037037036785E-002,
764:   4.72222222222222140E-002,
765:  -1.38888888888889055E-002,
766:   4.90740740740740478E-002,
767:  -4.72222222222221863E-002,
768:  -5.27777777777777499E-002,
769:   0.18981481481481469     ,
770:  -5.27777777777777499E-002,
771:   1.38888888888889072E-002,
772:  -5.64814814814814659E-002,
773:   1.38888888888889003E-002,
774:   5.27777777777777429E-002,
775:  -8.24074074074074153E-002,
776:  -1.38888888888888812E-002,
777:  -5.27777777777777429E-002,
778:  -1.38888888888888742E-002,
779:  -8.24074074074073876E-002,
780:  -1.38888888888889089E-002,
781:   1.38888888888888933E-002,
782:  -5.64814814814814589E-002,
783:   1.38888888888888812E-002,
784:   5.27777777777777429E-002,
785:  -8.24074074074073737E-002,
786:  -4.72222222222222071E-002,
787:   4.72222222222222002E-002,
788:  -7.87037037037036646E-002,
789:   1.38888888888889055E-002,
790:  -4.72222222222221932E-002,
791:   4.90740740740740339E-002,
792:   5.27777777777777499E-002,
793:  -5.27777777777777499E-002,
794:   0.18981481481481474     ,
795:   4.72222222222222002E-002,
796:  -1.38888888888888985E-002,
797:   4.90740740740740339E-002,
798:  -1.38888888888888846E-002,
799:   1.38888888888888812E-002,
800:   1.20370370370370284E-002,
801:  -7.87037037037036785E-002,
802:  -4.72222222222221932E-002,
803:  -4.72222222222222002E-002,
804:   1.20370370370370162E-002,
805:  -1.38888888888888812E-002,
806:  -1.38888888888888812E-002,
807:   4.90740740740740408E-002,
808:   4.72222222222222002E-002,
809:   1.38888888888888933E-002,
810:  -8.24074074074073876E-002,
811:   1.38888888888888708E-002,
812:  -5.27777777777777707E-002,
813:  -8.24074074074074153E-002,
814:  -5.27777777777777568E-002,
815:   1.38888888888888829E-002,
816:   4.90740740740740339E-002,
817:   1.38888888888889072E-002,
818:   4.72222222222222002E-002,
819:   0.18981481481481477     ,
820:   5.27777777777777429E-002,
821:   5.27777777777777568E-002,
822:  -5.64814814814814659E-002,
823:  -1.38888888888888846E-002,
824:  -1.38888888888888881E-002,
825:  -4.72222222222221932E-002,
826:  -7.87037037037036785E-002,
827:  -4.72222222222221932E-002,
828:   1.38888888888888673E-002,
829:  -8.24074074074073876E-002,
830:  -5.27777777777777568E-002,
831:   4.72222222222222002E-002,
832:   4.90740740740740269E-002,
833:   1.38888888888889020E-002,
834:  -1.38888888888888742E-002,
835:   1.20370370370370128E-002,
836:  -1.38888888888888829E-002,
837:  -5.27777777777777429E-002,
838:  -8.24074074074074153E-002,
839:   1.38888888888888777E-002,
840:  -1.38888888888889055E-002,
841:  -5.64814814814814659E-002,
842:  -1.38888888888888985E-002,
843:   5.27777777777777429E-002,
844:   0.18981481481481469     ,
845:   5.27777777777777429E-002,
846:   1.38888888888888933E-002,
847:   4.90740740740740339E-002,
848:   4.72222222222222071E-002,
849:  -4.72222222222222071E-002,
850:  -4.72222222222222002E-002,
851:  -7.87037037037036646E-002,
852:   1.38888888888888742E-002,
853:  -5.27777777777777568E-002,
854:  -8.24074074074073737E-002,
855:  -1.38888888888888951E-002,
856:  -1.38888888888888951E-002,
857:  -5.64814814814814589E-002,
858:  -5.27777777777777568E-002,
859:   1.38888888888888760E-002,
860:  -8.24074074074073876E-002,
861:  -1.38888888888888760E-002,
862:  -1.38888888888888812E-002,
863:   1.20370370370370249E-002,
864:   4.72222222222221932E-002,
865:   1.38888888888889003E-002,
866:   4.90740740740740339E-002,
867:   5.27777777777777568E-002,
868:   5.27777777777777429E-002,
869:   0.18981481481481474     ,
870:   1.38888888888888933E-002,
871:   4.72222222222222071E-002,
872:   4.90740740740740339E-002,
873:   1.20370370370370180E-002,
874:   1.38888888888888742E-002,
875:   1.38888888888888794E-002,
876:  -7.87037037037036785E-002,
877:   4.72222222222222071E-002,
878:   4.72222222222222002E-002,
879:  -8.24074074074073876E-002,
880:  -1.38888888888888846E-002,
881:   5.27777777777777568E-002,
882:   4.90740740740740616E-002,
883:  -4.72222222222222002E-002,
884:  -1.38888888888888881E-002,
885:   4.90740740740740408E-002,
886:  -1.38888888888888846E-002,
887:  -4.72222222222222002E-002,
888:  -8.24074074074074153E-002,
889:   5.27777777777777429E-002,
890:  -1.38888888888888846E-002,
891:  -5.64814814814814659E-002,
892:   1.38888888888888933E-002,
893:   1.38888888888888933E-002,
894:   0.18981481481481477     ,
895:  -5.27777777777777568E-002,
896:  -5.27777777777777637E-002,
897:  -1.38888888888888742E-002,
898:  -8.24074074074073598E-002,
899:  -5.27777777777777568E-002,
900:   4.72222222222222002E-002,
901:  -7.87037037037036924E-002,
902:  -4.72222222222222002E-002,
903:   1.38888888888888812E-002,
904:   1.20370370370370267E-002,
905:  -1.38888888888888794E-002,
906:  -4.72222222222222002E-002,
907:   4.90740740740740478E-002,
908:   1.38888888888888881E-002,
909:   1.38888888888888968E-002,
910:  -5.64814814814814659E-002,
911:  -1.38888888888888933E-002,
912:   5.27777777777777499E-002,
913:  -8.24074074074074153E-002,
914:   1.38888888888888812E-002,
915:  -1.38888888888888846E-002,
916:   4.90740740740740339E-002,
917:   4.72222222222222071E-002,
918:  -5.27777777777777568E-002,
919:   0.18981481481481477     ,
920:   5.27777777777777637E-002,
921:  -1.38888888888888829E-002,
922:  -5.27777777777777568E-002,
923:  -8.24074074074073598E-002,
924:   4.72222222222222071E-002,
925:  -4.72222222222222140E-002,
926:  -7.87037037037036924E-002,
927:   5.27777777777777637E-002,
928:   1.38888888888888916E-002,
929:  -8.24074074074073876E-002,
930:   1.38888888888888846E-002,
931:  -1.38888888888888951E-002,
932:  -5.64814814814814589E-002,
933:  -4.72222222222222071E-002,
934:   1.38888888888888812E-002,
935:   4.90740740740740339E-002,
936:   1.38888888888888829E-002,
937:  -1.38888888888888812E-002,
938:   1.20370370370370284E-002,
939:  -1.38888888888888881E-002,
940:   4.72222222222222071E-002,
941:   4.90740740740740339E-002,
942:  -5.27777777777777637E-002,
943:   5.27777777777777637E-002,
944:   0.18981481481481477     ,
945:   };

948:   PetscArraycpy(dd,DD,576);
949:   return(0);
950: }


953: /*TEST

955:    test:
956:       nsize: 8
957:       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 -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 -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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short
958:       filter: grep -v variant

960:    test:
961:       suffix: 2
962:       nsize: 8
963:       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 -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 -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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg
964:       filter: grep -v variant

966:    test:
967:       suffix: latebs
968:       filter: grep -v variant
969:       nsize: 8
970:       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 -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 -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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short -ksp_view

972:    test:
973:       suffix: latebs-2
974:       filter: grep -v variant
975:       nsize: 8
976:       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 -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 -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_use_parallel_coarse_grid_solver -mg_coarse_pc_type jacobi -mg_coarse_ksp_type cg -ksp_monitor_short -ksp_view

978:    test:
979:       suffix: ml
980:       nsize: 8
981:       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
982:       requires: ml

984:    test:
985:       suffix: nns
986:       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_max_it 10

988:    test:
989:       suffix: nns_telescope
990:       nsize: 2
991:       args: -use_mat_nearnullspace -ksp_monitor_short -pc_type telescope -pc_telescope_reduction_factor 2 -telescope_pc_type gamg

993:    test:
994:       suffix: seqaijmkl
995:       nsize: 8
996:       requires: mkl_sparse
997:       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 -pc_gamg_square_graph 1 -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

999: TEST*/