Actual source code: ex49.c

  1: static char help[] =  "   Solves the compressible plane strain elasticity equations in 2d on the unit domain using Q1 finite elements. \n\
  2:    Material properties E (Youngs modulus) and nu (Poisson ratio) may vary as a function of space. \n\
  3:    The model utilises boundary conditions which produce compression in the x direction. \n\
  4: Options: \n"
  5: "\
  6:      -mx : number of elements in x-direction \n\
  7:      -my : number of elements in y-direction \n\
  8:      -c_str : structure of the coefficients to use. \n"
  9: "\
 10:           -c_str 0 => isotropic material with constant coefficients. \n\
 11:                          Parameters: \n\
 12:                              -iso_E  : Youngs modulus \n\
 13:                              -iso_nu : Poisson ratio \n\
 14:           -c_str 1 => step function in the material properties in x. \n\
 15:                          Parameters: \n\
 16:                               -step_E0  : Youngs modulus to the left of the step \n\
 17:                               -step_nu0 : Poisson ratio to the left of the step \n\
 18:                               -step_E1  : Youngs modulus to the right of the step \n\
 19:                               -step_n1  : Poisson ratio to the right of the step \n\
 20:                               -step_xc  : x coordinate of the step \n"
 21: "\
 22:           -c_str 2 => checkerboard material with alternating properties. \n\
 23:                       Repeats the following pattern throughout the domain. For example with 4 materials specified, we would heve \n\
 24:                       -------------------------\n\
 25:                       |  D  |  A  |  B  |  C  |\n\
 26:                       ------|-----|-----|------\n\
 27:                       |  C  |  D  |  A  |  B  |\n\
 28:                       ------|-----|-----|------\n\
 29:                       |  B  |  C  |  D  |  A  |\n\
 30:                       ------|-----|-----|------\n\
 31:                       |  A  |  B  |  C  |  D  |\n\
 32:                       -------------------------\n\
 33:                       \n\
 34:                          Parameters: \n\
 35:                               -brick_E    : a comma separated list of Young's modulii \n\
 36:                               -brick_nu   : a comma separated list of Poisson ratios  \n\
 37:                               -brick_span : the number of elements in x and y each brick will span \n\
 38:           -c_str 3 => sponge-like material with alternating properties. \n\
 39:                       Repeats the following pattern throughout the domain \n"
 40: "\
 41:                       -----------------------------\n\
 42:                       |       [background]        |\n\
 43:                       |          E0,nu0           |\n\
 44:                       |     -----------------     |\n\
 45:                       |     |  [inclusion]  |     |\n\
 46:                       |     |    E1,nu1     |     |\n\
 47:                       |     |               |     |\n\
 48:                       |     | <---- w ----> |     |\n\
 49:                       |     |               |     |\n\
 50:                       |     |               |     |\n\
 51:                       |     -----------------     |\n\
 52:                       |                           |\n\
 53:                       |                           |\n\
 54:                       -----------------------------\n\
 55:                       <--------  t + w + t ------->\n\
 56:                       \n\
 57:                          Parameters: \n\
 58:                               -sponge_E0  : Youngs modulus of the surrounding material \n\
 59:                               -sponge_E1  : Youngs modulus of the inclusion \n\
 60:                               -sponge_nu0 : Poisson ratio of the surrounding material \n\
 61:                               -sponge_nu1 : Poisson ratio of the inclusion \n\
 62:                               -sponge_t   : the number of elements defining the border around each inclusion \n\
 63:                               -sponge_w   : the number of elements in x and y each inclusion will span\n\
 64:      -use_gp_coords : Evaluate the Youngs modulus, Poisson ratio and the body force at the global coordinates of the quadrature points.\n\
 65:      By default, E, nu and the body force are evaulated at the element center and applied as a constant over the entire element.\n\
 66:      -use_nonsymbc : Option to use non-symmetric boundary condition imposition. This choice will use less memory.";

 68: /* Contributed by Dave May */

 70: #include <petscksp.h>
 71: #include <petscdm.h>
 72: #include <petscdmda.h>

 74: static PetscErrorCode DMDABCApplyCompression(DM,Mat,Vec);
 75: static PetscErrorCode DMDABCApplySymmetricCompression(DM elas_da,Mat A,Vec f,IS *dofs,Mat *AA,Vec *ff);

 77: #define NSD            2 /* number of spatial dimensions */
 78: #define NODES_PER_EL   4 /* nodes per element */
 79: #define U_DOFS         2 /* degrees of freedom per displacement node */
 80: #define GAUSS_POINTS   4

 82: /* cell based evaluation */
 83: typedef struct {
 84:   PetscScalar E,nu,fx,fy;
 85: } Coefficients;

 87: /* Gauss point based evaluation 8+4+4+4 = 20 */
 88: typedef struct {
 89:   PetscScalar gp_coords[2*GAUSS_POINTS];
 90:   PetscScalar E[GAUSS_POINTS];
 91:   PetscScalar nu[GAUSS_POINTS];
 92:   PetscScalar fx[GAUSS_POINTS];
 93:   PetscScalar fy[GAUSS_POINTS];
 94: } GaussPointCoefficients;

 96: typedef struct {
 97:   PetscScalar ux_dof;
 98:   PetscScalar uy_dof;
 99: } ElasticityDOF;

101: /*

103:  D = E/((1+nu)(1-2nu)) * [ 1-nu   nu        0     ]
104:                          [  nu   1-nu       0     ]
105:                          [  0     0   0.5*(1-2nu) ]

107:  B = [ d_dx   0   ]
108:      [  0    d_dy ]
109:      [ d_dy  d_dx ]

111:  */

113: /* FEM routines */
114: /*
115:  Element: Local basis function ordering
116:  1-----2
117:  |     |
118:  |     |
119:  0-----3
120:  */
121: static void ConstructQ12D_Ni(PetscScalar _xi[],PetscScalar Ni[])
122: {
123:   PetscScalar xi  = _xi[0];
124:   PetscScalar eta = _xi[1];

126:   Ni[0] = 0.25*(1.0-xi)*(1.0-eta);
127:   Ni[1] = 0.25*(1.0-xi)*(1.0+eta);
128:   Ni[2] = 0.25*(1.0+xi)*(1.0+eta);
129:   Ni[3] = 0.25*(1.0+xi)*(1.0-eta);
130: }

132: static void ConstructQ12D_GNi(PetscScalar _xi[],PetscScalar GNi[][NODES_PER_EL])
133: {
134:   PetscScalar xi  = _xi[0];
135:   PetscScalar eta = _xi[1];

137:   GNi[0][0] = -0.25*(1.0-eta);
138:   GNi[0][1] = -0.25*(1.0+eta);
139:   GNi[0][2] =   0.25*(1.0+eta);
140:   GNi[0][3] =   0.25*(1.0-eta);

142:   GNi[1][0] = -0.25*(1.0-xi);
143:   GNi[1][1] =   0.25*(1.0-xi);
144:   GNi[1][2] =   0.25*(1.0+xi);
145:   GNi[1][3] = -0.25*(1.0+xi);
146: }

148: static void ConstructQ12D_GNx(PetscScalar GNi[][NODES_PER_EL],PetscScalar GNx[][NODES_PER_EL],PetscScalar coords[],PetscScalar *det_J)
149: {
150:   PetscScalar J00,J01,J10,J11,J;
151:   PetscScalar iJ00,iJ01,iJ10,iJ11;
152:   PetscInt    i;

154:   J00 = J01 = J10 = J11 = 0.0;
155:   for (i = 0; i < NODES_PER_EL; i++) {
156:     PetscScalar cx = coords[2*i+0];
157:     PetscScalar cy = coords[2*i+1];

159:     J00 = J00+GNi[0][i]*cx;      /* J_xx = dx/dxi */
160:     J01 = J01+GNi[0][i]*cy;      /* J_xy = dy/dxi */
161:     J10 = J10+GNi[1][i]*cx;      /* J_yx = dx/deta */
162:     J11 = J11+GNi[1][i]*cy;      /* J_yy = dy/deta */
163:   }
164:   J = (J00*J11)-(J01*J10);

166:   iJ00 =  J11/J;
167:   iJ01 = -J01/J;
168:   iJ10 = -J10/J;
169:   iJ11 =  J00/J;

171:   for (i = 0; i < NODES_PER_EL; i++) {
172:     GNx[0][i] = GNi[0][i]*iJ00+GNi[1][i]*iJ01;
173:     GNx[1][i] = GNi[0][i]*iJ10+GNi[1][i]*iJ11;
174:   }

176:   if (det_J) *det_J = J;
177: }

