Actual source code: ex56.c

petsc-master 2015-07-31
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 *);

 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:     PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
248:     MatView(Amat,viewer);
249:     PetscViewerDestroy(&viewer);
250:   }

252:   /* finish KSP/PC setup */
253:   KSPSetOperators(ksp, Amat, Amat);
254:   if (use_nearnullspace) {
255:     MatNullSpace matnull;
256:     Vec          vec_coords;
257:     PetscScalar  *c;

259:     VecCreate(MPI_COMM_WORLD,&vec_coords);
260:     VecSetBlockSize(vec_coords,3);
261:     VecSetSizes(vec_coords,m,PETSC_DECIDE);
262:     VecSetUp(vec_coords);
263:     VecGetArray(vec_coords,&c);
264:     for (i=0; i<m; i++) c[i] = coords[i]; /* Copy since Scalar type might be Complex */
265:     VecRestoreArray(vec_coords,&c);
266:     MatNullSpaceCreateRigidBody(vec_coords,&matnull);
267:     MatSetNearNullSpace(Amat,matnull);
268:     MatNullSpaceDestroy(&matnull);
269:     VecDestroy(&vec_coords);
270:   } else {
271:     PC             pc;
272:     KSPGetPC(ksp,&pc);
273:     PCSetCoordinates(pc, 3, m/3, coords);
274:   }

276:   MaybeLogStagePush(stage[0]);

278:   /* PC setup basically */
279:   KSPSetUp(ksp);

281:   MaybeLogStagePop();
282:   MaybeLogStagePush(stage[1]);

284:   /* test BCs */
285:   if (test_nonzero_cols) {
286:     VecZeroEntries(xx);
287:     if (mype==0) VecSetValue(xx,0,1.0,INSERT_VALUES);
288:     VecAssemblyBegin(xx);
289:     VecAssemblyEnd(xx);
290:     KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
291:   }

293:   /* 1st solve */
294:   KSPSolve(ksp, bb, xx);

296:   MaybeLogStagePop();

298:   /* 2nd solve */
299:   if (two_solves) {
300:     PetscReal emax, emin;
301:     PetscReal norm,norm2;
302:     Vec       res;

304:     MaybeLogStagePush(stage[2]);
305:     /* PC setup basically */
306:     MatScale(Amat, 100000.0);
307:     KSPSetOperators(ksp, Amat, Amat);
308:     KSPSetUp(ksp);

310:     MaybeLogStagePop();
311:     MaybeLogStagePush(stage[3]);
312:     KSPSolve(ksp, bb, xx);
313:     KSPComputeExtremeSingularValues(ksp, &emax, &emin);

315:     MaybeLogStagePop();
316:     MaybeLogStagePush(stage[4]);

318:     /* 3rd solve */
319:     MatScale(Amat, 100000.0);
320:     KSPSetOperators(ksp, Amat, Amat);
321:     KSPSetUp(ksp);

323:     MaybeLogStagePop();
324:     MaybeLogStagePush(stage[5]);

326:     KSPSolve(ksp, bb, xx);

328:     MaybeLogStagePop();


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

333:     VecDuplicate(xx, &res);
334:     MatMult(Amat, xx, res);
335:     VecAXPY(bb, -1.0, res);
336:     VecDestroy(&res);
337:     VecNorm(bb, NORM_2, &norm);
338:     PetscPrintf(PETSC_COMM_WORLD,"[%d]%s |b-Ax|/|b|=%e, |b|=%e, emax=%e\n",0,__FUNCT__,(double)(norm/norm2),(double)norm2,(double)emax);
339:     /*PetscViewerASCIIOpen(comm, "residual.m", &viewer);
340:      PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
341:      VecView(bb,viewer);
342:      PetscViewerDestroy(&viewer);*/


345:     /* PetscViewerASCIIOpen(comm, "rhs.m", &viewer); */
346:     /* PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB); */
347:     /*  */
348:     /* VecView(bb,viewer); */
349:     /* PetscViewerDestroy(&viewer); */

351:     /* PetscViewerASCIIOpen(comm, "solution.m", &viewer); */
352:     /* PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB); */
353:     /*  */
354:     /* VecView(xx, viewer); */
355:     /* PetscViewerDestroy(&viewer); */
356:   }

358:   /* Free work space */
359:   KSPDestroy(&ksp);
360:   VecDestroy(&xx);
361:   VecDestroy(&bb);
362:   MatDestroy(&Amat);
363:   PetscFree(coords);

365:   PetscFinalize();
366:   return 0;
367: }

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

955:   PetscMemcpy(dd,DD,sizeof(PetscScalar)*576);
956:   return(0);
957: }