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