Actual source code: ex49.c
petsc-3.4.5 2014-06-29
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 moduls) and nu (Poisson ratio) may vary as a function of space. \n\
3: The model utilisse boundary conditions which produce compression in the x direction. \n\
4: Options: \n\
5: -mx : number elements in x-direciton \n\
6: -my : number elements in y-direciton \n\
7: -c_str : indicates the structure of the coefficients to use. \n\
8: -c_str 0 => Setup for an isotropic material with constant coefficients. \n\
9: Parameters: \n\
10: -iso_E : Youngs modulus \n\
11: -iso_nu : Poisson ratio \n\
12: -c_str 1 => Setup for a step function in the material properties in x. \n\
13: Parameters: \n\
14: -step_E0 : Youngs modulus to the left of the step \n\
15: -step_nu0 : Poisson ratio to the left of the step \n\
16: -step_E1 : Youngs modulus to the right of the step \n\
17: -step_n1 : Poisson ratio to the right of the step \n\
18: -step_xc : x coordinate of the step \n\
19: -c_str 2 => Setup for a checkerboard material with alternating properties. \n\
20: Repeats the following pattern throughout the domain. For example with 4 materials specified, we would heve \n\
21: -------------------------\n\
22: | D | A | B | C |\n\
23: ------|-----|-----|------\n\
24: | C | D | A | B |\n\
25: ------|-----|-----|------\n\
26: | B | C | D | A |\n\
27: ------|-----|-----|------\n\
28: | A | B | C | D |\n\
29: -------------------------\n\
30: \n\
31: Parameters: \n\
32: -brick_E : a comma seperated list of Young's modulii \n\
33: -brick_nu : a comma seperated list of Poisson ratio's \n\
34: -brick_span : the number of elements in x and y each brick will span \n\
35: -c_str 3 => Setup for a sponge-like material with alternating properties. \n\
36: Repeats the following pattern throughout the domain \n\
37: -----------------------------\n\
38: | [background] |\n\
39: | E0,nu0 |\n\
40: | ----------------- |\n\
41: | | [inclusion] | |\n\
42: | | E1,nu1 | |\n\
43: | | | |\n\
44: | | <---- w ----> | |\n\
45: | | | |\n\
46: | | | |\n\
47: | ----------------- |\n\
48: | |\n\
49: | |\n\
50: -----------------------------\n\
51: <-------- t + w + t ------->\n\
52: \n\
53: Parameters: \n\
54: -sponge_E0 : Youngs moduls of the surrounding material \n\
55: -sponge_E1 : Youngs moduls of the inclusio \n\
56: -sponge_nu0 : Poisson ratio of the surrounding material \n\
57: -sponge_nu1 : Poisson ratio of the inclusio \n\
58: -sponge_t : the number of elements defining the border around each inclusion \n\
59: -sponge_w : the number of elements in x and y each inclusion will span\n\
60: -use_gp_coords : Evaluate the Youngs modulus, Poisson ratio and the body force at the global coordinates of the quadrature points.\n\
61: By default, E, nu and the body force are evaulated at the element center and applied as a constant over the entire element.\n\
62: -use_nonsymbc : Option to use non-symmetric boundary condition imposition. This choice will use less memory.";
64: /* Contributed by Dave May */
66: #include <petscksp.h>
67: #include <petscdmda.h>
69: static PetscErrorCode DMDABCApplyCompression(DM,Mat,Vec);
70: static PetscErrorCode DMDABCApplySymmetricCompression(DM elas_da,Mat A,Vec f,IS *dofs,Mat *AA,Vec *ff);
73: #define NSD 2 /* number of spatial dimensions */
74: #define NODES_PER_EL 4 /* nodes per element */
75: #define U_DOFS 2 /* degrees of freedom per displacement node */
76: #define GAUSS_POINTS 4
78: /* cell based evaluation */
79: typedef struct {
80: PetscScalar E,nu,fx,fy;
81: } Coefficients;
83: /* Gauss point based evaluation 8+4+4+4 = 20 */
84: typedef struct {
85: PetscScalar gp_coords[2*GAUSS_POINTS];
86: PetscScalar E[GAUSS_POINTS];
87: PetscScalar nu[GAUSS_POINTS];
88: PetscScalar fx[GAUSS_POINTS];
89: PetscScalar fy[GAUSS_POINTS];
90: } GaussPointCoefficients;
92: typedef struct {
93: PetscScalar ux_dof;
94: PetscScalar uy_dof;
95: } ElasticityDOF;
98: /*
100: D = E/((1+nu)(1-2nu)) * [ 1-nu nu 0 ]
101: [ nu 1-nu 0 ]
102: [ 0 0 0.5*(1-2nu) ]
104: B = [ d_dx 0 ]
105: [ 0 d_dy ]
106: [ d_dy d_dx ]
108: */
110: /* FEM routines */
111: /*
112: Element: Local basis function ordering
113: 1-----2
114: | |
115: | |
116: 0-----3
117: */
118: static void ConstructQ12D_Ni(PetscScalar _xi[],PetscScalar Ni[])
119: {
120: PetscScalar xi = _xi[0];
121: PetscScalar eta = _xi[1];
123: Ni[0] = 0.25*(1.0-xi)*(1.0-eta);
124: Ni[1] = 0.25*(1.0-xi)*(1.0+eta);
125: Ni[2] = 0.25*(1.0+xi)*(1.0+eta);
126: Ni[3] = 0.25*(1.0+xi)*(1.0-eta);
127: }
129: static void ConstructQ12D_GNi(PetscScalar _xi[],PetscScalar GNi[][NODES_PER_EL])
130: {
131: PetscScalar xi = _xi[0];
132: PetscScalar eta = _xi[1];
134: GNi[0][0] = -0.25*(1.0-eta);
135: GNi[0][1] = -0.25*(1.0+eta);
136: GNi[0][2] = 0.25*(1.0+eta);
137: GNi[0][3] = 0.25*(1.0-eta);
139: GNi[1][0] = -0.25*(1.0-xi);
140: GNi[1][1] = 0.25*(1.0-xi);
141: GNi[1][2] = 0.25*(1.0+xi);
142: GNi[1][3] = -0.25*(1.0+xi);
143: }
145: static void ConstructQ12D_GNx(PetscScalar GNi[][NODES_PER_EL],PetscScalar GNx[][NODES_PER_EL],PetscScalar coords[],PetscScalar *det_J)
146: {
147: PetscScalar J00,J01,J10,J11,J;
148: PetscScalar iJ00,iJ01,iJ10,iJ11;
149: PetscInt i;
151: J00 = J01 = J10 = J11 = 0.0;
152: for (i = 0; i < NODES_PER_EL; i++) {
153: PetscScalar cx = coords[2*i+0];
154: PetscScalar cy = coords[2*i+1];
156: J00 = J00+GNi[0][i]*cx; /* J_xx = dx/dxi */
157: J01 = J01+GNi[0][i]*cy; /* J_xy = dy/dxi */
158: J10 = J10+GNi[1][i]*cx; /* J_yx = dx/deta */
159: J11 = J11+GNi[1][i]*cy; /* J_yy = dy/deta */
160: }
161: J = (J00*J11)-(J01*J10);
163: iJ00 = J11/J;
164: iJ01 = -J01/J;
165: iJ10 = -J10/J;
166: iJ11 = J00/J;
169: for (i = 0; i < NODES_PER_EL; i++) {
170: GNx[0][i] = GNi[0][i]*iJ00+GNi[1][i]*iJ01;
171: GNx[1][i] = GNi[0][i]*iJ10+GNi[1][i]*iJ11;
172: }
174: if (det_J != NULL) *det_J = J;
175: }
177: static void ConstructGaussQuadrature(PetscInt *ngp,PetscScalar gp_xi[][2],PetscScalar gp_weight[])
178: {
179: *ngp = 4;
180: gp_xi[0][0] = -0.57735026919;gp_xi[0][1] = -0.57735026919;
181: gp_xi[1][0] = -0.57735026919;gp_xi[1][1] = 0.57735026919;
182: gp_xi[2][0] = 0.57735026919;gp_xi[2][1] = 0.57735026919;
183: gp_xi[3][0] = 0.57735026919;gp_xi[3][1] = -0.57735026919;
184: gp_weight[0] = 1.0;
185: gp_weight[1] = 1.0;
186: gp_weight[2] = 1.0;
187: gp_weight[3] = 1.0;
188: }
191: /* procs to the left claim the ghost node as their element */
194: static PetscErrorCode DMDAGetLocalElementSize(DM da,PetscInt *mxl,PetscInt *myl,PetscInt *mzl)
195: {
197: PetscInt m,n,p,M,N,P;
198: PetscInt sx,sy,sz;
201: DMDAGetInfo(da,0,&M,&N,&P,0,0,0,0,0,0,0,0,0);
202: DMDAGetCorners(da,&sx,&sy,&sz,&m,&n,&p);
204: if (mxl != NULL) {
205: *mxl = m;
206: if ((sx+m) == M) *mxl = m-1; /* last proc */
207: }
208: if (myl != NULL) {
209: *myl = n;
210: if ((sy+n) == N) *myl = n-1; /* last proc */
211: }
212: if (mzl != NULL) {
213: *mzl = p;
214: if ((sz+p) == P) *mzl = p-1; /* last proc */
215: }
216: return(0);
217: }
221: static PetscErrorCode DMDAGetElementCorners(DM da,PetscInt *sx,PetscInt *sy,PetscInt *sz,PetscInt *mx,PetscInt *my,PetscInt *mz)
222: {
224: PetscInt si,sj,sk;
227: DMDAGetGhostCorners(da,&si,&sj,&sk,0,0,0);
229: if (sx) {
230: *sx = si;
231: if (si != 0) *sx = si+1;
232: }
233: if (sy) {
234: *sy = sj;
235: if (sj != 0) *sy = sj+1;
236: }
238: if (sk) {
239: *sz = sk;
240: if (sk != 0) *sz = sk+1;
241: }
243: DMDAGetLocalElementSize(da,mx,my,mz);
244: return(0);
245: }
249: static PetscErrorCode DMDAGetElementOwnershipRanges2d(DM da,PetscInt **_lx,PetscInt **_ly)
250: {
252: PetscMPIInt rank;
253: PetscInt proc_I,proc_J;
254: PetscInt cpu_x,cpu_y;
255: PetscInt local_mx,local_my;
256: Vec vlx,vly;
257: PetscInt *LX,*LY,i;
258: PetscScalar *_a;
259: Vec V_SEQ;
260: VecScatter ctx;
263: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
265: DMDAGetInfo(da,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);
267: proc_J = rank/cpu_x;
268: proc_I = rank-cpu_x*proc_J;
270: PetscMalloc(sizeof(PetscInt)*cpu_x,&LX);
271: PetscMalloc(sizeof(PetscInt)*cpu_y,&LY);
273: DMDAGetLocalElementSize(da,&local_mx,&local_my,NULL);
274: VecCreate(PETSC_COMM_WORLD,&vlx);
275: VecSetSizes(vlx,PETSC_DECIDE,cpu_x);
276: VecSetFromOptions(vlx);
278: VecCreate(PETSC_COMM_WORLD,&vly);
279: VecSetSizes(vly,PETSC_DECIDE,cpu_y);
280: VecSetFromOptions(vly);
282: VecSetValue(vlx,proc_I,(PetscScalar)(local_mx+1.0e-9),INSERT_VALUES);
283: VecSetValue(vly,proc_J,(PetscScalar)(local_my+1.0e-9),INSERT_VALUES);
284: VecAssemblyBegin(vlx);VecAssemblyEnd(vlx);
285: VecAssemblyBegin(vly);VecAssemblyEnd(vly);
289: VecScatterCreateToAll(vlx,&ctx,&V_SEQ);
290: VecScatterBegin(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
291: VecScatterEnd(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
292: VecGetArray(V_SEQ,&_a);
293: for (i = 0; i < cpu_x; i++) LX[i] = (PetscInt)PetscRealPart(_a[i]);
294: VecRestoreArray(V_SEQ,&_a);
295: VecScatterDestroy(&ctx);
296: VecDestroy(&V_SEQ);
298: VecScatterCreateToAll(vly,&ctx,&V_SEQ);
299: VecScatterBegin(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
300: VecScatterEnd(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
301: VecGetArray(V_SEQ,&_a);
302: for (i = 0; i < cpu_y; i++) LY[i] = (PetscInt)PetscRealPart(_a[i]);
303: VecRestoreArray(V_SEQ,&_a);
304: VecScatterDestroy(&ctx);
305: VecDestroy(&V_SEQ);
307: *_lx = LX;
308: *_ly = LY;
310: VecDestroy(&vlx);
311: VecDestroy(&vly);
312: return(0);
313: }
317: static PetscErrorCode DMDACoordViewGnuplot2d(DM da,const char prefix[])
318: {
319: DM cda;
320: Vec coords;
321: DMDACoor2d **_coords;
322: PetscInt si,sj,nx,ny,i,j;
323: FILE *fp;
324: char fname[PETSC_MAX_PATH_LEN];
325: PetscMPIInt rank;
329: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
330: PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
331: PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
332: if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");
333: PetscFPrintf(PETSC_COMM_SELF,fp,"### Element geometry for processor %1.4d ### \n",rank);
335: DMGetCoordinateDM(da,&cda);
336: DMGetCoordinatesLocal(da,&coords);
337: DMDAVecGetArray(cda,coords,&_coords);
338: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
339: for (j = sj; j < sj+ny-1; j++) {
340: for (i = si; i < si+nx-1; i++) {
341: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j][i].x),PetscRealPart(_coords[j][i].y));
342: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j+1][i].x),PetscRealPart(_coords[j+1][i].y));
343: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j+1][i+1].x),PetscRealPart(_coords[j+1][i+1].y));
344: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j][i+1].x),PetscRealPart(_coords[j][i+1].y));
345: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n\n",PetscRealPart(_coords[j][i].x),PetscRealPart(_coords[j][i].y));
346: }
347: }
348: DMDAVecRestoreArray(cda,coords,&_coords);
350: PetscFClose(PETSC_COMM_SELF,fp);
351: return(0);
352: }
356: static PetscErrorCode DMDAViewGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
357: {
358: DM cda;
359: Vec coords,local_fields;
360: DMDACoor2d **_coords;
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;
366: PetscInt n_dofs,d;
367: PetscScalar *_fields;
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");
386: DMGetCoordinateDM(da,&cda);
387: DMGetCoordinatesLocal(da,&coords);
388: DMDAVecGetArray(cda,coords,&_coords);
389: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
391: DMCreateLocalVector(da,&local_fields);
392: DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
393: DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
394: VecGetArray(local_fields,&_fields);
396: for (j = sj; j < sj+ny; j++) {
397: for (i = si; i < si+nx; i++) {
398: PetscScalar coord_x,coord_y;
399: PetscScalar field_d;
401: coord_x = _coords[j][i].x;
402: coord_y = _coords[j][i].y;
404: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",PetscRealPart(coord_x),PetscRealPart(coord_y));
405: for (d = 0; d < n_dofs; d++) {
406: field_d = _fields[n_dofs*((i-si)+(j-sj)*(nx))+d];
407: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e ",PetscRealPart(field_d));
408: }
409: PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
410: }
411: }
412: VecRestoreArray(local_fields,&_fields);
413: VecDestroy(&local_fields);
415: DMDAVecRestoreArray(cda,coords,&_coords);
417: PetscFClose(PETSC_COMM_SELF,fp);
418: return(0);
419: }
423: static PetscErrorCode DMDAViewCoefficientsGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
424: {
425: DM cda;
426: Vec local_fields;
427: FILE *fp;
428: char fname[PETSC_MAX_PATH_LEN];
429: const char *field_name;
430: PetscMPIInt rank;
431: PetscInt si,sj,nx,ny,i,j,p;
432: PetscInt n_dofs,d;
433: GaussPointCoefficients **_coefficients;
434: PetscErrorCode ierr;
437: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
438: PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
439: PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
440: if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");
442: PetscFPrintf(PETSC_COMM_SELF,fp,"### %s (processor %1.4d) ### \n",comment,rank);
443: DMDAGetInfo(da,0,0,0,0,0,0,0,&n_dofs,0,0,0,0,0);
444: PetscFPrintf(PETSC_COMM_SELF,fp,"### x y ");
445: for (d = 0; d < n_dofs; d++) {
446: DMDAGetFieldName(da,d,&field_name);
447: PetscFPrintf(PETSC_COMM_SELF,fp,"%s ",field_name);
448: }
449: PetscFPrintf(PETSC_COMM_SELF,fp,"###\n");
452: DMGetCoordinateDM(da,&cda);
453: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
455: DMCreateLocalVector(da,&local_fields);
456: DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
457: DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
458: DMDAVecGetArray(da,local_fields,&_coefficients);
460: for (j = sj; j < sj+ny; j++) {
461: for (i = si; i < si+nx; i++) {
462: PetscScalar coord_x,coord_y;
464: for (p = 0; p < GAUSS_POINTS; p++) {
465: coord_x = _coefficients[j][i].gp_coords[2*p];
466: coord_y = _coefficients[j][i].gp_coords[2*p+1];
468: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",PetscRealPart(coord_x),PetscRealPart(coord_y));
470: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e %1.6e %1.6e",
471: PetscRealPart(_coefficients[j][i].E[p]),PetscRealPart(_coefficients[j][i].nu[p]),
472: PetscRealPart(_coefficients[j][i].fx[p]),PetscRealPart(_coefficients[j][i].fy[p]));
473: PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
474: }
475: }
476: }
477: DMDAVecRestoreArray(da,local_fields,&_coefficients);
478: VecDestroy(&local_fields);
480: PetscFClose(PETSC_COMM_SELF,fp);
481: return(0);
482: }
484: static void FormStressOperatorQ1(PetscScalar Ke[],PetscScalar coords[],PetscScalar E[],PetscScalar nu[])
485: {
486: PetscInt ngp;
487: PetscScalar gp_xi[GAUSS_POINTS][2];
488: PetscScalar gp_weight[GAUSS_POINTS];
489: PetscInt p,i,j,k,l;
490: PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
491: PetscScalar J_p;
492: PetscScalar B[3][U_DOFS*NODES_PER_EL];
493: PetscScalar prop_E,prop_nu,factor,constit_D[3][3];
495: /* define quadrature rule */
496: ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);
498: /* evaluate integral */
499: for (p = 0; p < ngp; p++) {
500: ConstructQ12D_GNi(gp_xi[p],GNi_p);
501: ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
503: for (i = 0; i < NODES_PER_EL; i++) {
504: PetscScalar d_dx_i = GNx_p[0][i];
505: PetscScalar d_dy_i = GNx_p[1][i];
507: B[0][2*i] = d_dx_i; B[0][2*i+1] = 0.0;
508: B[1][2*i] = 0.0; B[1][2*i+1] = d_dy_i;
509: B[2][2*i] = d_dy_i; B[2][2*i+1] = d_dx_i;
510: }
512: /* form D for the quadrature point */
513: prop_E = E[p];
514: prop_nu = nu[p];
515: factor = prop_E / ((1.0+prop_nu)*(1.0-2.0*prop_nu));
516: constit_D[0][0] = 1.0-prop_nu; constit_D[0][1] = prop_nu; constit_D[0][2] = 0.0;
517: constit_D[1][0] = prop_nu; constit_D[1][1] = 1.0-prop_nu; constit_D[1][2] = 0.0;
518: constit_D[2][0] = 0.0; constit_D[2][1] = 0.0; constit_D[2][2] = 0.5*(1.0-2.0*prop_nu);
519: for (i = 0; i < 3; i++) {
520: for (j = 0; j < 3; j++) {
521: constit_D[i][j] = factor * constit_D[i][j] * gp_weight[p] * J_p;
522: }
523: }
525: /* form Bt tildeD B */
526: /*
527: Ke_ij = Bt_ik . D_kl . B_lj
528: = B_ki . D_kl . B_lj
529: */
530: for (i = 0; i < 8; i++) {
531: for (j = 0; j < 8; j++) {
532: for (k = 0; k < 3; k++) {
533: for (l = 0; l < 3; l++) {
534: Ke[8*i+j] = Ke[8*i+j] + B[k][i] * constit_D[k][l] * B[l][j];
535: }
536: }
537: }
538: }
540: } /* end quadrature */
541: }
543: static void FormMomentumRhsQ1(PetscScalar Fe[],PetscScalar coords[],PetscScalar fx[],PetscScalar fy[])
544: {
545: PetscInt ngp;
546: PetscScalar gp_xi[GAUSS_POINTS][2];
547: PetscScalar gp_weight[GAUSS_POINTS];
548: PetscInt p,i;
549: PetscScalar Ni_p[NODES_PER_EL];
550: PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
551: PetscScalar J_p,fac;
553: /* define quadrature rule */
554: ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);
556: /* evaluate integral */
557: for (p = 0; p < ngp; p++) {
558: ConstructQ12D_Ni(gp_xi[p],Ni_p);
559: ConstructQ12D_GNi(gp_xi[p],GNi_p);
560: ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
561: fac = gp_weight[p]*J_p;
563: for (i = 0; i < NODES_PER_EL; i++) {
564: Fe[NSD*i] += fac*Ni_p[i]*fx[p];
565: Fe[NSD*i+1] += fac*Ni_p[i]*fy[p];
566: }
567: }
568: }
570: /*
571: i,j are the element indices
572: The unknown is a vector quantity.
573: The s[].c is used to indicate the degree of freedom.
574: */
577: static PetscErrorCode DMDAGetElementEqnums_u(MatStencil s_u[],PetscInt i,PetscInt j)
578: {
580: /* displacement */
581: /* node 0 */
582: s_u[0].i = i;s_u[0].j = j;s_u[0].c = 0; /* Ux0 */
583: s_u[1].i = i;s_u[1].j = j;s_u[1].c = 1; /* Uy0 */
585: /* node 1 */
586: s_u[2].i = i;s_u[2].j = j+1;s_u[2].c = 0; /* Ux1 */
587: s_u[3].i = i;s_u[3].j = j+1;s_u[3].c = 1; /* Uy1 */
589: /* node 2 */
590: s_u[4].i = i+1;s_u[4].j = j+1;s_u[4].c = 0; /* Ux2 */
591: s_u[5].i = i+1;s_u[5].j = j+1;s_u[5].c = 1; /* Uy2 */
593: /* node 3 */
594: s_u[6].i = i+1;s_u[6].j = j;s_u[6].c = 0; /* Ux3 */
595: s_u[7].i = i+1;s_u[7].j = j;s_u[7].c = 1; /* Uy3 */
596: return(0);
597: }
601: static PetscErrorCode GetElementCoords(DMDACoor2d **_coords,PetscInt ei,PetscInt ej,PetscScalar el_coords[])
602: {
604: /* get coords for the element */
605: el_coords[NSD*0+0] = _coords[ej][ei].x; el_coords[NSD*0+1] = _coords[ej][ei].y;
606: el_coords[NSD*1+0] = _coords[ej+1][ei].x; el_coords[NSD*1+1] = _coords[ej+1][ei].y;
607: el_coords[NSD*2+0] = _coords[ej+1][ei+1].x; el_coords[NSD*2+1] = _coords[ej+1][ei+1].y;
608: el_coords[NSD*3+0] = _coords[ej][ei+1].x; el_coords[NSD*3+1] = _coords[ej][ei+1].y;
609: return(0);
610: }
614: static PetscErrorCode AssembleA_Elasticity(Mat A,DM elas_da,DM properties_da,Vec properties)
615: {
616: DM cda;
617: Vec coords;
618: DMDACoor2d **_coords;
619: MatStencil u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
620: PetscInt sex,sey,mx,my;
621: PetscInt ei,ej;
622: PetscScalar Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
623: PetscScalar el_coords[NODES_PER_EL*NSD];
624: Vec local_properties;
625: GaussPointCoefficients **props;
626: PetscScalar *prop_E,*prop_nu;
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: DMCreateLocalVector(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: DMDAGetElementCorners(elas_da,&sex,&sey,0,&mx,&my,0);
642: for (ej = sey; ej < sey+my; ej++) {
643: for (ei = sex; ei < sex+mx; ei++) {
644: /* get coords for the element */
645: GetElementCoords(_coords,ei,ej,el_coords);
647: /* get coefficients for the element */
648: prop_E = props[ej][ei].E;
649: prop_nu = props[ej][ei].nu;
651: /* initialise element stiffness matrix */
652: PetscMemzero(Ae,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS);
654: /* form element stiffness matrix */
655: FormStressOperatorQ1(Ae,el_coords,prop_E,prop_nu);
657: /* insert element matrix into global matrix */
658: DMDAGetElementEqnums_u(u_eqn,ei,ej);
659: MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
660: }
661: }
662: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
663: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
665: DMDAVecRestoreArray(cda,coords,&_coords);
667: DMDAVecRestoreArray(properties_da,local_properties,&props);
668: VecDestroy(&local_properties);
669: return(0);
670: }
675: static PetscErrorCode DMDASetValuesLocalStencil_ADD_VALUES(ElasticityDOF **fields_F,MatStencil u_eqn[],PetscScalar Fe_u[])
676: {
677: PetscInt n;
680: for (n = 0; n < 4; n++) {
681: 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];
682: 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];
683: }
684: return(0);
685: }
689: static PetscErrorCode AssembleF_Elasticity(Vec F,DM elas_da,DM properties_da,Vec properties)
690: {
691: DM cda;
692: Vec coords;
693: DMDACoor2d **_coords;
694: MatStencil u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
695: PetscInt sex,sey,mx,my;
696: PetscInt ei,ej;
697: PetscScalar Fe[NODES_PER_EL*U_DOFS];
698: PetscScalar el_coords[NODES_PER_EL*NSD];
699: Vec local_properties;
700: GaussPointCoefficients **props;
701: PetscScalar *prop_fx,*prop_fy;
702: Vec local_F;
703: ElasticityDOF **ff;
704: PetscErrorCode ierr;
707: /* setup for coords */
708: DMGetCoordinateDM(elas_da,&cda);
709: DMGetCoordinatesLocal(elas_da,&coords);
710: DMDAVecGetArray(cda,coords,&_coords);
712: /* setup for coefficients */
713: DMGetLocalVector(properties_da,&local_properties);
714: DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
715: DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
716: DMDAVecGetArray(properties_da,local_properties,&props);
718: /* get acces to the vector */
719: DMGetLocalVector(elas_da,&local_F);
720: VecZeroEntries(local_F);
721: DMDAVecGetArray(elas_da,local_F,&ff);
724: DMDAGetElementCorners(elas_da,&sex,&sey,0,&mx,&my,0);
725: for (ej = sey; ej < sey+my; ej++) {
726: for (ei = sex; ei < sex+mx; ei++) {
727: /* get coords for the element */
728: GetElementCoords(_coords,ei,ej,el_coords);
730: /* get coefficients for the element */
731: prop_fx = props[ej][ei].fx;
732: prop_fy = props[ej][ei].fy;
734: /* initialise element stiffness matrix */
735: PetscMemzero(Fe,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS);
737: /* form element stiffness matrix */
738: FormMomentumRhsQ1(Fe,el_coords,prop_fx,prop_fy);
740: /* insert element matrix into global matrix */
741: DMDAGetElementEqnums_u(u_eqn,ei,ej);
743: DMDASetValuesLocalStencil_ADD_VALUES(ff,u_eqn,Fe);
744: }
745: }
747: DMDAVecRestoreArray(elas_da,local_F,&ff);
748: DMLocalToGlobalBegin(elas_da,local_F,ADD_VALUES,F);
749: DMLocalToGlobalEnd(elas_da,local_F,ADD_VALUES,F);
750: DMRestoreLocalVector(elas_da,&local_F);
752: DMDAVecRestoreArray(cda,coords,&_coords);
754: DMDAVecRestoreArray(properties_da,local_properties,&props);
755: DMRestoreLocalVector(properties_da,&local_properties);
756: return(0);
757: }
761: static PetscErrorCode solve_elasticity_2d(PetscInt mx,PetscInt my)
762: {
763: DM elas_da,da_prop;
764: PetscInt u_dof,dof,stencil_width;
765: Mat A;
766: PetscInt mxl,myl;
767: DM prop_cda,vel_cda;
768: Vec prop_coords,vel_coords;
769: PetscInt si,sj,nx,ny,i,j,p;
770: Vec f,X;
771: PetscInt prop_dof,prop_stencil_width;
772: Vec properties,l_properties;
773: MatNullSpace matnull;
774: PetscReal dx,dy;
775: PetscInt M,N;
776: DMDACoor2d **_prop_coords,**_vel_coords;
777: GaussPointCoefficients **element_props;
778: KSP ksp_E;
779: PetscInt coefficient_structure = 0;
780: PetscInt cpu_x,cpu_y,*lx = NULL,*ly = NULL;
781: PetscBool use_gp_coords = PETSC_FALSE;
782: PetscBool use_nonsymbc = PETSC_FALSE;
783: PetscBool no_view = PETSC_FALSE;
784: PetscBool flg;
785: PetscErrorCode ierr;
788: /* Generate the da for velocity and pressure */
789: /*
790: We use Q1 elements for the temperature.
791: FEM has a 9-point stencil (BOX) or connectivity pattern
792: Num nodes in each direction is mx+1, my+1
793: */
794: u_dof = U_DOFS; /* Vx, Vy - velocities */
795: dof = u_dof;
796: stencil_width = 1;
797: DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,
798: mx+1,my+1,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&elas_da);
799: DMDASetFieldName(elas_da,0,"Ux");
800: DMDASetFieldName(elas_da,1,"Uy");
802: /* unit box [0,1] x [0,1] */
803: DMDASetUniformCoordinates(elas_da,0.0,1.0,0.0,1.0,0.0,1.0);
806: /* Generate element properties, we will assume all material properties are constant over the element */
807: /* local number of elements */
808: DMDAGetLocalElementSize(elas_da,&mxl,&myl,NULL);
810: /* !!! IN PARALLEL WE MUST MAKE SURE THE TWO DMDA's ALIGN !!! */
811: DMDAGetInfo(elas_da,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);
812: DMDAGetElementOwnershipRanges2d(elas_da,&lx,&ly);
814: prop_dof = (PetscInt)(sizeof(GaussPointCoefficients)/sizeof(PetscScalar)); /* gauss point setup */
815: prop_stencil_width = 0;
816: DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,
817: mx,my,cpu_x,cpu_y,prop_dof,prop_stencil_width,lx,ly,&da_prop);
818: PetscFree(lx);
819: PetscFree(ly);
821: /* define centroid positions */
822: DMDAGetInfo(da_prop,0,&M,&N,0,0,0,0,0,0,0,0,0,0);
823: dx = 1.0/((PetscReal)(M));
824: dy = 1.0/((PetscReal)(N));
826: 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);
828: /* define coefficients */
829: PetscOptionsGetInt(NULL,"-c_str",&coefficient_structure,NULL);
831: DMCreateGlobalVector(da_prop,&properties);
832: DMCreateLocalVector(da_prop,&l_properties);
833: DMDAVecGetArray(da_prop,l_properties,&element_props);
835: DMGetCoordinateDM(da_prop,&prop_cda);
836: DMGetCoordinatesLocal(da_prop,&prop_coords);
837: DMDAVecGetArray(prop_cda,prop_coords,&_prop_coords);
839: DMDAGetGhostCorners(prop_cda,&si,&sj,0,&nx,&ny,0);
841: DMGetCoordinateDM(elas_da,&vel_cda);
842: DMGetCoordinatesLocal(elas_da,&vel_coords);
843: DMDAVecGetArray(vel_cda,vel_coords,&_vel_coords);
846: /* interpolate the coordinates */
847: for (j = sj; j < sj+ny; j++) {
848: for (i = si; i < si+nx; i++) {
849: PetscInt ngp;
850: PetscScalar gp_xi[GAUSS_POINTS][2],gp_weight[GAUSS_POINTS];
851: PetscScalar el_coords[8];
853: GetElementCoords(_vel_coords,i,j,el_coords);
854: ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);
856: for (p = 0; p < GAUSS_POINTS; p++) {
857: PetscScalar gp_x,gp_y;
858: PetscInt n;
859: PetscScalar xi_p[2],Ni_p[4];
861: xi_p[0] = gp_xi[p][0];
862: xi_p[1] = gp_xi[p][1];
863: ConstructQ12D_Ni(xi_p,Ni_p);
865: gp_x = 0.0;
866: gp_y = 0.0;
867: for (n = 0; n < NODES_PER_EL; n++) {
868: gp_x = gp_x+Ni_p[n]*el_coords[2*n];
869: gp_y = gp_y+Ni_p[n]*el_coords[2*n+1];
870: }
871: element_props[j][i].gp_coords[2*p] = gp_x;
872: element_props[j][i].gp_coords[2*p+1] = gp_y;
873: }
874: }
875: }
877: /* define the coefficients */
878: PetscOptionsGetBool(NULL,"-use_gp_coords",&use_gp_coords,&flg);
880: for (j = sj; j < sj+ny; j++) {
881: for (i = si; i < si+nx; i++) {
882: PetscScalar centroid_x = _prop_coords[j][i].x; /* centroids of cell */
883: PetscScalar centroid_y = _prop_coords[j][i].y;
884: PETSC_UNUSED PetscScalar coord_x,coord_y;
887: if (coefficient_structure == 0) { /* isotropic */
888: PetscScalar opts_E,opts_nu;
890: opts_E = 1.0;
891: opts_nu = 0.33;
892: PetscOptionsGetScalar(NULL,"-iso_E",&opts_E,&flg);
893: PetscOptionsGetScalar(NULL,"-iso_nu",&opts_nu,&flg);
895: for (p = 0; p < GAUSS_POINTS; p++) {
896: element_props[j][i].E[p] = opts_E;
897: element_props[j][i].nu[p] = opts_nu;
899: element_props[j][i].fx[p] = 0.0;
900: element_props[j][i].fy[p] = 0.0;
901: }
902: } else if (coefficient_structure == 1) { /* step */
903: PetscScalar opts_E0,opts_nu0,opts_xc;
904: PetscScalar opts_E1,opts_nu1;
906: opts_E0 = opts_E1 = 1.0;
907: opts_nu0 = opts_nu1 = 0.333;
908: opts_xc = 0.5;
909: PetscOptionsGetScalar(NULL,"-step_E0",&opts_E0,&flg);
910: PetscOptionsGetScalar(NULL,"-step_nu0",&opts_nu0,&flg);
911: PetscOptionsGetScalar(NULL,"-step_E1",&opts_E1,&flg);
912: PetscOptionsGetScalar(NULL,"-step_nu1",&opts_nu1,&flg);
913: PetscOptionsGetScalar(NULL,"-step_xc",&opts_xc,&flg);
915: for (p = 0; p < GAUSS_POINTS; p++) {
916: coord_x = centroid_x;
917: coord_y = centroid_y;
918: if (use_gp_coords) {
919: coord_x = element_props[j][i].gp_coords[2*p];
920: coord_y = element_props[j][i].gp_coords[2*p+1];
921: }
923: element_props[j][i].E[p] = opts_E0;
924: element_props[j][i].nu[p] = opts_nu0;
925: if (PetscRealPart(coord_x) > PetscRealPart(opts_xc)) {
926: element_props[j][i].E[p] = opts_E1;
927: element_props[j][i].nu[p] = opts_nu1;
928: }
930: element_props[j][i].fx[p] = 0.0;
931: element_props[j][i].fy[p] = 0.0;
932: }
933: } else if (coefficient_structure == 2) { /* brick */
934: PetscReal values_E[10];
935: PetscReal values_nu[10];
936: PetscInt nbricks,maxnbricks;
937: PetscInt index,span;
938: PetscInt jj;
940: flg = PETSC_FALSE;
941: maxnbricks = 10;
942: PetscOptionsGetRealArray(NULL, "-brick_E",values_E,&maxnbricks,&flg);
943: nbricks = maxnbricks;
944: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply a list of E values for each brick");
946: flg = PETSC_FALSE;
947: maxnbricks = 10;
948: PetscOptionsGetRealArray(NULL, "-brick_nu",values_nu,&maxnbricks,&flg);
949: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply a list of nu values for each brick");
950: if (maxnbricks != nbricks) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply equal numbers of values for E and nu");
952: span = 1;
953: PetscOptionsGetInt(NULL,"-brick_span",&span,&flg);
955: /* cycle through the indices so that no two material properties are repeated in lines of x or y */
956: jj = (j/span)%nbricks;
957: index = (jj+i/span)%nbricks;
958: /*printf("j=%d: index = %d \n", j,index); */
960: for (p = 0; p < GAUSS_POINTS; p++) {
961: element_props[j][i].E[p] = values_E[index];
962: element_props[j][i].nu[p] = values_nu[index];
963: }
964: } else if (coefficient_structure == 3) { /* sponge */
965: PetscScalar opts_E0,opts_nu0;
966: PetscScalar opts_E1,opts_nu1;
967: PetscInt opts_t,opts_w;
968: PetscInt ii,jj,ci,cj;
970: opts_E0 = opts_E1 = 1.0;
971: opts_nu0 = opts_nu1 = 0.333;
972: PetscOptionsGetScalar(NULL,"-sponge_E0",&opts_E0,&flg);
973: PetscOptionsGetScalar(NULL,"-sponge_nu0",&opts_nu0,&flg);
974: PetscOptionsGetScalar(NULL,"-sponge_E1",&opts_E1,&flg);
975: PetscOptionsGetScalar(NULL,"-sponge_nu1",&opts_nu1,&flg);
977: opts_t = opts_w = 1;
978: PetscOptionsGetInt(NULL,"-sponge_t",&opts_t,&flg);
979: PetscOptionsGetInt(NULL,"-sponge_w",&opts_w,&flg);
981: ii = (i)/(opts_t+opts_w+opts_t);
982: jj = (j)/(opts_t+opts_w+opts_t);
984: ci = i - ii*(opts_t+opts_w+opts_t);
985: cj = j - jj*(opts_t+opts_w+opts_t);
987: for (p = 0; p < GAUSS_POINTS; p++) {
988: element_props[j][i].E[p] = opts_E0;
989: element_props[j][i].nu[p] = opts_nu0;
990: }
991: if ((ci >= opts_t) && (ci < opts_t+opts_w)) {
992: if ((cj >= opts_t) && (cj < opts_t+opts_w)) {
993: for (p = 0; p < GAUSS_POINTS; p++) {
994: element_props[j][i].E[p] = opts_E1;
995: element_props[j][i].nu[p] = opts_nu1;
996: }
997: }
998: }
1000: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown coefficient_structure");
1001: }
1002: }
1003: DMDAVecRestoreArray(prop_cda,prop_coords,&_prop_coords);
1005: DMDAVecRestoreArray(vel_cda,vel_coords,&_vel_coords);
1007: DMDAVecRestoreArray(da_prop,l_properties,&element_props);
1008: DMLocalToGlobalBegin(da_prop,l_properties,ADD_VALUES,properties);
1009: DMLocalToGlobalEnd(da_prop,l_properties,ADD_VALUES,properties);
1011: PetscOptionsGetBool(NULL,"-no_view",&no_view,NULL);
1012: if (!no_view) {
1013: DMDAViewCoefficientsGnuplot2d(da_prop,properties,"Coeffcients for elasticity eqn.","properties");
1014: DMDACoordViewGnuplot2d(elas_da,"mesh");
1015: }
1017: /* Generate a matrix with the correct non-zero pattern of type AIJ. This will work in parallel and serial */
1018: DMCreateMatrix(elas_da,MATAIJ,&A);
1019: DMGetCoordinates(elas_da,&vel_coords);
1020: MatNullSpaceCreateRigidBody(vel_coords,&matnull);
1021: MatSetNearNullSpace(A,matnull);
1022: MatNullSpaceDestroy(&matnull);
1023: MatGetVecs(A,&f,&X);
1025: /* assemble A11 */
1026: MatZeroEntries(A);
1027: VecZeroEntries(f);
1029: AssembleA_Elasticity(A,elas_da,da_prop,properties);
1030: /* build force vector */
1031: AssembleF_Elasticity(f,elas_da,da_prop,properties);
1034: KSPCreate(PETSC_COMM_WORLD,&ksp_E);
1035: KSPSetOptionsPrefix(ksp_E,"elas_"); /* elasticity */
1037: PetscOptionsGetBool(NULL,"-use_nonsymbc",&use_nonsymbc,&flg);
1038: /* solve */
1039: if (!use_nonsymbc) {
1040: Mat AA;
1041: Vec ff,XX;
1042: IS is;
1043: VecScatter scat;
1045: DMDABCApplySymmetricCompression(elas_da,A,f,&is,&AA,&ff);
1046: VecDuplicate(ff,&XX);
1048: KSPSetOperators(ksp_E,AA,AA,SAME_NONZERO_PATTERN);
1049: KSPSetFromOptions(ksp_E);
1051: KSPSolve(ksp_E,ff,XX);
1053: /* push XX back into X */
1054: DMDABCApplyCompression(elas_da,NULL,X);
1056: VecScatterCreate(XX,NULL,X,is,&scat);
1057: VecScatterBegin(scat,XX,X,INSERT_VALUES,SCATTER_FORWARD);
1058: VecScatterEnd(scat,XX,X,INSERT_VALUES,SCATTER_FORWARD);
1059: VecScatterDestroy(&scat);
1061: MatDestroy(&AA);
1062: VecDestroy(&ff);
1063: VecDestroy(&XX);
1064: ISDestroy(&is);
1065: } else {
1066: DMDABCApplyCompression(elas_da,A,f);
1068: KSPSetOperators(ksp_E,A,A,SAME_NONZERO_PATTERN);
1069: KSPSetFromOptions(ksp_E);
1071: KSPSolve(ksp_E,f,X);
1072: }
1074: if (!no_view) {DMDAViewGnuplot2d(elas_da,X,"Displacement solution for elasticity eqn.","X");}
1075: KSPDestroy(&ksp_E);
1077: VecDestroy(&X);
1078: VecDestroy(&f);
1079: MatDestroy(&A);
1081: DMDestroy(&elas_da);
1082: DMDestroy(&da_prop);
1084: VecDestroy(&properties);
1085: VecDestroy(&l_properties);
1086: return(0);
1087: }
1091: int main(int argc,char **args)
1092: {
1094: PetscInt mx,my;
1096: PetscInitialize(&argc,&args,(char*)0,help);
1098: mx = my = 10;
1099: PetscOptionsGetInt(NULL,"-mx",&mx,NULL);
1100: PetscOptionsGetInt(NULL,"-my",&my,NULL);
1102: solve_elasticity_2d(mx,my);
1104: PetscFinalize();
1105: return 0;
1106: }
1108: /* -------------------------- helpers for boundary conditions -------------------------------- */
1112: static PetscErrorCode BCApply_EAST(DM da,PetscInt d_idx,PetscScalar bc_val,Mat A,Vec b)
1113: {
1114: DM cda;
1115: Vec coords;
1116: PetscInt si,sj,nx,ny,i,j;
1117: PetscInt M,N;
1118: DMDACoor2d **_coords;
1119: PetscInt *g_idx;
1120: PetscInt *bc_global_ids;
1121: PetscScalar *bc_vals;
1122: PetscInt nbcs;
1123: PetscInt n_dofs;
1127: /* enforce bc's */
1128: DMDAGetGlobalIndices(da,NULL,&g_idx);
1130: DMGetCoordinateDM(da,&cda);
1131: DMGetCoordinatesLocal(da,&coords);
1132: DMDAVecGetArray(cda,coords,&_coords);
1133: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1134: DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);
1136: /* --- */
1138: PetscMalloc(sizeof(PetscInt)*ny*n_dofs,&bc_global_ids);
1139: PetscMalloc(sizeof(PetscScalar)*ny*n_dofs,&bc_vals);
1141: /* init the entries to -1 so VecSetValues will ignore them */
1142: for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;
1144: i = nx-1;
1145: for (j = 0; j < ny; j++) {
1146: PetscInt local_id;
1147: PETSC_UNUSED PetscScalar coordx,coordy;
1149: local_id = i+j*nx;
1151: bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];
1153: coordx = _coords[j+sj][i+si].x;
1154: coordy = _coords[j+sj][i+si].y;
1156: bc_vals[j] = bc_val;
1157: }
1158: nbcs = 0;
1159: if ((si+nx) == (M)) nbcs = ny;
1161: if (b) {
1162: VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1163: VecAssemblyBegin(b);
1164: VecAssemblyEnd(b);
1165: }
1166: if (A) {
1167: MatZeroRows(A,nbcs,bc_global_ids,1.0,0,0);
1168: }
1170: PetscFree(bc_vals);
1171: PetscFree(bc_global_ids);
1173: DMDAVecRestoreArray(cda,coords,&_coords);
1174: return(0);
1175: }
1179: static PetscErrorCode BCApply_WEST(DM da,PetscInt d_idx,PetscScalar bc_val,Mat A,Vec b)
1180: {
1181: DM cda;
1182: Vec coords;
1183: PetscInt si,sj,nx,ny,i,j;
1184: PetscInt M,N;
1185: DMDACoor2d **_coords;
1186: PetscInt *g_idx;
1187: PetscInt *bc_global_ids;
1188: PetscScalar *bc_vals;
1189: PetscInt nbcs;
1190: PetscInt n_dofs;
1194: /* enforce bc's */
1195: DMDAGetGlobalIndices(da,NULL,&g_idx);
1197: DMGetCoordinateDM(da,&cda);
1198: DMGetCoordinatesLocal(da,&coords);
1199: DMDAVecGetArray(cda,coords,&_coords);
1200: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1201: DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);
1203: /* --- */
1205: PetscMalloc(sizeof(PetscInt)*ny*n_dofs,&bc_global_ids);
1206: PetscMalloc(sizeof(PetscScalar)*ny*n_dofs,&bc_vals);
1208: /* init the entries to -1 so VecSetValues will ignore them */
1209: for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;
1211: i = 0;
1212: for (j = 0; j < ny; j++) {
1213: PetscInt local_id;
1214: PETSC_UNUSED PetscScalar coordx,coordy;
1216: local_id = i+j*nx;
1218: bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];
1220: coordx = _coords[j+sj][i+si].x;
1221: coordy = _coords[j+sj][i+si].y;
1223: bc_vals[j] = bc_val;
1224: }
1225: nbcs = 0;
1226: if (si == 0) nbcs = ny;
1228: if (b) {
1229: VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1230: VecAssemblyBegin(b);
1231: VecAssemblyEnd(b);
1232: }
1233: if (A) {
1234: MatZeroRows(A,nbcs,bc_global_ids,1.0,0,0);
1235: }
1237: PetscFree(bc_vals);
1238: PetscFree(bc_global_ids);
1240: DMDAVecRestoreArray(cda,coords,&_coords);
1241: return(0);
1242: }
1246: static PetscErrorCode DMDABCApplyCompression(DM elas_da,Mat A,Vec f)
1247: {
1251: BCApply_EAST(elas_da,0,-1.0,A,f);
1252: BCApply_EAST(elas_da,1, 0.0,A,f);
1253: BCApply_WEST(elas_da,0,1.0,A,f);
1254: BCApply_WEST(elas_da,1,0.0,A,f);
1255: return(0);
1256: }
1260: static PetscErrorCode DMDABCApplySymmetricCompression(DM elas_da,Mat A,Vec f,IS *dofs,Mat *AA,Vec *ff)
1261: {
1263: PetscInt start,end,m;
1264: PetscInt *unconstrained;
1265: PetscInt cnt,i;
1266: Vec x;
1267: PetscScalar *_x;
1268: IS is;
1269: VecScatter scat;
1272: /* push bc's into f and A */
1273: VecDuplicate(f,&x);
1274: BCApply_EAST(elas_da,0,-1.0,A,x);
1275: BCApply_EAST(elas_da,1, 0.0,A,x);
1276: BCApply_WEST(elas_da,0,1.0,A,x);
1277: BCApply_WEST(elas_da,1,0.0,A,x);
1279: /* define which dofs are not constrained */
1280: VecGetLocalSize(x,&m);
1281: PetscMalloc(sizeof(PetscInt)*m,&unconstrained);
1282: VecGetOwnershipRange(x,&start,&end);
1283: VecGetArray(x,&_x);
1284: cnt = 0;
1285: for (i = 0; i < m; i++) {
1286: PetscReal val;
1288: val = PetscRealPart(_x[i]);
1289: if (fabs(val) < 0.1) {
1290: unconstrained[cnt] = start + i;
1291: cnt++;
1292: }
1293: }
1294: VecRestoreArray(x,&_x);
1296: ISCreateGeneral(PETSC_COMM_WORLD,cnt,unconstrained,PETSC_COPY_VALUES,&is);
1297: PetscFree(unconstrained);
1299: /* define correction for dirichlet in the rhs */
1300: MatMult(A,x,f);
1301: VecScale(f,-1.0);
1303: /* get new matrix */
1304: MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,AA);
1305: /* get new vector */
1306: MatGetVecs(*AA,NULL,ff);
1308: VecScatterCreate(f,is,*ff,NULL,&scat);
1309: VecScatterBegin(scat,f,*ff,INSERT_VALUES,SCATTER_FORWARD);
1310: VecScatterEnd(scat,f,*ff,INSERT_VALUES,SCATTER_FORWARD);
1311: VecScatterDestroy(&scat);
1313: *dofs = is;
1314: VecDestroy(&x);
1315: return(0);
1316: }