179: static void ConstructGaussQuadrature(PetscInt *ngp,PetscScalar gp_xi[][2],PetscScalar gp_weight[])
180: {
181:   *ngp         = 4;
182:   gp_xi[0][0]  = -0.57735026919;gp_xi[0][1] = -0.57735026919;
183:   gp_xi[1][0]  = -0.57735026919;gp_xi[1][1] =  0.57735026919;
184:   gp_xi[2][0]  =  0.57735026919;gp_xi[2][1] =  0.57735026919;
185:   gp_xi[3][0]  =  0.57735026919;gp_xi[3][1] = -0.57735026919;
186:   gp_weight[0] = 1.0;
187:   gp_weight[1] = 1.0;
188:   gp_weight[2] = 1.0;
189:   gp_weight[3] = 1.0;
190: }

192: static PetscErrorCode DMDAGetElementOwnershipRanges2d(DM da,PetscInt **_lx,PetscInt **_ly)
193: {
195:   PetscMPIInt    rank;
196:   PetscInt       proc_I,proc_J;
197:   PetscInt       cpu_x,cpu_y;
198:   PetscInt       local_mx,local_my;
199:   Vec            vlx,vly;
200:   PetscInt       *LX,*LY,i;
201:   PetscScalar    *_a;
202:   Vec            V_SEQ;
203:   VecScatter     ctx;

206:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

208:   DMDAGetInfo(da,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);

210:   proc_J = rank/cpu_x;
211:   proc_I = rank-cpu_x*proc_J;

213:   PetscMalloc1(cpu_x,&LX);
214:   PetscMalloc1(cpu_y,&LY);

216:   DMDAGetElementsSizes(da,&local_mx,&local_my,NULL);
217:   VecCreate(PETSC_COMM_WORLD,&vlx);
218:   VecSetSizes(vlx,PETSC_DECIDE,cpu_x);
219:   VecSetFromOptions(vlx);

221:   VecCreate(PETSC_COMM_WORLD,&vly);
222:   VecSetSizes(vly,PETSC_DECIDE,cpu_y);
223:   VecSetFromOptions(vly);

225:   VecSetValue(vlx,proc_I,(PetscScalar)(local_mx+1.0e-9),INSERT_VALUES);
226:   VecSetValue(vly,proc_J,(PetscScalar)(local_my+1.0e-9),INSERT_VALUES);
227:   VecAssemblyBegin(vlx);VecAssemblyEnd(vlx);
228:   VecAssemblyBegin(vly);VecAssemblyEnd(vly);

230:   VecScatterCreateToAll(vlx,&ctx,&V_SEQ);
231:   VecScatterBegin(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
232:   VecScatterEnd(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
233:   VecGetArray(V_SEQ,&_a);
234:   for (i = 0; i < cpu_x; i++) LX[i] = (PetscInt)PetscRealPart(_a[i]);
235:   VecRestoreArray(V_SEQ,&_a);
236:   VecScatterDestroy(&ctx);
237:   VecDestroy(&V_SEQ);

239:   VecScatterCreateToAll(vly,&ctx,&V_SEQ);
240:   VecScatterBegin(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
241:   VecScatterEnd(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
242:   VecGetArray(V_SEQ,&_a);
243:   for (i = 0; i < cpu_y; i++) LY[i] = (PetscInt)PetscRealPart(_a[i]);
244:   VecRestoreArray(V_SEQ,&_a);
245:   VecScatterDestroy(&ctx);
246:   VecDestroy(&V_SEQ);

248:   *_lx = LX;
249:   *_ly = LY;

251:   VecDestroy(&vlx);
252:   VecDestroy(&vly);
253:   return(0);
254: }

256: static PetscErrorCode DMDACoordViewGnuplot2d(DM da,const char prefix[])
257: {
258:   DM             cda;
259:   Vec            coords;
260:   DMDACoor2d     **_coords;
261:   PetscInt       si,sj,nx,ny,i,j;
262:   FILE           *fp;
263:   char           fname[PETSC_MAX_PATH_LEN];
264:   PetscMPIInt    rank;

268:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
269:   PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
270:   PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
271:   if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");
272:   PetscFPrintf(PETSC_COMM_SELF,fp,"### Element geometry for processor %1.4d ### \n",rank);

274:   DMGetCoordinateDM(da,&cda);
275:   DMGetCoordinatesLocal(da,&coords);
276:   DMDAVecGetArray(cda,coords,&_coords);
277:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
278:   for (j = sj; j < sj+ny-1; j++) {
279:     for (i = si; i < si+nx-1; i++) {
280:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j][i].x),PetscRealPart(_coords[j][i].y));
281:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j+1][i].x),PetscRealPart(_coords[j+1][i].y));
282:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j+1][i+1].x),PetscRealPart(_coords[j+1][i+1].y));
283:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j][i+1].x),PetscRealPart(_coords[j][i+1].y));
284:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n\n",PetscRealPart(_coords[j][i].x),PetscRealPart(_coords[j][i].y));
285:     }
286:   }
287:   DMDAVecRestoreArray(cda,coords,&_coords);

289:   PetscFClose(PETSC_COMM_SELF,fp);
290:   return(0);
291: }

293: static PetscErrorCode DMDAViewGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
294: {
295:   DM             cda;
296:   Vec            coords,local_fields;
297:   DMDACoor2d     **_coords;
298:   FILE           *fp;
299:   char           fname[PETSC_MAX_PATH_LEN];
300:   const char     *field_name;
301:   PetscMPIInt    rank;
302:   PetscInt       si,sj,nx,ny,i,j;
303:   PetscInt       n_dofs,d;
304:   PetscScalar    *_fields;

308:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
309:   PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
310:   PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
311:   if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");

313:   PetscFPrintf(PETSC_COMM_SELF,fp,"### %s (processor %1.4d) ### \n",comment,rank);
314:   DMDAGetInfo(da,0,0,0,0,0,0,0,&n_dofs,0,0,0,0,0);
315:   PetscFPrintf(PETSC_COMM_SELF,fp,"### x y ");
316:   for (d = 0; d < n_dofs; d++) {
317:     DMDAGetFieldName(da,d,&field_name);
318:     PetscFPrintf(PETSC_COMM_SELF,fp,"%s ",field_name);
319:   }
320:   PetscFPrintf(PETSC_COMM_SELF,fp,"###\n");

322:   DMGetCoordinateDM(da,&cda);
323:   DMGetCoordinatesLocal(da,&coords);
324:   DMDAVecGetArray(cda,coords,&_coords);
325:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);

327:   DMCreateLocalVector(da,&local_fields);
328:   DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
329:   DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
330:   VecGetArray(local_fields,&_fields);

332:   for (j = sj; j < sj+ny; j++) {
333:     for (i = si; i < si+nx; i++) {
334:       PetscScalar coord_x,coord_y;
335:       PetscScalar field_d;

337:       coord_x = _coords[j][i].x;
338:       coord_y = _coords[j][i].y;

340:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",PetscRealPart(coord_x),PetscRealPart(coord_y));
341:       for (d = 0; d < n_dofs; d++) {
342:         field_d = _fields[n_dofs*((i-si)+(j-sj)*(nx))+d];
343:         PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e ",PetscRealPart(field_d));
344:       }
345:       PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
346:     }
347:   }
348:   VecRestoreArray(local_fields,&_fields);
349:   VecDestroy(&local_fields);

351:   DMDAVecRestoreArray(cda,coords,&_coords);

353:   PetscFClose(PETSC_COMM_SELF,fp);
354:   return(0);
355: }

