Actual source code: ex56.c
petsc-3.10.5 2019-03-28
1: static char help[] = "3D, tri-quadratic hexahedra (Q1), displacement finite element formulation\n\
2: of linear elasticity. E=1.0, nu=1/3.\n\
3: Unit cube domain with Dirichlet boundary\n\n";
5: #include <petscdmplex.h>
6: #include <petscsnes.h>
7: #include <petscds.h>
8: #include <petscdmforest.h>
10: static PetscReal s_soft_alpha=1.e-3;
11: static PetscReal s_mu=0.4;
12: static PetscReal s_lambda=0.4;
14: static void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
15: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
16: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17: PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
18: {
19: f0[0] = 1; /* x direction pull */
20: f0[1] = -x[2]; /* add a twist around x-axis */
21: f0[2] = x[1];
22: }
24: static void f1_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
25: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
26: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
27: PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
28: {
29: const PetscInt Ncomp = dim;
30: PetscInt comp, d;
31: for (comp = 0; comp < Ncomp; ++comp) {
32: for (d = 0; d < dim; ++d) {
33: f1[comp*dim+d] = 0.0;
34: }
35: }
36: }
38: /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
39: static void f1_u_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux,
40: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
41: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
42: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
43: {
44: PetscReal trace,mu=s_mu,lambda=s_lambda,rad;
45: PetscInt i,j;
46: for (i=0,rad=0.;i<dim;i++) {
47: PetscReal t=x[i];
48: rad += t*t;
49: }
50: rad = PetscSqrtReal(rad);
51: if (rad>0.25) {
52: mu *= s_soft_alpha;
53: lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
54: }
55: for (i=0,trace=0; i < dim; ++i) {
56: trace += PetscRealPart(u_x[i*dim+i]);
57: }
58: for (i=0; i < dim; ++i) {
59: for (j=0; j < dim; ++j) {
60: f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]);
61: }
62: f1[i*dim+i] += lambda*trace;
63: }
64: }
66: /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
67: static void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
68: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
69: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
70: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
71: {
72: PetscReal trace,mu=s_mu,lambda=s_lambda;
73: PetscInt i,j;
74: for (i=0,trace=0; i < dim; ++i) {
75: trace += PetscRealPart(u_x[i*dim+i]);
76: }
77: for (i=0; i < dim; ++i) {
78: for (j=0; j < dim; ++j) {
79: f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]);
80: }
81: f1[i*dim+i] += lambda*trace;
82: }
83: }
85: /* 3D elasticity */
86: #define IDX(ii,jj,kk,ll) (27*ii+9*jj+3*kk+ll)
88: void g3_uu_3d_private( PetscScalar g3[], const PetscReal mu, const PetscReal lambda)
89: {
90: if (1) {
91: g3[0] += lambda;
92: g3[0] += mu;
93: g3[0] += mu;
94: g3[4] += lambda;
95: g3[8] += lambda;
96: g3[10] += mu;
97: g3[12] += mu;
98: g3[20] += mu;
99: g3[24] += mu;
100: g3[28] += mu;
101: g3[30] += mu;
102: g3[36] += lambda;
103: g3[40] += lambda;
104: g3[40] += mu;
105: g3[40] += mu;
106: g3[44] += lambda;
107: g3[50] += mu;
108: g3[52] += mu;
109: g3[56] += mu;
110: g3[60] += mu;
111: g3[68] += mu;
112: g3[70] += mu;
113: g3[72] += lambda;
114: g3[76] += lambda;
115: g3[80] += lambda;
116: g3[80] += mu;
117: g3[80] += mu;
118: } else {
119: int i,j,k,l;
120: static int cc=-1;
121: cc++;
122: for (i = 0; i < 3; ++i) {
123: for (j = 0; j < 3; ++j) {
124: for (k = 0; k < 3; ++k) {
125: for (l = 0; l < 3; ++l) {
126: if (k==l && i==j) g3[IDX(i,j,k,l)] += lambda;
127: if (i==k && j==l) g3[IDX(i,j,k,l)] += mu;
128: if (i==l && j==k) g3[IDX(i,j,k,l)] += mu;
129: if (k==l && i==j && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += lambda;\n",IDX(i,j,k,l));
130: if (i==k && j==l && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l));
131: if (i==l && j==k && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l));
132: }
133: }
134: }
135: }
136: }
137: }
139: static void g3_uu_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux,
140: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
141: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
142: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
143: {
144: PetscReal mu=s_mu, lambda=s_lambda,rad;
145: PetscInt i;
146: for (i=0,rad=0.;i<dim;i++) {
147: PetscReal t=x[i];
148: rad += t*t;
149: }
150: rad = PetscSqrtReal(rad);
151: if (rad>0.25) {
152: mu *= s_soft_alpha;
153: lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
154: }
155: g3_uu_3d_private(g3,mu,lambda);
156: }
158: static void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
159: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
160: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
161: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
162: {
163: g3_uu_3d_private(g3,s_mu,s_lambda);
164: }
166: static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
167: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
168: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
169: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
170: {
171: const PetscInt Ncomp = dim;
172: PetscInt comp;
174: for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0;
175: }
177: /* PI_i (x_i^4 - x_i^2) */
178: static void f0_u_x4(PetscInt dim, PetscInt Nf, PetscInt NfAux,
179: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
180: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
181: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
182: {
183: const PetscInt Ncomp = dim;
184: PetscInt comp,i;
186: for (comp = 0; comp < Ncomp; ++comp) {
187: f0[comp] = 1e5;
188: for (i = 0; i < Ncomp; ++i) {
189: f0[comp] *= /* (comp+1)* */(x[i]*x[i]*x[i]*x[i] - x[i]*x[i]); /* assumes (0,1]^D domain */
190: }
191: }
192: }
194: PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
195: {
196: const PetscInt Ncomp = dim;
197: PetscInt comp;
199: for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0;
200: return 0;
201: }
203: int main(int argc,char **args)
204: {
205: Mat Amat;
207: SNES snes;
208: KSP ksp;
209: MPI_Comm comm;
210: PetscMPIInt rank;
211: PetscLogStage stage[7];
212: PetscBool test_nonzero_cols=PETSC_FALSE,use_nearnullspace=PETSC_TRUE;
213: Vec xx,bb;
214: PetscInt iter,i,N,dim=3,cells[3]={1,1,1},max_conv_its,local_sizes[7],run_type=1;
215: DM dm,distdm,basedm;
216: PetscBool flg;
217: char convType[256];
218: PetscReal Lx,mdisp[10],err[10];
219: const char * const options[10] = {"-ex56_dm_refine 0",
220: "-ex56_dm_refine 1",
221: "-ex56_dm_refine 2",
222: "-ex56_dm_refine 3",
223: "-ex56_dm_refine 4",
224: "-ex56_dm_refine 5",
225: "-ex56_dm_refine 6",
226: "-ex56_dm_refine 7",
227: "-ex56_dm_refine 8",
228: "-ex56_dm_refine 9"};
230: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
231: comm = PETSC_COMM_WORLD;
232: MPI_Comm_rank(comm, &rank);
233: /* options */
234: PetscOptionsBegin(comm,NULL,"3D bilinear Q1 elasticity options","");
235: {
236: i = 3;
237: PetscOptionsIntArray("-cells", "Number of (flux tube) processor in each dimension", "ex56.c", cells, &i, NULL);
239: Lx = 1.; /* or ne for rod */
240: max_conv_its = 3;
241: PetscOptionsInt("-max_conv_its","Number of iterations in convergence study","",max_conv_its,&max_conv_its,NULL);
242: if (max_conv_its<=0 || max_conv_its>7) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_USER, "Bad number of iterations for convergence test (%D)",max_conv_its);
243: PetscOptionsReal("-lx","Length of domain","",Lx,&Lx,NULL);
244: PetscOptionsReal("-alpha","material coefficient inside circle","",s_soft_alpha,&s_soft_alpha,NULL);
245: PetscOptionsBool("-test_nonzero_cols","nonzero test","",test_nonzero_cols,&test_nonzero_cols,NULL);
246: PetscOptionsBool("-use_mat_nearnullspace","MatNearNullSpace API test","",use_nearnullspace,&use_nearnullspace,NULL);
247: PetscOptionsInt("-run_type","0: twisting load on cantalever, 1: 3rd order accurate convergence test","",run_type,&run_type,NULL);
248: i = 3;
249: PetscOptionsInt("-mat_block_size","","",i,&i,&flg);
250: if (!flg || i!=3) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_USER, "'-mat_block_size 3' must be set (%D) and = 3 (%D)",flg,flg? i : 3);
251: }
252: PetscOptionsEnd();
253: PetscLogStageRegister("Mesh Setup", &stage[6]);
254: PetscLogStageRegister("1st Setup", &stage[0]);
255: PetscLogStageRegister("1st Solve", &stage[1]);
257: /* create DM, Plex calls DMSetup */
258: PetscLogStagePush(stage[6]);
259: DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm);
260: {
261: DMLabel label;
262: IS is;
263: DMCreateLabel(dm, "boundary");
264: DMGetLabel(dm, "boundary", &label);
265: DMPlexMarkBoundaryFaces(dm, 1, label);
266: if (run_type==0) {
267: DMGetStratumIS(dm, "boundary", 1, &is);
268: DMCreateLabel(dm,"Faces");
269: if (is) {
270: PetscInt d, f, Nf;
271: const PetscInt *faces;
272: PetscInt csize;
273: PetscSection cs;
274: Vec coordinates ;
275: DM cdm;
276: ISGetLocalSize(is, &Nf);
277: ISGetIndices(is, &faces);
278: DMGetCoordinatesLocal(dm, &coordinates);
279: DMGetCoordinateDM(dm, &cdm);
280: DMGetSection(cdm, &cs);
281: /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
282: for (f = 0; f < Nf; ++f) {
283: PetscReal faceCoord;
284: PetscInt b,v;
285: PetscScalar *coords = NULL;
286: PetscInt Nv;
287: DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords);
288: Nv = csize/dim; /* Calculate mean coordinate vector */
289: for (d = 0; d < dim; ++d) {
290: faceCoord = 0.0;
291: for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]);
292: faceCoord /= Nv;
293: for (b = 0; b < 2; ++b) {
294: if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */
295: DMSetLabelValue(dm, "Faces", faces[f], d*2+b+1);
296: }
297: }
298: }
299: DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords);
300: }
301: ISRestoreIndices(is, &faces);
302: }
303: ISDestroy(&is);
304: DMGetLabel(dm, "Faces", &label);
305: DMPlexLabelComplete(dm, label);
306: }
307: }
308: {
309: PetscInt dimEmbed, i;
310: PetscInt nCoords;
311: PetscScalar *coords,bounds[] = {0,1,-.5,.5,-.5,.5,}; /* x_min,x_max,y_min,y_max */
312: Vec coordinates;
313: bounds[1] = Lx;
314: if (run_type==1) {
315: for (i = 0; i < 2*dim; i++) bounds[i] = (i%2) ? 1 : 0;
316: }
317: DMGetCoordinatesLocal(dm,&coordinates);
318: DMGetCoordinateDim(dm,&dimEmbed);
319: if (dimEmbed != dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"dimEmbed != dim %D",dimEmbed);
320: VecGetLocalSize(coordinates,&nCoords);
321: if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
322: VecGetArray(coordinates,&coords);
323: for (i = 0; i < nCoords; i += dimEmbed) {
324: PetscInt j;
325: PetscScalar *coord = &coords[i];
326: for (j = 0; j < dimEmbed; j++) {
327: coord[j] = bounds[2 * j] + coord[j] * (bounds[2 * j + 1] - bounds[2 * j]);
328: }
329: }
330: VecRestoreArray(coordinates,&coords);
331: DMSetCoordinatesLocal(dm,coordinates);
332: }
334: /* convert to p4est, and distribute */
336: PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
337: PetscOptionsFList("-dm_type","Convert DMPlex to another format (should not be Plex!)","ex56.c",DMList,DMPLEX,convType,256,&flg);
338: PetscOptionsEnd();
339: if (flg) {
340: DM newdm;
341: DMConvert(dm,convType,&newdm);
342: if (newdm) {
343: const char *prefix;
344: PetscBool isForest;
345: PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix);
346: PetscObjectSetOptionsPrefix((PetscObject)newdm,prefix);
347: DMIsForest(newdm,&isForest);
348: if (isForest) {
349: } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Converted to non Forest?");
350: DMDestroy(&dm);
351: dm = newdm;
352: } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Convert failed?");
353: } else {
354: PetscPartitioner part;
355: /* Plex Distribute mesh over processes */
356: #if 1
357: DMPlexGetPartitioner(dm,&part);
358: PetscPartitionerSetFromOptions(part);
359: #endif
360: DMPlexDistribute(dm, 0, NULL, &distdm);
361: if (distdm) {
362: const char *prefix;
363: PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix);
364: PetscObjectSetOptionsPrefix((PetscObject)distdm,prefix);
365: DMDestroy(&dm);
366: dm = distdm;
367: }
368: }
369: PetscLogStagePop();
370: basedm = dm; dm = NULL;
372: for (iter=0 ; iter<max_conv_its ; iter++) {
373: PetscLogStagePush(stage[6]);
374: /* make new DM */
375: DMClone(basedm, &dm);
376: PetscObjectSetOptionsPrefix((PetscObject) dm, "ex56_");
377: PetscObjectSetName( (PetscObject)dm,"Mesh");
378: PetscOptionsClearValue(NULL,"-ex56_dm_refine");
379: PetscOptionsInsertString(NULL,options[iter]);
380: DMSetFromOptions(dm); /* refinement done here in Plex, p4est */
381: /* snes */
382: SNESCreate(comm, &snes);
383: SNESSetDM(snes, dm);
384: /* fem */
385: {
386: const PetscInt Ncomp = dim;
387: const PetscInt components[] = {0,1,2};
388: const PetscInt Nfid = 1, Npid = 1;
389: const PetscInt fid[] = {1}; /* The fixed faces (x=0) */
390: const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */
391: PetscFE fe;
392: PetscDS prob;
393: DM cdm = dm;
395: PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, PETSC_FALSE, NULL, PETSC_DECIDE, &fe); /* elasticity */
396: PetscObjectSetName((PetscObject) fe, "deformation");
397: /* FEM prob */
398: DMGetDS(dm, &prob);
399: PetscDSSetDiscretization(prob, 0, (PetscObject) fe);
400: /* setup problem */
401: if (run_type==1) {
402: PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d);
403: PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d);
404: } else {
405: PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha);
406: PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha);
407: PetscDSSetBdResidual(prob, 0, f0_bd_u_3d, f1_bd_u);
408: }
409: /* bcs */
410: if (run_type==1) {
411: PetscInt id = 1;
412: DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "boundary", 0, 0, NULL, (void (*)(void)) zero, 1, &id, NULL);
413: } else {
414: PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "fixed", "Faces", 0, Ncomp, components, (void (*)(void)) zero, Nfid, fid, NULL);
415: PetscDSAddBoundary(prob, DM_BC_NATURAL, "traction", "Faces", 0, Ncomp, components, NULL, Npid, pid, NULL);
416: }
417: while (cdm) {
418: DMSetDS(cdm,prob);
419: DMGetCoarseDM(cdm, &cdm);
420: }
421: PetscFEDestroy(&fe);
422: }
423: /* vecs & mat */
424: DMCreateGlobalVector(dm,&xx);
425: VecDuplicate(xx, &bb);
426: PetscObjectSetName((PetscObject) bb, "b");
427: PetscObjectSetName((PetscObject) xx, "u");
428: DMCreateMatrix(dm, &Amat);
429: VecGetSize(bb,&N);
430: local_sizes[iter] = N;
431: PetscPrintf(PETSC_COMM_WORLD,"[%d] %D global equations, %D vertices\n",rank,N,N/dim);
432: if (use_nearnullspace && N/dim > 1) {
433: /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
434: DM subdm;
435: MatNullSpace nearNullSpace;
436: PetscInt fields = 0;
437: PetscObject deformation;
438: DMCreateSubDM(dm, 1, &fields, NULL, &subdm);
439: DMPlexCreateRigidBody(subdm, &nearNullSpace);
440: DMGetField(dm, 0, &deformation);
441: PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace);
442: DMDestroy(&subdm);
443: MatNullSpaceDestroy(&nearNullSpace); /* created by DM and destroyed by Mat */
444: }
445: DMPlexSetSNESLocalFEM(dm,NULL,NULL,NULL);
446: SNESSetJacobian(snes, Amat, Amat, NULL, NULL);
447: SNESSetFromOptions(snes);
448: DMSetUp(dm);
449: PetscLogStagePop();
450: PetscLogStagePush(stage[0]);
451: /* ksp */
452: SNESGetKSP(snes, &ksp);
453: KSPSetComputeSingularValues(ksp,PETSC_TRUE);
454: /* test BCs */
455: VecZeroEntries(xx);
456: if (test_nonzero_cols) {
457: if (rank==0) VecSetValue(xx,0,1.0,INSERT_VALUES);
458: VecAssemblyBegin(xx);
459: VecAssemblyEnd(xx);
460: }
461: VecZeroEntries(bb);
462: VecGetSize(bb,&i);
463: local_sizes[iter] = i;
464: PetscPrintf(PETSC_COMM_WORLD,"[%d] %D equations in vector, %D vertices\n",rank,i,i/dim);
465: /* setup solver, dummy solve to really setup */
466: if (0) {
467: KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
468: SNESSolve(snes, bb, xx);
469: KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,50);
470: VecZeroEntries(xx);
471: }
472: PetscLogStagePop();
473: /* solve */
474: PetscLogStagePush(stage[1]);
475: SNESSolve(snes, bb, xx);
476: PetscLogStagePop();
477: VecNorm(xx,NORM_INFINITY,&mdisp[iter]);
478: DMViewFromOptions(dm, NULL, "-dm_view");
479: {
480: PetscViewer viewer = NULL;
481: PetscViewerFormat fmt;
482: PetscOptionsGetViewer(comm,"ex56_","-vec_view",&viewer,&fmt,&flg);
483: if (flg) {
484: PetscViewerPushFormat(viewer,fmt);
485: VecView(xx,viewer);
486: VecView(bb,viewer);
487: PetscViewerPopFormat(viewer);
488: }
489: PetscViewerDestroy(&viewer);
490: }
491: /* Free work space */
492: DMDestroy(&dm);
493: SNESDestroy(&snes);
494: VecDestroy(&xx);
495: VecDestroy(&bb);
496: MatDestroy(&Amat);
497: }
498: DMDestroy(&basedm);
499: if (run_type==1) {
500: err[0] = 59.975208 - mdisp[0]; /* error with what I think is the exact solution */
501: } else {
502: err[0] = 171.038 - mdisp[0];
503: }
504: for (iter=1 ; iter<max_conv_its ; iter++) {
505: if (run_type==1) {
506: err[iter] = 59.975208 - mdisp[iter];
507: } else {
508: err[iter] = 171.038 - mdisp[iter];
509: }
510: PetscPrintf(PETSC_COMM_WORLD,"[%d] %D) N=%12D, max displ=%9.7e, disp diff=%9.2e, error=%4.3e, rate=%3.2g\n",rank,iter,local_sizes[iter],mdisp[iter],
511: mdisp[iter]-mdisp[iter-1],err[iter],PetscLogReal(err[iter-1]/err[iter])/PetscLogReal(2.));
512: }
514: PetscFinalize();
515: return ierr;
516: }
518: /*TEST
520: test:
521: suffix: 0
522: nsize: 4
523: requires: !single
524: args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -matptap_via scalable -ex56_dm_view -run_type 1
525: timeoutfactor: 2
527: test:
528: suffix: hypre
529: requires: hypre !single
530: nsize: 4
531: args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -pc_type hypre -pc_hypre_type boomeramg -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -ksp_converged_reason -use_mat_nearnullspace true -mat_block_size 3 -petscpartitioner_type simple
533: test:
534: suffix: ml
535: requires: ml !single
536: nsize: 4
537: args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type ml -mg_levels_ksp_type chebyshev -mg_levels_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type sor -mg_levels_esteig_ksp_type cg -mat_block_size 3 -petscpartitioner_type simple -use_mat_nearnullspace
539: test:
540: nsize: 4
541: requires: parmetis !single
542: args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -pc_gamg_mat_partitioning_type parmetis -mat_block_size 3 -run_type 1 -pc_gamg_repartition true
544: TEST*/