357: static PetscErrorCode DMDAViewCoefficientsGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
358: {
359:   DM                     cda;
360:   Vec                    local_fields;
361:   FILE                   *fp;
362:   char                   fname[PETSC_MAX_PATH_LEN];
363:   const char             *field_name;
364:   PetscMPIInt            rank;
365:   PetscInt               si,sj,nx,ny,i,j,p;
366:   PetscInt               n_dofs,d;
367:   GaussPointCoefficients **_coefficients;
368:   PetscErrorCode         ierr;

371:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
372:   PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
373:   PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
374:   if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");

376:   PetscFPrintf(PETSC_COMM_SELF,fp,"### %s (processor %1.4d) ### \n",comment,rank);
377:   DMDAGetInfo(da,0,0,0,0,0,0,0,&n_dofs,0,0,0,0,0);
378:   PetscFPrintf(PETSC_COMM_SELF,fp,"### x y ");
379:   for (d = 0; d < n_dofs; d++) {
380:     DMDAGetFieldName(da,d,&field_name);
381:     PetscFPrintf(PETSC_COMM_SELF,fp,"%s ",field_name);
382:   }
383:   PetscFPrintf(PETSC_COMM_SELF,fp,"###\n");

385:   DMGetCoordinateDM(da,&cda);
386:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);

388:   DMCreateLocalVector(da,&local_fields);
389:   DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
390:   DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
391:   DMDAVecGetArray(da,local_fields,&_coefficients);

393:   for (j = sj; j < sj+ny; j++) {
394:     for (i = si; i < si+nx; i++) {
395:       PetscScalar coord_x,coord_y;

397:       for (p = 0; p < GAUSS_POINTS; p++) {
398:         coord_x = _coefficients[j][i].gp_coords[2*p];
399:         coord_y = _coefficients[j][i].gp_coords[2*p+1];

401:         PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",PetscRealPart(coord_x),PetscRealPart(coord_y));

403:         PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e %1.6e %1.6e",
404:                             PetscRealPart(_coefficients[j][i].E[p]),PetscRealPart(_coefficients[j][i].nu[p]),
405:                             PetscRealPart(_coefficients[j][i].fx[p]),PetscRealPart(_coefficients[j][i].fy[p]));
406:         PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
407:       }
408:     }
409:   }
410:   DMDAVecRestoreArray(da,local_fields,&_coefficients);
411:   VecDestroy(&local_fields);

413:   PetscFClose(PETSC_COMM_SELF,fp);
414:   return(0);
415: }

417: static void FormStressOperatorQ1(PetscScalar Ke[],PetscScalar coords[],PetscScalar E[],PetscScalar nu[])
418: {
419:   PetscInt    ngp;
420:   PetscScalar gp_xi[GAUSS_POINTS][2];
421:   PetscScalar gp_weight[GAUSS_POINTS];
422:   PetscInt    p,i,j,k,l;
423:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
424:   PetscScalar J_p;
425:   PetscScalar B[3][U_DOFS*NODES_PER_EL];
426:   PetscScalar prop_E,prop_nu,factor,constit_D[3][3];

428:   /* define quadrature rule */
429:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

431:   /* evaluate integral */
432:   for (p = 0; p < ngp; p++) {
433:     ConstructQ12D_GNi(gp_xi[p],GNi_p);
434:     ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);

436:     for (i = 0; i < NODES_PER_EL; i++) {
437:       PetscScalar d_dx_i = GNx_p[0][i];
438:       PetscScalar d_dy_i = GNx_p[1][i];

440:       B[0][2*i] = d_dx_i;  B[0][2*i+1] = 0.0;
441:       B[1][2*i] = 0.0;     B[1][2*i+1] = d_dy_i;
442:       B[2][2*i] = d_dy_i;  B[2][2*i+1] = d_dx_i;
443:     }

445:     /* form D for the quadrature point */
446:     prop_E          = E[p];
447:     prop_nu         = nu[p];
448:     factor          = prop_E / ((1.0+prop_nu)*(1.0-2.0*prop_nu));
449:     constit_D[0][0] = 1.0-prop_nu;  constit_D[0][1] = prop_nu;      constit_D[0][2] = 0.0;
450:     constit_D[1][0] = prop_nu;      constit_D[1][1] = 1.0-prop_nu;  constit_D[1][2] = 0.0;
451:     constit_D[2][0] = 0.0;          constit_D[2][1] = 0.0;          constit_D[2][2] = 0.5*(1.0-2.0*prop_nu);
452:     for (i = 0; i < 3; i++) {
453:       for (j = 0; j < 3; j++) {
454:         constit_D[i][j] = factor * constit_D[i][j] * gp_weight[p] * J_p;
455:       }
456:     }

458:     /* form Bt tildeD B */
459:     /*
460:      Ke_ij = Bt_ik . D_kl . B_lj
461:      = B_ki . D_kl . B_lj
462:      */
463:     for (i = 0; i < 8; i++) {
464:       for (j = 0; j < 8; j++) {
465:         for (k = 0; k < 3; k++) {
466:           for (l = 0; l < 3; l++) {
467:             Ke[8*i+j] = Ke[8*i+j] + B[k][i] * constit_D[k][l] * B[l][j];
468:           }
469:         }
470:       }
471:     }

473:   } /* end quadrature */
474: }

476: static void FormMomentumRhsQ1(PetscScalar Fe[],PetscScalar coords[],PetscScalar fx[],PetscScalar fy[])
477: {
478:   PetscInt    ngp;
479:   PetscScalar gp_xi[GAUSS_POINTS][2];
480:   PetscScalar gp_weight[GAUSS_POINTS];
481:   PetscInt    p,i;
482:   PetscScalar Ni_p[NODES_PER_EL];
483:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
484:   PetscScalar J_p,fac;

486:   /* define quadrature rule */
487:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

489:   /* evaluate integral */
490:   for (p = 0; p < ngp; p++) {
491:     ConstructQ12D_Ni(gp_xi[p],Ni_p);
492:     ConstructQ12D_GNi(gp_xi[p],GNi_p);
493:     ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
494:     fac = gp_weight[p]*J_p;

496:     for (i = 0; i < NODES_PER_EL; i++) {
497:       Fe[NSD*i]   += fac*Ni_p[i]*fx[p];
498:       Fe[NSD*i+1] += fac*Ni_p[i]*fy[p];
499:     }
500:   }
501: }

503: /*
504:  i,j are the element indices
505:  The unknown is a vector quantity.
506:  The s[].c is used to indicate the degree of freedom.
507:  */
508: static PetscErrorCode DMDAGetElementEqnums_u(MatStencil s_u[],PetscInt i,PetscInt j)
509: {
511:   /* displacement */
512:   /* node 0 */
513:   s_u[0].i = i;s_u[0].j = j;s_u[0].c = 0;          /* Ux0 */
514:   s_u[1].i = i;s_u[1].j = j;s_u[1].c = 1;          /* Uy0 */

516:   /* node 1 */
517:   s_u[2].i = i;s_u[2].j = j+1;s_u[2].c = 0;        /* Ux1 */
518:   s_u[3].i = i;s_u[3].j = j+1;s_u[3].c = 1;        /* Uy1 */

520:   /* node 2 */
521:   s_u[4].i = i+1;s_u[4].j = j+1;s_u[4].c = 0;      /* Ux2 */
522:   s_u[5].i = i+1;s_u[5].j = j+1;s_u[5].c = 1;      /* Uy2 */

524:   /* node 3 */
525:   s_u[6].i = i+1;s_u[6].j = j;s_u[6].c = 0;        /* Ux3 */
526:   s_u[7].i = i+1;s_u[7].j = j;s_u[7].c = 1;        /* Uy3 */
527:   return(0);
528: }

530: static PetscErrorCode GetElementCoords(DMDACoor2d **_coords,PetscInt ei,PetscInt ej,PetscScalar el_coords[])
531: {
533:   /* get coords for the element */
534:   el_coords[NSD*0+0] = _coords[ej][ei].x;      el_coords[NSD*0+1] = _coords[ej][ei].y;
535:   el_coords[NSD*1+0] = _coords[ej+1][ei].x;    el_coords[NSD*1+1] = _coords[ej+1][ei].y;
536:   el_coords[NSD*2+0] = _coords[ej+1][ei+1].x;  el_coords[NSD*2+1] = _coords[ej+1][ei+1].y;
537:   el_coords[NSD*3+0] = _coords[ej][ei+1].x;    el_coords[NSD*3+1] = _coords[ej][ei+1].y;
538:   return(0);
539: }

541: static PetscErrorCode AssembleA_Elasticity(Mat A,DM elas_da,DM properties_da,Vec properties)
542: {
543:   DM                     cda;
544:   Vec                    coords;
545:   DMDACoor2d             **_coords;
546:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
547:   PetscInt               sex,sey,mx,my;
548:   PetscInt               ei,ej;
549:   PetscScalar            Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
550:   PetscScalar            el_coords[NODES_PER_EL*NSD];
551:   Vec                    local_properties;
552:   GaussPointCoefficients **props;
553:   PetscScalar            *prop_E,*prop_nu;
554:   PetscErrorCode         ierr;

557:   /* setup for coords */
558:   DMGetCoordinateDM(elas_da,&cda);
559:   DMGetCoordinatesLocal(elas_da,&coords);
560:   DMDAVecGetArray(cda,coords,&_coords);

562:   /* setup for coefficients */
563:   DMCreateLocalVector(properties_da,&local_properties);
564:   DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
565:   DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
566:   DMDAVecGetArray(properties_da,local_properties,&props);

568:   DMDAGetElementsCorners(elas_da,&sex,&sey,0);
569:   DMDAGetElementsSizes(elas_da,&mx,&my,0);
570:   for (ej = sey; ej < sey+my; ej++) {
571:     for (ei = sex; ei < sex+mx; ei++) {
572:       /* get coords for the element */
573:       GetElementCoords(_coords,ei,ej,el_coords);

575:       /* get coefficients for the element */
576:       prop_E  = props[ej][ei].E;
577:       prop_nu = props[ej][ei].nu;

579:       /* initialise element stiffness matrix */
580:       PetscMemzero(Ae,sizeof(Ae));

582:       /* form element stiffness matrix */
583:       FormStressOperatorQ1(Ae,el_coords,prop_E,prop_nu);

585:       /* insert element matrix into global matrix */
586:       DMDAGetElementEqnums_u(u_eqn,ei,ej);
587:       MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
588:     }
589:   }
590:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
591:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

593:   DMDAVecRestoreArray(cda,coords,&_coords);

595:   DMDAVecRestoreArray(properties_da,local_properties,&props);
596:   VecDestroy(&local_properties);
597:   return(0);
598: }

600: static PetscErrorCode DMDASetValuesLocalStencil_ADD_VALUES(ElasticityDOF **fields_F,MatStencil u_eqn[],PetscScalar Fe_u[])
601: {
602:   PetscInt n;

605:   for (n = 0; n < 4; n++) {
606:     fields_F[u_eqn[2*n].j][u_eqn[2*n].i].ux_dof     = fields_F[u_eqn[2*n].j][u_eqn[2*n].i].ux_dof+Fe_u[2*n];
607:     fields_F[u_eqn[2*n+1].j][u_eqn[2*n+1].i].uy_dof = fields_F[u_eqn[2*n+1].j][u_eqn[2*n+1].i].uy_dof+Fe_u[2*n+1];
608:   }
609:   return(0);
610: }

612: static PetscErrorCode AssembleF_Elasticity(Vec F,DM elas_da,DM properties_da,Vec properties)
613: {
614:   DM                     cda;
615:   Vec                    coords;
616:   DMDACoor2d             **_coords;
617:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
618:   PetscInt               sex,sey,mx,my;
619:   PetscInt               ei,ej;
620:   PetscScalar            Fe[NODES_PER_EL*U_DOFS];
621:   PetscScalar            el_coords[NODES_PER_EL*NSD];
622:   Vec                    local_properties;
623:   GaussPointCoefficients **props;
624:   PetscScalar            *prop_fx,*prop_fy;
625:   Vec                    local_F;
626:   ElasticityDOF          **ff;
627:   PetscErrorCode         ierr;

630:   /* setup for coords */
631:   DMGetCoordinateDM(elas_da,&cda);
632:   DMGetCoordinatesLocal(elas_da,&coords);
633:   DMDAVecGetArray(cda,coords,&_coords);

635:   /* setup for coefficients */
636:   DMGetLocalVector(properties_da,&local_properties);
637:   DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
638:   DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
639:   DMDAVecGetArray(properties_da,local_properties,&props);

641:   /* get access to the vector */
642:   DMGetLocalVector(elas_da,&local_F);
643:   VecZeroEntries(local_F);
644:   DMDAVecGetArray(elas_da,local_F,&ff);

646:   DMDAGetElementsCorners(elas_da,&sex,&sey,0);
647:   DMDAGetElementsSizes(elas_da,&mx,&my,0);
648:   for (ej = sey; ej < sey+my; ej++) {
649:     for (ei = sex; ei < sex+mx; ei++) {
650:       /* get coords for the element */
651:       GetElementCoords(_coords,ei,ej,el_coords);

653:       /* get coefficients for the element */
654:       prop_fx = props[ej][ei].fx;
655:       prop_fy = props[ej][ei].fy;

657:       /* initialise element stiffness matrix */
658:       PetscMemzero(Fe,sizeof(Fe));

660:       /* form element stiffness matrix */
661:       FormMomentumRhsQ1(Fe,el_coords,prop_fx,prop_fy);

663:       /* insert element matrix into global matrix */
664:       DMDAGetElementEqnums_u(u_eqn,ei,ej);

666:       DMDASetValuesLocalStencil_ADD_VALUES(ff,u_eqn,Fe);
667:     }
668:   }

670:   DMDAVecRestoreArray(elas_da,local_F,&ff);
671:   DMLocalToGlobalBegin(elas_da,local_F,ADD_VALUES,F);
672:   DMLocalToGlobalEnd(elas_da,local_F,ADD_VALUES,F);
673:   DMRestoreLocalVector(elas_da,&local_F);

675:   DMDAVecRestoreArray(cda,coords,&_coords);

677:   DMDAVecRestoreArray(properties_da,local_properties,&props);
678:   DMRestoreLocalVector(properties_da,&local_properties);
679:   return(0);
680: }

682: static PetscErrorCode solve_elasticity_2d(PetscInt mx,PetscInt my)
683: {
684:   DM                     elas_da,da_prop;
685:   PetscInt               u_dof,dof,stencil_width;
686:   Mat                    A;
687:   PetscInt               mxl,myl;
688:   DM                     prop_cda,vel_cda;
689:   Vec                    prop_coords,vel_coords;
690:   PetscInt               si,sj,nx,ny,i,j,p;
691:   Vec                    f,X;
692:   PetscInt               prop_dof,prop_stencil_width;
693:   Vec                    properties,l_properties;
694:   MatNullSpace           matnull;
695:   PetscReal              dx,dy;
696:   PetscInt               M,N;
697:   DMDACoor2d             **_prop_coords,**_vel_coords;
698:   GaussPointCoefficients **element_props;
699:   KSP                    ksp_E;
700:   PetscInt               coefficient_structure = 0;
701:   PetscInt               cpu_x,cpu_y,*lx = NULL,*ly = NULL;
702:   PetscBool              use_gp_coords = PETSC_FALSE;
703:   PetscBool              use_nonsymbc  = PETSC_FALSE;
704:   PetscBool              no_view       = PETSC_FALSE;
705:   PetscBool              flg;
706:   PetscErrorCode         ierr;

709:   /* Generate the da for velocity and pressure */
710:   /*
711:    We use Q1 elements for the temperature.
712:    FEM has a 9-point stencil (BOX) or connectivity pattern
713:    Num nodes in each direction is mx+1, my+1
714:    */
715:   u_dof         = U_DOFS; /* Vx, Vy - velocities */
716:   dof           = u_dof;
717:   stencil_width = 1;
718:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx+1,my+1,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&elas_da);

720:   DMSetMatType(elas_da,MATAIJ);
721:   DMSetFromOptions(elas_da);
722:   DMSetUp(elas_da);

724:   DMDASetFieldName(elas_da,0,"Ux");
725:   DMDASetFieldName(elas_da,1,"Uy");

727:   /* unit box [0,1] x [0,1] */
728:   DMDASetUniformCoordinates(elas_da,0.0,1.0,0.0,1.0,0.0,1.0);

730:   /* Generate element properties, we will assume all material properties are constant over the element */
731:   /* local number of elements */
732:   DMDAGetElementsSizes(elas_da,&mxl,&myl,NULL);

734:   /* !!! IN PARALLEL WE MUST MAKE SURE THE TWO DMDA's ALIGN !!! */
735:   DMDAGetInfo(elas_da,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);
736:   DMDAGetElementOwnershipRanges2d(elas_da,&lx,&ly);

738:   prop_dof           = (PetscInt)(sizeof(GaussPointCoefficients)/sizeof(PetscScalar)); /* gauss point setup */
739:   prop_stencil_width = 0;
740:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx,my,cpu_x,cpu_y,prop_dof,prop_stencil_width,lx,ly,&da_prop);
741:   DMSetFromOptions(da_prop);
742:   DMSetUp(da_prop);

744:   PetscFree(lx);
745:   PetscFree(ly);

747:   /* define centroid positions */
748:   DMDAGetInfo(da_prop,0,&M,&N,0,0,0,0,0,0,0,0,0,0);
749:   dx   = 1.0/((PetscReal)(M));
750:   dy   = 1.0/((PetscReal)(N));

752:   DMDASetUniformCoordinates(da_prop,0.0+0.5*dx,1.0-0.5*dx,0.0+0.5*dy,1.0-0.5*dy,0.0,1.0);

754:   /* define coefficients */
755:   PetscOptionsGetInt(NULL,NULL,"-c_str",&coefficient_structure,NULL);

757:   DMCreateGlobalVector(da_prop,&properties);
758:   DMCreateLocalVector(da_prop,&l_properties);
759:   DMDAVecGetArray(da_prop,l_properties,&element_props);

761:   DMGetCoordinateDM(da_prop,&prop_cda);
762:   DMGetCoordinatesLocal(da_prop,&prop_coords);
763:   DMDAVecGetArray(prop_cda,prop_coords,&_prop_coords);

765:   DMDAGetGhostCorners(prop_cda,&si,&sj,0,&nx,&ny,0);

767:   DMGetCoordinateDM(elas_da,&vel_cda);
768:   DMGetCoordinatesLocal(elas_da,&vel_coords);
769:   DMDAVecGetArray(vel_cda,vel_coords,&_vel_coords);

771:   /* interpolate the coordinates */
772:   for (j = sj; j < sj+ny; j++) {
773:     for (i = si; i < si+nx; i++) {
774:       PetscInt    ngp;
775:       PetscScalar gp_xi[GAUSS_POINTS][2],gp_weight[GAUSS_POINTS];
776:       PetscScalar el_coords[8];

778:       GetElementCoords(_vel_coords,i,j,el_coords);
779:       ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

781:       for (p = 0; p < GAUSS_POINTS; p++) {
782:         PetscScalar gp_x,gp_y;
783:         PetscInt    n;
784:         PetscScalar xi_p[2],Ni_p[4];

786:         xi_p[0] = gp_xi[p][0];
787:         xi_p[1] = gp_xi[p][1];
788:         ConstructQ12D_Ni(xi_p,Ni_p);

790:         gp_x = 0.0;
791:         gp_y = 0.0;
792:         for (n = 0; n < NODES_PER_EL; n++) {
793:           gp_x = gp_x+Ni_p[n]*el_coords[2*n];
794:           gp_y = gp_y+Ni_p[n]*el_coords[2*n+1];
795:         }
796:         element_props[j][i].gp_coords[2*p]   = gp_x;
797:         element_props[j][i].gp_coords[2*p+1] = gp_y;
798:       }
799:     }
800:   }

802:   /* define the coefficients */
803:   PetscOptionsGetBool(NULL,NULL,"-use_gp_coords",&use_gp_coords,&flg);

805:   for (j = sj; j < sj+ny; j++) {
806:     for (i = si; i < si+nx; i++) {
807:       PetscScalar              centroid_x = _prop_coords[j][i].x; /* centroids of cell */
808:       PetscScalar              centroid_y = _prop_coords[j][i].y;
809:       PETSC_UNUSED PetscScalar coord_x,coord_y;

811:       if (coefficient_structure == 0) { /* isotropic */
812:         PetscScalar opts_E,opts_nu;

814:         opts_E  = 1.0;
815:         opts_nu = 0.33;
816:         PetscOptionsGetScalar(NULL,NULL,"-iso_E",&opts_E,&flg);
817:         PetscOptionsGetScalar(NULL,NULL,"-iso_nu",&opts_nu,&flg);

819:         for (p = 0; p < GAUSS_POINTS; p++) {
820:           element_props[j][i].E[p]  = opts_E;
821:           element_props[j][i].nu[p] = opts_nu;

823:           element_props[j][i].fx[p] = 0.0;
824:           element_props[j][i].fy[p] = 0.0;
825:         }
826:       } else if (coefficient_structure == 1) { /* step */
827:         PetscScalar opts_E0,opts_nu0,opts_xc;
828:         PetscScalar opts_E1,opts_nu1;

830:         opts_E0  = opts_E1  = 1.0;
831:         opts_nu0 = opts_nu1 = 0.333;
832:         opts_xc  = 0.5;
833:         PetscOptionsGetScalar(NULL,NULL,"-step_E0",&opts_E0,&flg);
834:         PetscOptionsGetScalar(NULL,NULL,"-step_nu0",&opts_nu0,&flg);
835:         PetscOptionsGetScalar(NULL,NULL,"-step_E1",&opts_E1,&flg);
836:         PetscOptionsGetScalar(NULL,NULL,"-step_nu1",&opts_nu1,&flg);
837:         PetscOptionsGetScalar(NULL,NULL,"-step_xc",&opts_xc,&flg);

839:         for (p = 0; p < GAUSS_POINTS; p++) {
840:           coord_x = centroid_x;
841:           coord_y = centroid_y;
842:           if (use_gp_coords) {
843:             coord_x = element_props[j][i].gp_coords[2*p];
844:             coord_y = element_props[j][i].gp_coords[2*p+1];
845:           }

847:           element_props[j][i].E[p]  = opts_E0;
848:           element_props[j][i].nu[p] = opts_nu0;
849:           if (PetscRealPart(coord_x) > PetscRealPart(opts_xc)) {
850:             element_props[j][i].E[p]  = opts_E1;
851:             element_props[j][i].nu[p] = opts_nu1;
852:           }

854:           element_props[j][i].fx[p] = 0.0;
855:           element_props[j][i].fy[p] = 0.0;
856:         }
857:       } else if (coefficient_structure == 2) { /* brick */
858:         PetscReal values_E[10];
859:         PetscReal values_nu[10];
860:         PetscInt  nbricks,maxnbricks;
861:         PetscInt  index,span;
862:         PetscInt  jj;

864:         flg        = PETSC_FALSE;
865:         maxnbricks = 10;
866:         PetscOptionsGetRealArray(NULL,NULL, "-brick_E",values_E,&maxnbricks,&flg);
867:         nbricks    = maxnbricks;
868:         if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply a list of E values for each brick");

870:         flg        = PETSC_FALSE;
871:         maxnbricks = 10;
872:         PetscOptionsGetRealArray(NULL,NULL, "-brick_nu",values_nu,&maxnbricks,&flg);
873:         if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply a list of nu values for each brick");
874:         if (maxnbricks != nbricks) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply equal numbers of values for E and nu");

876:         span = 1;
877:         PetscOptionsGetInt(NULL,NULL,"-brick_span",&span,&flg);

879:         /* cycle through the indices so that no two material properties are repeated in lines of x or y */
880:         jj    = (j/span)%nbricks;
881:         index = (jj+i/span)%nbricks;
882:         /*printf("j=%d: index = %d \n", j,index); */

884:         for (p = 0; p < GAUSS_POINTS; p++) {
885:           element_props[j][i].E[p]  = values_E[index];
886:           element_props[j][i].nu[p] = values_nu[index];
887:         }
888:       } else if (coefficient_structure == 3) { /* sponge */
889:         PetscScalar opts_E0,opts_nu0;
890:         PetscScalar opts_E1,opts_nu1;
891:         PetscInt    opts_t,opts_w;
892:         PetscInt    ii,jj,ci,cj;

894:         opts_E0  = opts_E1  = 1.0;
895:         opts_nu0 = opts_nu1 = 0.333;
896:         PetscOptionsGetScalar(NULL,NULL,"-sponge_E0",&opts_E0,&flg);
897:         PetscOptionsGetScalar(NULL,NULL,"-sponge_nu0",&opts_nu0,&flg);
898:         PetscOptionsGetScalar(NULL,NULL,"-sponge_E1",&opts_E1,&flg);
899:         PetscOptionsGetScalar(NULL,NULL,"-sponge_nu1",&opts_nu1,&flg);

901:         opts_t = opts_w = 1;
902:         PetscOptionsGetInt(NULL,NULL,"-sponge_t",&opts_t,&flg);
903:         PetscOptionsGetInt(NULL,NULL,"-sponge_w",&opts_w,&flg);

905:         ii = (i)/(opts_t+opts_w+opts_t);
906:         jj = (j)/(opts_t+opts_w+opts_t);

908:         ci = i - ii*(opts_t+opts_w+opts_t);
909:         cj = j - jj*(opts_t+opts_w+opts_t);

911:         for (p = 0; p < GAUSS_POINTS; p++) {
912:           element_props[j][i].E[p]  = opts_E0;
913:           element_props[j][i].nu[p] = opts_nu0;
914:         }
915:         if ((ci >= opts_t) && (ci < opts_t+opts_w)) {
916:           if ((cj >= opts_t) && (cj < opts_t+opts_w)) {
917:             for (p = 0; p < GAUSS_POINTS; p++) {
918:               element_props[j][i].E[p]  = opts_E1;
919:               element_props[j][i].nu[p] = opts_nu1;
920:             }
921:           }
922:         }
923:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown coefficient_structure");
924:     }
925:   }
926:   DMDAVecRestoreArray(prop_cda,prop_coords,&_prop_coords);

928:   DMDAVecRestoreArray(vel_cda,vel_coords,&_vel_coords);

930:   DMDAVecRestoreArray(da_prop,l_properties,&element_props);
931:   DMLocalToGlobalBegin(da_prop,l_properties,ADD_VALUES,properties);
932:   DMLocalToGlobalEnd(da_prop,l_properties,ADD_VALUES,properties);

934:   PetscOptionsGetBool(NULL,NULL,"-no_view",&no_view,NULL);
935:   if (!no_view) {
936:     DMDAViewCoefficientsGnuplot2d(da_prop,properties,"Coeffcients for elasticity eqn.","properties");
937:     DMDACoordViewGnuplot2d(elas_da,"mesh");
938:   }

940:   /* Generate a matrix with the correct non-zero pattern of type AIJ. This will work in parallel and serial */
941:   DMCreateMatrix(elas_da,&A);
942:   DMGetCoordinates(elas_da,&vel_coords);
943:   MatNullSpaceCreateRigidBody(vel_coords,&matnull);
944:   MatSetNearNullSpace(A,matnull);
945:   MatNullSpaceDestroy(&matnull);
946:   MatCreateVecs(A,&f,&X);

948:   /* assemble A11 */
949:   MatZeroEntries(A);
950:   VecZeroEntries(f);

952:   AssembleA_Elasticity(A,elas_da,da_prop,properties);
953:   /* build force vector */
954:   AssembleF_Elasticity(f,elas_da,da_prop,properties);

956:   KSPCreate(PETSC_COMM_WORLD,&ksp_E);
957:   KSPSetOptionsPrefix(ksp_E,"elas_");  /* elasticity */

959:   PetscOptionsGetBool(NULL,NULL,"-use_nonsymbc",&use_nonsymbc,&flg);
960:   /* solve */
961:   if (!use_nonsymbc) {
962:     Mat        AA;
963:     Vec        ff,XX;
964:     IS         is;
965:     VecScatter scat;

967:     DMDABCApplySymmetricCompression(elas_da,A,f,&is,&AA,&ff);
968:     VecDuplicate(ff,&XX);

970:     KSPSetOperators(ksp_E,AA,AA);
971:     KSPSetFromOptions(ksp_E);

973:     KSPSolve(ksp_E,ff,XX);

975:     /* push XX back into X */
976:     DMDABCApplyCompression(elas_da,NULL,X);

978:     VecScatterCreate(XX,NULL,X,is,&scat);
979:     VecScatterBegin(scat,XX,X,INSERT_VALUES,SCATTER_FORWARD);
980:     VecScatterEnd(scat,XX,X,INSERT_VALUES,SCATTER_FORWARD);
981:     VecScatterDestroy(&scat);

983:     MatDestroy(&AA);
984:     VecDestroy(&ff);
985:     VecDestroy(&XX);
986:     ISDestroy(&is);
987:   } else {
988:     DMDABCApplyCompression(elas_da,A,f);

990:     KSPSetOperators(ksp_E,A,A);
991:     KSPSetFromOptions(ksp_E);

993:     KSPSolve(ksp_E,f,X);
994:   }

996:   if (!no_view) {DMDAViewGnuplot2d(elas_da,X,"Displacement solution for elasticity eqn.","X");}
997:   KSPDestroy(&ksp_E);

999:   VecDestroy(&X);
1000:   VecDestroy(&f);
1001:   MatDestroy(&A);

1003:   DMDestroy(&elas_da);
1004:   DMDestroy(&da_prop);

1006:   VecDestroy(&properties);
1007:   VecDestroy(&l_properties);
1008:   return(0);
1009: }

1011: int main(int argc,char **args)
1012: {
1014:   PetscInt       mx,my;

1016:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
1017:   mx   = my = 10;
1018:   PetscOptionsGetInt(NULL,NULL,"-mx",&mx,NULL);
1019:   PetscOptionsGetInt(NULL,NULL,"-my",&my,NULL);
1020:   solve_elasticity_2d(mx,my);
1021:   PetscFinalize();
1022:   return ierr;
1023: }

1025: /* -------------------------- helpers for boundary conditions -------------------------------- */

1027: static PetscErrorCode BCApply_EAST(DM da,PetscInt d_idx,PetscScalar bc_val,Mat A,Vec b)
1028: {
1029:   DM                     cda;
1030:   Vec                    coords;
1031:   PetscInt               si,sj,nx,ny,i,j;
1032:   PetscInt               M,N;
1033:   DMDACoor2d             **_coords;
1034:   const PetscInt         *g_idx;
1035:   PetscInt               *bc_global_ids;
1036:   PetscScalar            *bc_vals;
1037:   PetscInt               nbcs;
1038:   PetscInt               n_dofs;
1039:   PetscErrorCode         ierr;
1040:   ISLocalToGlobalMapping ltogm;

1043:   /* enforce bc's */
1044:   DMGetLocalToGlobalMapping(da,&ltogm);
1045:   ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);

1047:   DMGetCoordinateDM(da,&cda);
1048:   DMGetCoordinatesLocal(da,&coords);
1049:   DMDAVecGetArray(cda,coords,&_coords);
1050:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1051:   DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);

1053:   /* --- */

1055:   PetscMalloc1(ny*n_dofs,&bc_global_ids);
1056:   PetscMalloc1(ny*n_dofs,&bc_vals);

1058:   /* init the entries to -1 so VecSetValues will ignore them */
1059:   for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;

1061:   i = nx-1;
1062:   for (j = 0; j < ny; j++) {
1063:     PetscInt                 local_id;
1064:     PETSC_UNUSED PetscScalar coordx,coordy;

1066:     local_id = i+j*nx;

1068:     bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];

1070:     coordx = _coords[j+sj][i+si].x;
1071:     coordy = _coords[j+sj][i+si].y;

1073:     bc_vals[j] =  bc_val;
1074:   }
1075:   ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1076:   nbcs = 0;
1077:   if ((si+nx) == (M)) nbcs = ny;

1079:   if (b) {
1080:     VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1081:     VecAssemblyBegin(b);
1082:     VecAssemblyEnd(b);
1083:   }
1084:   if (A) {
1085:     MatZeroRows(A,nbcs,bc_global_ids,1.0,0,0);
1086:   }

1088:   PetscFree(bc_vals);
1089:   PetscFree(bc_global_ids);

1091:   DMDAVecRestoreArray(cda,coords,&_coords);
1092:   return(0);
1093: }

1095: static PetscErrorCode BCApply_WEST(DM da,PetscInt d_idx,PetscScalar bc_val,Mat A,Vec b)
1096: {
1097:   DM                     cda;
1098:   Vec                    coords;
1099:   PetscInt               si,sj,nx,ny,i,j;
1100:   PetscInt               M,N;
1101:   DMDACoor2d             **_coords;
1102:   const PetscInt         *g_idx;
1103:   PetscInt               *bc_global_ids;
1104:   PetscScalar            *bc_vals;
1105:   PetscInt               nbcs;
1106:   PetscInt               n_dofs;
1107:   PetscErrorCode         ierr;
1108:   ISLocalToGlobalMapping ltogm;

1111:   /* enforce bc's */
1112:   DMGetLocalToGlobalMapping(da,&ltogm);
1113:   ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);

1115:   DMGetCoordinateDM(da,&cda);
1116:   DMGetCoordinatesLocal(da,&coords);
1117:   DMDAVecGetArray(cda,coords,&_coords);
1118:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1119:   DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);

1121:   /* --- */

1123:   PetscMalloc1(ny*n_dofs,&bc_global_ids);
1124:   PetscMalloc1(ny*n_dofs,&bc_vals);

1126:   /* init the entries to -1 so VecSetValues will ignore them */
1127:   for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;

1129:   i = 0;
1130:   for (j = 0; j < ny; j++) {
1131:     PetscInt                 local_id;
1132:     PETSC_UNUSED PetscScalar coordx,coordy;

1134:     local_id = i+j*nx;

1136:     bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];

1138:     coordx = _coords[j+sj][i+si].x;
1139:     coordy = _coords[j+sj][i+si].y;

1141:     bc_vals[j] =  bc_val;
1142:   }
1143:   ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1144:   nbcs = 0;
1145:   if (si == 0) nbcs = ny;

1147:   if (b) {
1148:     VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1149:     VecAssemblyBegin(b);
1150:     VecAssemblyEnd(b);
1151:   }
1152:   if (A) {
1153:     MatZeroRows(A,nbcs,bc_global_ids,1.0,0,0);
1154:   }

1156:   PetscFree(bc_vals);
1157:   PetscFree(bc_global_ids);

1159:   DMDAVecRestoreArray(cda,coords,&_coords);
1160:   return(0);
1161: }

1163: static PetscErrorCode DMDABCApplyCompression(DM elas_da,Mat A,Vec f)
1164: {

1168:   BCApply_EAST(elas_da,0,-1.0,A,f);
1169:   BCApply_EAST(elas_da,1, 0.0,A,f);
1170:   BCApply_WEST(elas_da,0,1.0,A,f);
1171:   BCApply_WEST(elas_da,1,0.0,A,f);
1172:   return(0);
1173: }

1175: static PetscErrorCode Orthogonalize(PetscInt n,Vec *vecs)
1176: {
1177:   PetscInt       i,j;
1178:   PetscScalar    dot;

1182:   for (i=0; i<n; i++) {
1183:   VecNormalize(vecs[i],NULL);
1184:      for (j=i+1; j<n; j++) {
1185:        VecDot(vecs[i],vecs[j],&dot);
1186:        VecAXPY(vecs[j],-dot,vecs[i]);
1187:      }
1188:   }
1189:   return(0);
1190: }

1192: static PetscErrorCode DMDABCApplySymmetricCompression(DM elas_da,Mat A,Vec f,IS *dofs,Mat *AA,Vec *ff)
1193: {
1195:   PetscInt       start,end,m;
1196:   PetscInt       *unconstrained;
1197:   PetscInt       cnt,i;
1198:   Vec            x;
1199:   PetscScalar    *_x;
1200:   IS             is;
1201:   VecScatter     scat;

1204:   /* push bc's into f and A */
1205:   VecDuplicate(f,&x);
1206:   BCApply_EAST(elas_da,0,-1.0,A,x);
1207:   BCApply_EAST(elas_da,1, 0.0,A,x);
1208:   BCApply_WEST(elas_da,0,1.0,A,x);
1209:   BCApply_WEST(elas_da,1,0.0,A,x);

1211:   /* define which dofs are not constrained */
1212:   VecGetLocalSize(x,&m);
1213:   PetscMalloc1(m,&unconstrained);
1214:   VecGetOwnershipRange(x,&start,&end);
1215:   VecGetArray(x,&_x);
1216:   cnt  = 0;
1217:   for (i = 0; i < m; i+=2) {
1218:     PetscReal val1,val2;

1220:     val1 = PetscRealPart(_x[i]);
1221:     val2 = PetscRealPart(_x[i+1]);
1222:     if (PetscAbs(val1) < 0.1 && PetscAbs(val2) < 0.1) {
1223:       unconstrained[cnt] = start + i;
1224:       cnt++;
1225:       unconstrained[cnt] = start + i + 1;
1226:       cnt++;
1227:     }
1228:   }
1229:   VecRestoreArray(x,&_x);

1231:   ISCreateGeneral(PETSC_COMM_WORLD,cnt,unconstrained,PETSC_COPY_VALUES,&is);
1232:   PetscFree(unconstrained);
1233:   ISSetBlockSize(is,2);

1235:   /* define correction for dirichlet in the rhs */
1236:   MatMult(A,x,f);
1237:   VecScale(f,-1.0);

1239:   /* get new matrix */
1240:   MatCreateSubMatrix(A,is,is,MAT_INITIAL_MATRIX,AA);
1241:   /* get new vector */
1242:   MatCreateVecs(*AA,NULL,ff);

1244:   VecScatterCreate(f,is,*ff,NULL,&scat);
1245:   VecScatterBegin(scat,f,*ff,INSERT_VALUES,SCATTER_FORWARD);
1246:   VecScatterEnd(scat,f,*ff,INSERT_VALUES,SCATTER_FORWARD);

1248:   {                             /* Constrain near-null space */
1249:     PetscInt     nvecs;
1250:     const        Vec *vecs;
1251:     Vec          *uvecs;
1252:     PetscBool    has_const;
1253:     MatNullSpace mnull,unull;

1255:     MatGetNearNullSpace(A,&mnull);
1256:     MatNullSpaceGetVecs(mnull,&has_const,&nvecs,&vecs);
1257:     VecDuplicateVecs(*ff,nvecs,&uvecs);
1258:     for (i=0; i<nvecs; i++) {
1259:       VecScatterBegin(scat,vecs[i],uvecs[i],INSERT_VALUES,SCATTER_FORWARD);
1260:       VecScatterEnd(scat,vecs[i],uvecs[i],INSERT_VALUES,SCATTER_FORWARD);
1261:     }
1262:     Orthogonalize(nvecs,uvecs);
1263:     MatNullSpaceCreate(PetscObjectComm((PetscObject)A),PETSC_FALSE,nvecs,uvecs,&unull);
1264:     MatSetNearNullSpace(*AA,unull);
1265:     MatNullSpaceDestroy(&unull);
1266:     VecDestroyVecs(nvecs,&uvecs);
1267:   }

1269:   VecScatterDestroy(&scat);

1271:   *dofs = is;
1272:   VecDestroy(&x);
1273:   return(0);
1274: }

1276: /*TEST

1278:    build:
1279:       requires: !complex !single

1281:    test:
1282:       args: -mx 20 -my 30 -elas_ksp_monitor_short -no_view -c_str 3 -sponge_E0 1 -sponge_E1 1000 -sponge_nu0 0.4 -sponge_nu1 0.2 -sponge_t 1 -sponge_w 8 -elas_ksp_rtol 5e-3 -elas_ksp_view
1283:       output_file: output/ex49_1.out

1285:    test:
1286:       suffix: 2
1287:       nsize: 4
1288:       args: -mx 20 -my 30 -elas_ksp_monitor_short -no_view -c_str 3 -sponge_E0 1 -sponge_E1 1000 -sponge_nu0 0.4 -sponge_nu1 0.2 -sponge_t 1 -sponge_w 8 -elas_ksp_type gcr -elas_pc_type asm -elas_sub_pc_type lu -elas_ksp_rtol 5e-3

1290:    test:
1291:       suffix: 3
1292:       nsize: 4
1293:       args: -mx 20 -my 30 -elas_ksp_monitor_short -no_view -c_str 2 -brick_E 1,10,1000,100 -brick_nu 0.4,0.2,0.3,0.1 -brick_span 3 -elas_pc_type asm -elas_sub_pc_type lu -elas_ksp_rtol 5e-3

1295:    test:
1296:       suffix: 4
1297:       nsize: 4
1298:       args: -elas_ksp_monitor_short -elas_ksp_converged_reason -elas_ksp_type cg -elas_ksp_norm_type unpreconditioned -mx 40 -my 40 -c_str 2 -brick_E 1,1e-6,1e-2 -brick_nu .3,.2,.4 -brick_span 8 -elas_mg_levels_ksp_type chebyshev -elas_pc_type ml -elas_mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.1 -elas_mg_levels_pc_type pbjacobi -elas_mg_levels_ksp_max_it 3 -use_nonsymbc -elas_pc_ml_nullspace user
1299:       requires: ml

1301:    test:
1302:       suffix: 5
1303:       nsize: 3
1304:       args: -elas_ksp_monitor_short -elas_ksp_converged_reason -elas_ksp_type cg -elas_ksp_norm_type natural -mx 22 -my 22 -c_str 2 -brick_E 1,1e-6,1e-2 -brick_nu .3,.2,.4 -brick_span 8 -elas_pc_type gamg -elas_mg_levels_ksp_type chebyshev -elas_mg_levels_ksp_max_it 1 -elas_mg_levels_ksp_chebyshev_esteig 0.2,1.1 -elas_mg_levels_pc_type jacobi

1306:    test:
1307:       suffix: 6
1308:       nsize: 4
1309:       args: -mx 20 -my 30 -elas_ksp_monitor_short -no_view -c_str 3 -sponge_E0 1 -sponge_E1 1000 -sponge_nu0 0.4 -sponge_nu1 0.2 -sponge_t 1 -sponge_w 8 -elas_ksp_type pipegcr -elas_pc_type asm -elas_sub_pc_type lu

1311:    test:
1312:       suffix: 7
1313:       nsize: 4
1314:       args: -mx 20 -my 30 -elas_ksp_monitor_short -no_view -c_str 3 -sponge_E0 1 -sponge_E1 1000 -sponge_nu0 0.4 -sponge_nu1 0.2 -sponge_t 1 -sponge_w 8 -elas_ksp_type pipegcr -elas_pc_type asm -elas_sub_pc_type ksp -elas_sub_ksp_ksp_type cg -elas_sub_ksp_ksp_max_it 15

1316:    test:
1317:       suffix: 8
1318:       nsize: 4
1319:       args: -mx 20 -my 30 -elas_ksp_monitor_short -no_view -c_str 3 -sponge_E0 1 -sponge_E1 1000 -sponge_nu0 0.4 -sponge_nu1 0.2 -sponge_t 1 -sponge_w 8 -elas_ksp_type pipefgmres -elas_pc_type asm -elas_sub_pc_type ksp -elas_sub_ksp_ksp_type cg -elas_sub_ksp_ksp_max_it 15

1321:    test:
1322:       suffix: hypre_nullspace
1323:       requires: hypre !defined(PETSC_HAVE_HYPRE_DEVICE)
1324:       args: -elas_ksp_monitor_short -elas_ksp_converged_reason -elas_ksp_type cg -elas_ksp_norm_type natural -mx 22 -my 22 -c_str 2 -brick_E 1,1e-6,1e-2 -brick_nu .3,.2,.4 -brick_span 8 -elas_pc_type hypre -elas_pc_hypre_boomeramg_nodal_coarsen 6 -elas_pc_hypre_boomeramg_vec_interp_variant 3 -elas_pc_hypre_boomeramg_interp_type ext+i -elas_ksp_view

1326:    test:
1327:       nsize: 4
1328:       suffix: bddc
1329:       args: -elas_ksp_monitor_short -no_view -elas_ksp_converged_reason -elas_ksp_type cg -elas_ksp_norm_type natural -mx 22 -my 22 -dm_mat_type is -elas_pc_type bddc -elas_pc_bddc_monolithic

1331:    test:
1332:       nsize: 4
1333:       suffix: bddc_unsym
1334:       args: -elas_ksp_monitor_short -no_view -elas_ksp_converged_reason -elas_ksp_type cg -elas_ksp_norm_type natural -mx 22 -my 22 -dm_mat_type is -elas_pc_type bddc -elas_pc_bddc_monolithic -use_nonsymbc -elas_pc_bddc_symmetric 0

1336:    test:
1337:       nsize: 4
1338:       suffix: bddc_unsym_deluxe
1339:       args: -elas_ksp_monitor_short -no_view -elas_ksp_converged_reason -elas_ksp_type cg -elas_ksp_norm_type natural -mx 22 -my 22 -dm_mat_type is -elas_pc_type bddc -elas_pc_bddc_monolithic -use_nonsymbc -elas_pc_bddc_symmetric 0 -elas_pc_bddc_use_deluxe_scaling -elas_sub_schurs_symmetric 0

1341:    test:
1342:       nsize: 4
1343:       suffix: fetidp_unsym_deluxe
1344:       args: -elas_ksp_monitor_short -no_view -elas_ksp_converged_reason -elas_ksp_type fetidp -elas_fetidp_ksp_type cg -elas_ksp_norm_type natural -mx 22 -my 22 -dm_mat_type is -elas_fetidp_bddc_pc_bddc_monolithic -use_nonsymbc -elas_fetidp_bddc_pc_bddc_use_deluxe_scaling -elas_fetidp_bddc_sub_schurs_symmetric 0 -elas_fetidp_bddc_pc_bddc_deluxe_singlemat

1346:    test:
1347:       nsize: 4
1348:       suffix: bddc_layerjump
1349:       args: -mx 40 -my 40 -elas_ksp_monitor_short -no_view -c_str 3 -sponge_E0 1 -sponge_E1 1000 -sponge_nu0 0.4 -sponge_nu1 0.2 -sponge_t 1 -sponge_w 8 -elas_ksp_type cg -elas_pc_type bddc -elas_pc_bddc_monolithic -dm_mat_type is -elas_ksp_norm_type natural

1351:    test:
1352:       nsize: 4
1353:       suffix: bddc_subdomainjump
1354:       args: -mx 40 -my 40 -elas_ksp_monitor_short -no_view -c_str 2 -brick_E 1,1000 -brick_nu 0.4,0.2 -brick_span 20  -elas_ksp_type cg -elas_pc_type bddc -elas_pc_bddc_monolithic -dm_mat_type is -elas_pc_is_use_stiffness_scaling -elas_ksp_norm_type natural

1356:    test:
1357:       nsize: 9
1358:       suffix: bddc_subdomainjump_deluxe
1359:       args: -mx 30 -my 30 -elas_ksp_monitor_short -no_view -c_str 2 -brick_E 1,1000 -brick_nu 0.4,0.2 -brick_span 10  -elas_ksp_type cg -elas_pc_type bddc -elas_pc_bddc_monolithic -dm_mat_type is -elas_pc_bddc_use_deluxe_scaling -elas_ksp_norm_type natural -elas_pc_bddc_schur_layers 1
1360: TEST*/