Actual source code: ex49.c
petsc-3.7.5 2017-01-01
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 : indicates the structure of the coefficients to use. \n"
9: "\
10: -c_str 0 => Setup for an 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 => Setup for a 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 => Setup for a 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 => Setup for a 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);
78: #define NSD 2 /* number of spatial dimensions */
79: #define NODES_PER_EL 4 /* nodes per element */
80: #define U_DOFS 2 /* degrees of freedom per displacement node */
81: #define GAUSS_POINTS 4
83: /* cell based evaluation */
84: typedef struct {
85: PetscScalar E,nu,fx,fy;
86: } Coefficients;
88: /* Gauss point based evaluation 8+4+4+4 = 20 */
89: typedef struct {
90: PetscScalar gp_coords[2*GAUSS_POINTS];
91: PetscScalar E[GAUSS_POINTS];
92: PetscScalar nu[GAUSS_POINTS];
93: PetscScalar fx[GAUSS_POINTS];
94: PetscScalar fy[GAUSS_POINTS];
95: } GaussPointCoefficients;
97: typedef struct {
98: PetscScalar ux_dof;
99: PetscScalar uy_dof;
100: } ElasticityDOF;
103: /*
105: D = E/((1+nu)(1-2nu)) * [ 1-nu nu 0 ]
106: [ nu 1-nu 0 ]
107: [ 0 0 0.5*(1-2nu) ]
109: B = [ d_dx 0 ]
110: [ 0 d_dy ]
111: [ d_dy d_dx ]
113: */
115: /* FEM routines */
116: /*
117: Element: Local basis function ordering
118: 1-----2
119: | |
120: | |
121: 0-----3
122: */
123: static void ConstructQ12D_Ni(PetscScalar _xi[],PetscScalar Ni[])
124: {
125: PetscScalar xi = _xi[0];
126: PetscScalar eta = _xi[1];
128: Ni[0] = 0.25*(1.0-xi)*(1.0-eta);
129: Ni[1] = 0.25*(1.0-xi)*(1.0+eta);
130: Ni[2] = 0.25*(1.0+xi)*(1.0+eta);
131: Ni[3] = 0.25*(1.0+xi)*(1.0-eta);
132: }
134: static void ConstructQ12D_GNi(PetscScalar _xi[],PetscScalar GNi[][NODES_PER_EL])
135: {
136: PetscScalar xi = _xi[0];
137: PetscScalar eta = _xi[1];
139: GNi[0][0] = -0.25*(1.0-eta);
140: GNi[0][1] = -0.25*(1.0+eta);
141: GNi[0][2] = 0.25*(1.0+eta);
142: GNi[0][3] = 0.25*(1.0-eta);
144: GNi[1][0] = -0.25*(1.0-xi);
145: GNi[1][1] = 0.25*(1.0-xi);
146: GNi[1][2] = 0.25*(1.0+xi);
147: GNi[1][3] = -0.25*(1.0+xi);
148: }
150: static void ConstructQ12D_GNx(PetscScalar GNi[][NODES_PER_EL],PetscScalar GNx[][NODES_PER_EL],PetscScalar coords[],PetscScalar *det_J)
151: {
152: PetscScalar J00,J01,J10,J11,J;
153: PetscScalar iJ00,iJ01,iJ10,iJ11;
154: PetscInt i;
156: J00 = J01 = J10 = J11 = 0.0;
157: for (i = 0; i < NODES_PER_EL; i++) {
158: PetscScalar cx = coords[2*i+0];
159: PetscScalar cy = coords[2*i+1];
161: J00 = J00+GNi[0][i]*cx; /* J_xx = dx/dxi */
162: J01 = J01+GNi[0][i]*cy; /* J_xy = dy/dxi */
163: J10 = J10+GNi[1][i]*cx; /* J_yx = dx/deta */
164: J11 = J11+GNi[1][i]*cy; /* J_yy = dy/deta */
165: }
166: J = (J00*J11)-(J01*J10);
168: iJ00 = J11/J;
169: iJ01 = -J01/J;
170: iJ10 = -J10/J;
171: iJ11 = J00/J;
174: for (i = 0; i < NODES_PER_EL; i++) {
175: GNx[0][i] = GNi[0][i]*iJ00+GNi[1][i]*iJ01;
176: GNx[1][i] = GNi[0][i]*iJ10+GNi[1][i]*iJ11;
177: }
179: if (det_J) *det_J = J;
180: }
182: static void ConstructGaussQuadrature(PetscInt *ngp,PetscScalar gp_xi[][2],PetscScalar gp_weight[])
183: {
184: *ngp = 4;
185: gp_xi[0][0] = -0.57735026919;gp_xi[0][1] = -0.57735026919;
186: gp_xi[1][0] = -0.57735026919;gp_xi[1][1] = 0.57735026919;
187: gp_xi[2][0] = 0.57735026919;gp_xi[2][1] = 0.57735026919;
188: gp_xi[3][0] = 0.57735026919;gp_xi[3][1] = -0.57735026919;
189: gp_weight[0] = 1.0;
190: gp_weight[1] = 1.0;
191: gp_weight[2] = 1.0;
192: gp_weight[3] = 1.0;
193: }
196: /* procs to the left claim the ghost node as their element */
199: static PetscErrorCode DMDAGetLocalElementSize(DM da,PetscInt *mxl,PetscInt *myl,PetscInt *mzl)
200: {
202: PetscInt m,n,p,M,N,P;
203: PetscInt sx,sy,sz;
206: DMDAGetInfo(da,0,&M,&N,&P,0,0,0,0,0,0,0,0,0);
207: DMDAGetCorners(da,&sx,&sy,&sz,&m,&n,&p);
209: if (mxl) {
210: *mxl = m;
211: if ((sx+m) == M) *mxl = m-1; /* last proc */
212: }
213: if (myl) {
214: *myl = n;
215: if ((sy+n) == N) *myl = n-1; /* last proc */
216: }
217: if (mzl) {
218: *mzl = p;
219: if ((sz+p) == P) *mzl = p-1; /* last proc */
220: }
221: return(0);
222: }
226: static PetscErrorCode DMDAGetElementCorners(DM da,PetscInt *sx,PetscInt *sy,PetscInt *sz,PetscInt *mx,PetscInt *my,PetscInt *mz)
227: {
229: PetscInt si,sj,sk;
232: DMDAGetGhostCorners(da,&si,&sj,&sk,0,0,0);
234: if (sx) {
235: *sx = si;
236: if (si != 0) *sx = si+1;
237: }
238: if (sy) {
239: *sy = sj;
240: if (sj != 0) *sy = sj+1;
241: }
243: if (sk) {
244: *sz = sk;
245: if (sk != 0) *sz = sk+1;
246: }
248: DMDAGetLocalElementSize(da,mx,my,mz);
249: return(0);
250: }
254: static PetscErrorCode DMDAGetElementOwnershipRanges2d(DM da,PetscInt **_lx,PetscInt **_ly)
255: {
257: PetscMPIInt rank;
258: PetscInt proc_I,proc_J;
259: PetscInt cpu_x,cpu_y;
260: PetscInt local_mx,local_my;
261: Vec vlx,vly;
262: PetscInt *LX,*LY,i;
263: PetscScalar *_a;
264: Vec V_SEQ;
265: VecScatter ctx;
268: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
270: DMDAGetInfo(da,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);
272: proc_J = rank/cpu_x;
273: proc_I = rank-cpu_x*proc_J;
275: PetscMalloc1(cpu_x,&LX);
276: PetscMalloc1(cpu_y,&LY);
278: DMDAGetLocalElementSize(da,&local_mx,&local_my,NULL);
279: VecCreate(PETSC_COMM_WORLD,&vlx);
280: VecSetSizes(vlx,PETSC_DECIDE,cpu_x);
281: VecSetFromOptions(vlx);
283: VecCreate(PETSC_COMM_WORLD,&vly);
284: VecSetSizes(vly,PETSC_DECIDE,cpu_y);
285: VecSetFromOptions(vly);
287: VecSetValue(vlx,proc_I,(PetscScalar)(local_mx+1.0e-9),INSERT_VALUES);
288: VecSetValue(vly,proc_J,(PetscScalar)(local_my+1.0e-9),INSERT_VALUES);
289: VecAssemblyBegin(vlx);VecAssemblyEnd(vlx);
290: VecAssemblyBegin(vly);VecAssemblyEnd(vly);
294: VecScatterCreateToAll(vlx,&ctx,&V_SEQ);
295: VecScatterBegin(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
296: VecScatterEnd(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
297: VecGetArray(V_SEQ,&_a);
298: for (i = 0; i < cpu_x; i++) LX[i] = (PetscInt)PetscRealPart(_a[i]);
299: VecRestoreArray(V_SEQ,&_a);
300: VecScatterDestroy(&ctx);
301: VecDestroy(&V_SEQ);
303: VecScatterCreateToAll(vly,&ctx,&V_SEQ);
304: VecScatterBegin(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
305: VecScatterEnd(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
306: VecGetArray(V_SEQ,&_a);
307: for (i = 0; i < cpu_y; i++) LY[i] = (PetscInt)PetscRealPart(_a[i]);
308: VecRestoreArray(V_SEQ,&_a);
309: VecScatterDestroy(&ctx);
310: VecDestroy(&V_SEQ);
312: *_lx = LX;
313: *_ly = LY;
315: VecDestroy(&vlx);
316: VecDestroy(&vly);
317: return(0);
318: }
322: static PetscErrorCode DMDACoordViewGnuplot2d(DM da,const char prefix[])
323: {
324: DM cda;
325: Vec coords;
326: DMDACoor2d **_coords;
327: PetscInt si,sj,nx,ny,i,j;
328: FILE *fp;
329: char fname[PETSC_MAX_PATH_LEN];
330: PetscMPIInt rank;
334: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
335: PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
336: PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
337: if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");
338: PetscFPrintf(PETSC_COMM_SELF,fp,"### Element geometry for processor %1.4d ### \n",rank);
340: DMGetCoordinateDM(da,&cda);
341: DMGetCoordinatesLocal(da,&coords);
342: DMDAVecGetArray(cda,coords,&_coords);
343: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
344: for (j = sj; j < sj+ny-1; j++) {
345: for (i = si; i < si+nx-1; i++) {
346: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j][i].x),PetscRealPart(_coords[j][i].y));
347: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j+1][i].x),PetscRealPart(_coords[j+1][i].y));
348: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j+1][i+1].x),PetscRealPart(_coords[j+1][i+1].y));
349: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",PetscRealPart(_coords[j][i+1].x),PetscRealPart(_coords[j][i+1].y));
350: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n\n",PetscRealPart(_coords[j][i].x),PetscRealPart(_coords[j][i].y));
351: }
352: }
353: DMDAVecRestoreArray(cda,coords,&_coords);
355: PetscFClose(PETSC_COMM_SELF,fp);
356: return(0);
357: }
361: static PetscErrorCode DMDAViewGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
362: {
363: DM cda;
364: Vec coords,local_fields;
365: DMDACoor2d **_coords;
366: FILE *fp;
367: char fname[PETSC_MAX_PATH_LEN];
368: const char *field_name;
369: PetscMPIInt rank;
370: PetscInt si,sj,nx,ny,i,j;
371: PetscInt n_dofs,d;
372: PetscScalar *_fields;
376: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
377: PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
378: PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
379: if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");
381: PetscFPrintf(PETSC_COMM_SELF,fp,"### %s (processor %1.4d) ### \n",comment,rank);
382: DMDAGetInfo(da,0,0,0,0,0,0,0,&n_dofs,0,0,0,0,0);
383: PetscFPrintf(PETSC_COMM_SELF,fp,"### x y ");
384: for (d = 0; d < n_dofs; d++) {
385: DMDAGetFieldName(da,d,&field_name);
386: PetscFPrintf(PETSC_COMM_SELF,fp,"%s ",field_name);
387: }
388: PetscFPrintf(PETSC_COMM_SELF,fp,"###\n");
391: DMGetCoordinateDM(da,&cda);
392: DMGetCoordinatesLocal(da,&coords);
393: DMDAVecGetArray(cda,coords,&_coords);
394: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
396: DMCreateLocalVector(da,&local_fields);
397: DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
398: DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
399: VecGetArray(local_fields,&_fields);
401: for (j = sj; j < sj+ny; j++) {
402: for (i = si; i < si+nx; i++) {
403: PetscScalar coord_x,coord_y;
404: PetscScalar field_d;
406: coord_x = _coords[j][i].x;
407: coord_y = _coords[j][i].y;
409: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",PetscRealPart(coord_x),PetscRealPart(coord_y));
410: for (d = 0; d < n_dofs; d++) {
411: field_d = _fields[n_dofs*((i-si)+(j-sj)*(nx))+d];
412: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e ",PetscRealPart(field_d));
413: }
414: PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
415: }
416: }
417: VecRestoreArray(local_fields,&_fields);
418: VecDestroy(&local_fields);
420: DMDAVecRestoreArray(cda,coords,&_coords);
422: PetscFClose(PETSC_COMM_SELF,fp);
423: return(0);
424: }
428: static PetscErrorCode DMDAViewCoefficientsGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
429: {
430: DM cda;
431: Vec local_fields;
432: FILE *fp;
433: char fname[PETSC_MAX_PATH_LEN];
434: const char *field_name;
435: PetscMPIInt rank;
436: PetscInt si,sj,nx,ny,i,j,p;
437: PetscInt n_dofs,d;
438: GaussPointCoefficients **_coefficients;
439: PetscErrorCode ierr;
442: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
443: PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
444: PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
445: if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");
447: PetscFPrintf(PETSC_COMM_SELF,fp,"### %s (processor %1.4d) ### \n",comment,rank);
448: DMDAGetInfo(da,0,0,0,0,0,0,0,&n_dofs,0,0,0,0,0);
449: PetscFPrintf(PETSC_COMM_SELF,fp,"### x y ");
450: for (d = 0; d < n_dofs; d++) {
451: DMDAGetFieldName(da,d,&field_name);
452: PetscFPrintf(PETSC_COMM_SELF,fp,"%s ",field_name);
453: }
454: PetscFPrintf(PETSC_COMM_SELF,fp,"###\n");
457: DMGetCoordinateDM(da,&cda);
458: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
460: DMCreateLocalVector(da,&local_fields);
461: DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
462: DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
463: DMDAVecGetArray(da,local_fields,&_coefficients);
465: for (j = sj; j < sj+ny; j++) {
466: for (i = si; i < si+nx; i++) {
467: PetscScalar coord_x,coord_y;
469: for (p = 0; p < GAUSS_POINTS; p++) {
470: coord_x = _coefficients[j][i].gp_coords[2*p];
471: coord_y = _coefficients[j][i].gp_coords[2*p+1];
473: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",PetscRealPart(coord_x),PetscRealPart(coord_y));
475: PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e %1.6e %1.6e",
476: PetscRealPart(_coefficients[j][i].E[p]),PetscRealPart(_coefficients[j][i].nu[p]),
477: PetscRealPart(_coefficients[j][i].fx[p]),PetscRealPart(_coefficients[j][i].fy[p]));
478: PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
479: }
480: }
481: }
482: DMDAVecRestoreArray(da,local_fields,&_coefficients);
483: VecDestroy(&local_fields);
485: PetscFClose(PETSC_COMM_SELF,fp);
486: return(0);
487: }
489: static void FormStressOperatorQ1(PetscScalar Ke[],PetscScalar coords[],PetscScalar E[],PetscScalar nu[])
490: {
491: PetscInt ngp;
492: PetscScalar gp_xi[GAUSS_POINTS][2];
493: PetscScalar gp_weight[GAUSS_POINTS];
494: PetscInt p,i,j,k,l;
495: PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
496: PetscScalar J_p;
497: PetscScalar B[3][U_DOFS*NODES_PER_EL];
498: PetscScalar prop_E,prop_nu,factor,constit_D[3][3];
500: /* define quadrature rule */
501: ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);
503: /* evaluate integral */
504: for (p = 0; p < ngp; p++) {
505: ConstructQ12D_GNi(gp_xi[p],GNi_p);
506: ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
508: for (i = 0; i < NODES_PER_EL; i++) {
509: PetscScalar d_dx_i = GNx_p[0][i];
510: PetscScalar d_dy_i = GNx_p[1][i];
512: B[0][2*i] = d_dx_i; B[0][2*i+1] = 0.0;
513: B[1][2*i] = 0.0; B[1][2*i+1] = d_dy_i;
514: B[2][2*i] = d_dy_i; B[2][2*i+1] = d_dx_i;
515: }
517: /* form D for the quadrature point */
518: prop_E = E[p];
519: prop_nu = nu[p];
520: factor = prop_E / ((1.0+prop_nu)*(1.0-2.0*prop_nu));
521: constit_D[0][0] = 1.0-prop_nu; constit_D[0][1] = prop_nu; constit_D[0][2] = 0.0;
522: constit_D[1][0] = prop_nu; constit_D[1][1] = 1.0-prop_nu; constit_D[1][2] = 0.0;
523: constit_D[2][0] = 0.0; constit_D[2][1] = 0.0; constit_D[2][2] = 0.5*(1.0-2.0*prop_nu);
524: for (i = 0; i < 3; i++) {
525: for (j = 0; j < 3; j++) {
526: constit_D[i][j] = factor * constit_D[i][j] * gp_weight[p] * J_p;
527: }
528: }
530: /* form Bt tildeD B */
531: /*
532: Ke_ij = Bt_ik . D_kl . B_lj
533: = B_ki . D_kl . B_lj
534: */
535: for (i = 0; i < 8; i++) {
536: for (j = 0; j < 8; j++) {
537: for (k = 0; k < 3; k++) {
538: for (l = 0; l < 3; l++) {
539: Ke[8*i+j] = Ke[8*i+j] + B[k][i] * constit_D[k][l] * B[l][j];
540: }
541: }
542: }
543: }
545: } /* end quadrature */
546: }
548: static void FormMomentumRhsQ1(PetscScalar Fe[],PetscScalar coords[],PetscScalar fx[],PetscScalar fy[])
549: {
550: PetscInt ngp;
551: PetscScalar gp_xi[GAUSS_POINTS][2];
552: PetscScalar gp_weight[GAUSS_POINTS];
553: PetscInt p,i;
554: PetscScalar Ni_p[NODES_PER_EL];
555: PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
556: PetscScalar J_p,fac;
558: /* define quadrature rule */
559: ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);
561: /* evaluate integral */
562: for (p = 0; p < ngp; p++) {
563: ConstructQ12D_Ni(gp_xi[p],Ni_p);
564: ConstructQ12D_GNi(gp_xi[p],GNi_p);
565: ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
566: fac = gp_weight[p]*J_p;
568: for (i = 0; i < NODES_PER_EL; i++) {
569: Fe[NSD*i] += fac*Ni_p[i]*fx[p];
570: Fe[NSD*i+1] += fac*Ni_p[i]*fy[p];
571: }
572: }
573: }
575: /*
576: i,j are the element indices
577: The unknown is a vector quantity.
578: The s[].c is used to indicate the degree of freedom.
579: */
582: static PetscErrorCode DMDAGetElementEqnums_u(MatStencil s_u[],PetscInt i,PetscInt j)
583: {
585: /* displacement */
586: /* node 0 */
587: s_u[0].i = i;s_u[0].j = j;s_u[0].c = 0; /* Ux0 */
588: s_u[1].i = i;s_u[1].j = j;s_u[1].c = 1; /* Uy0 */
590: /* node 1 */
591: s_u[2].i = i;s_u[2].j = j+1;s_u[2].c = 0; /* Ux1 */
592: s_u[3].i = i;s_u[3].j = j+1;s_u[3].c = 1; /* Uy1 */
594: /* node 2 */
595: s_u[4].i = i+1;s_u[4].j = j+1;s_u[4].c = 0; /* Ux2 */
596: s_u[5].i = i+1;s_u[5].j = j+1;s_u[5].c = 1; /* Uy2 */
598: /* node 3 */
599: s_u[6].i = i+1;s_u[6].j = j;s_u[6].c = 0; /* Ux3 */
600: s_u[7].i = i+1;s_u[7].j = j;s_u[7].c = 1; /* Uy3 */
601: return(0);
602: }
606: static PetscErrorCode GetElementCoords(DMDACoor2d **_coords,PetscInt ei,PetscInt ej,PetscScalar el_coords[])
607: {
609: /* get coords for the element */
610: el_coords[NSD*0+0] = _coords[ej][ei].x; el_coords[NSD*0+1] = _coords[ej][ei].y;
611: el_coords[NSD*1+0] = _coords[ej+1][ei].x; el_coords[NSD*1+1] = _coords[ej+1][ei].y;
612: el_coords[NSD*2+0] = _coords[ej+1][ei+1].x; el_coords[NSD*2+1] = _coords[ej+1][ei+1].y;
613: el_coords[NSD*3+0] = _coords[ej][ei+1].x; el_coords[NSD*3+1] = _coords[ej][ei+1].y;
614: return(0);
615: }
619: static PetscErrorCode AssembleA_Elasticity(Mat A,DM elas_da,DM properties_da,Vec properties)
620: {
621: DM cda;
622: Vec coords;
623: DMDACoor2d **_coords;
624: MatStencil u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
625: PetscInt sex,sey,mx,my;
626: PetscInt ei,ej;
627: PetscScalar Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
628: PetscScalar el_coords[NODES_PER_EL*NSD];
629: Vec local_properties;
630: GaussPointCoefficients **props;
631: PetscScalar *prop_E,*prop_nu;
632: PetscErrorCode ierr;
635: /* setup for coords */
636: DMGetCoordinateDM(elas_da,&cda);
637: DMGetCoordinatesLocal(elas_da,&coords);
638: DMDAVecGetArray(cda,coords,&_coords);
640: /* setup for coefficients */
641: DMCreateLocalVector(properties_da,&local_properties);
642: DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
643: DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
644: DMDAVecGetArray(properties_da,local_properties,&props);
646: DMDAGetElementCorners(elas_da,&sex,&sey,0,&mx,&my,0);
647: for (ej = sey; ej < sey+my; ej++) {
648: for (ei = sex; ei < sex+mx; ei++) {
649: /* get coords for the element */
650: GetElementCoords(_coords,ei,ej,el_coords);
652: /* get coefficients for the element */
653: prop_E = props[ej][ei].E;
654: prop_nu = props[ej][ei].nu;
656: /* initialise element stiffness matrix */
657: PetscMemzero(Ae,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS);
659: /* form element stiffness matrix */
660: FormStressOperatorQ1(Ae,el_coords,prop_E,prop_nu);
662: /* insert element matrix into global matrix */
663: DMDAGetElementEqnums_u(u_eqn,ei,ej);
664: MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
665: }
666: }
667: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
668: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
670: DMDAVecRestoreArray(cda,coords,&_coords);
672: DMDAVecRestoreArray(properties_da,local_properties,&props);
673: VecDestroy(&local_properties);
674: return(0);
675: }
680: static PetscErrorCode DMDASetValuesLocalStencil_ADD_VALUES(ElasticityDOF **fields_F,MatStencil u_eqn[],PetscScalar Fe_u[])
681: {
682: PetscInt n;
685: for (n = 0; n < 4; n++) {
686: 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];
687: 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];
688: }
689: return(0);
690: }
694: static PetscErrorCode AssembleF_Elasticity(Vec F,DM elas_da,DM properties_da,Vec properties)
695: {
696: DM cda;
697: Vec coords;
698: DMDACoor2d **_coords;
699: MatStencil u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
700: PetscInt sex,sey,mx,my;
701: PetscInt ei,ej;
702: PetscScalar Fe[NODES_PER_EL*U_DOFS];
703: PetscScalar el_coords[NODES_PER_EL*NSD];
704: Vec local_properties;
705: GaussPointCoefficients **props;
706: PetscScalar *prop_fx,*prop_fy;
707: Vec local_F;
708: ElasticityDOF **ff;
709: PetscErrorCode ierr;
712: /* setup for coords */
713: DMGetCoordinateDM(elas_da,&cda);
714: DMGetCoordinatesLocal(elas_da,&coords);
715: DMDAVecGetArray(cda,coords,&_coords);
717: /* setup for coefficients */
718: DMGetLocalVector(properties_da,&local_properties);
719: DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
720: DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
721: DMDAVecGetArray(properties_da,local_properties,&props);
723: /* get acces to the vector */
724: DMGetLocalVector(elas_da,&local_F);
725: VecZeroEntries(local_F);
726: DMDAVecGetArray(elas_da,local_F,&ff);
729: DMDAGetElementCorners(elas_da,&sex,&sey,0,&mx,&my,0);
730: for (ej = sey; ej < sey+my; ej++) {
731: for (ei = sex; ei < sex+mx; ei++) {
732: /* get coords for the element */
733: GetElementCoords(_coords,ei,ej,el_coords);
735: /* get coefficients for the element */
736: prop_fx = props[ej][ei].fx;
737: prop_fy = props[ej][ei].fy;
739: /* initialise element stiffness matrix */
740: PetscMemzero(Fe,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS);
742: /* form element stiffness matrix */
743: FormMomentumRhsQ1(Fe,el_coords,prop_fx,prop_fy);
745: /* insert element matrix into global matrix */
746: DMDAGetElementEqnums_u(u_eqn,ei,ej);
748: DMDASetValuesLocalStencil_ADD_VALUES(ff,u_eqn,Fe);
749: }
750: }
752: DMDAVecRestoreArray(elas_da,local_F,&ff);
753: DMLocalToGlobalBegin(elas_da,local_F,ADD_VALUES,F);
754: DMLocalToGlobalEnd(elas_da,local_F,ADD_VALUES,F);
755: DMRestoreLocalVector(elas_da,&local_F);
757: DMDAVecRestoreArray(cda,coords,&_coords);
759: DMDAVecRestoreArray(properties_da,local_properties,&props);
760: DMRestoreLocalVector(properties_da,&local_properties);
761: return(0);
762: }
766: static PetscErrorCode solve_elasticity_2d(PetscInt mx,PetscInt my)
767: {
768: DM elas_da,da_prop;
769: PetscInt u_dof,dof,stencil_width;
770: Mat A;
771: PetscInt mxl,myl;
772: DM prop_cda,vel_cda;
773: Vec prop_coords,vel_coords;
774: PetscInt si,sj,nx,ny,i,j,p;
775: Vec f,X;
776: PetscInt prop_dof,prop_stencil_width;
777: Vec properties,l_properties;
778: MatNullSpace matnull;
779: PetscReal dx,dy;
780: PetscInt M,N;
781: DMDACoor2d **_prop_coords,**_vel_coords;
782: GaussPointCoefficients **element_props;
783: KSP ksp_E;
784: PetscInt coefficient_structure = 0;
785: PetscInt cpu_x,cpu_y,*lx = NULL,*ly = NULL;
786: PetscBool use_gp_coords = PETSC_FALSE;
787: PetscBool use_nonsymbc = PETSC_FALSE;
788: PetscBool no_view = PETSC_FALSE;
789: PetscBool flg;
790: PetscErrorCode ierr;
793: /* Generate the da for velocity and pressure */
794: /*
795: We use Q1 elements for the temperature.
796: FEM has a 9-point stencil (BOX) or connectivity pattern
797: Num nodes in each direction is mx+1, my+1
798: */
799: u_dof = U_DOFS; /* Vx, Vy - velocities */
800: dof = u_dof;
801: stencil_width = 1;
802: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,
803: mx+1,my+1,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&elas_da);
804: DMDASetFieldName(elas_da,0,"Ux");
805: DMDASetFieldName(elas_da,1,"Uy");
807: /* unit box [0,1] x [0,1] */
808: DMDASetUniformCoordinates(elas_da,0.0,1.0,0.0,1.0,0.0,1.0);
811: /* Generate element properties, we will assume all material properties are constant over the element */
812: /* local number of elements */
813: DMDAGetLocalElementSize(elas_da,&mxl,&myl,NULL);
815: /* !!! IN PARALLEL WE MUST MAKE SURE THE TWO DMDA's ALIGN !!! */
816: DMDAGetInfo(elas_da,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);
817: DMDAGetElementOwnershipRanges2d(elas_da,&lx,&ly);
819: prop_dof = (PetscInt)(sizeof(GaussPointCoefficients)/sizeof(PetscScalar)); /* gauss point setup */
820: prop_stencil_width = 0;
821: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,
822: mx,my,cpu_x,cpu_y,prop_dof,prop_stencil_width,lx,ly,&da_prop);
823: PetscFree(lx);
824: PetscFree(ly);
826: /* define centroid positions */
827: DMDAGetInfo(da_prop,0,&M,&N,0,0,0,0,0,0,0,0,0,0);
828: dx = 1.0/((PetscReal)(M));
829: dy = 1.0/((PetscReal)(N));
831: 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);
833: /* define coefficients */
834: PetscOptionsGetInt(NULL,NULL,"-c_str",&coefficient_structure,NULL);
836: DMCreateGlobalVector(da_prop,&properties);
837: DMCreateLocalVector(da_prop,&l_properties);
838: DMDAVecGetArray(da_prop,l_properties,&element_props);
840: DMGetCoordinateDM(da_prop,&prop_cda);
841: DMGetCoordinatesLocal(da_prop,&prop_coords);
842: DMDAVecGetArray(prop_cda,prop_coords,&_prop_coords);
844: DMDAGetGhostCorners(prop_cda,&si,&sj,0,&nx,&ny,0);
846: DMGetCoordinateDM(elas_da,&vel_cda);
847: DMGetCoordinatesLocal(elas_da,&vel_coords);
848: DMDAVecGetArray(vel_cda,vel_coords,&_vel_coords);
851: /* interpolate the coordinates */
852: for (j = sj; j < sj+ny; j++) {
853: for (i = si; i < si+nx; i++) {
854: PetscInt ngp;
855: PetscScalar gp_xi[GAUSS_POINTS][2],gp_weight[GAUSS_POINTS];
856: PetscScalar el_coords[8];
858: GetElementCoords(_vel_coords,i,j,el_coords);
859: ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);
861: for (p = 0; p < GAUSS_POINTS; p++) {
862: PetscScalar gp_x,gp_y;
863: PetscInt n;
864: PetscScalar xi_p[2],Ni_p[4];
866: xi_p[0] = gp_xi[p][0];
867: xi_p[1] = gp_xi[p][1];
868: ConstructQ12D_Ni(xi_p,Ni_p);
870: gp_x = 0.0;
871: gp_y = 0.0;
872: for (n = 0; n < NODES_PER_EL; n++) {
873: gp_x = gp_x+Ni_p[n]*el_coords[2*n];
874: gp_y = gp_y+Ni_p[n]*el_coords[2*n+1];
875: }
876: element_props[j][i].gp_coords[2*p] = gp_x;
877: element_props[j][i].gp_coords[2*p+1] = gp_y;
878: }
879: }
880: }
882: /* define the coefficients */
883: PetscOptionsGetBool(NULL,NULL,"-use_gp_coords",&use_gp_coords,&flg);
885: for (j = sj; j < sj+ny; j++) {
886: for (i = si; i < si+nx; i++) {
887: PetscScalar centroid_x = _prop_coords[j][i].x; /* centroids of cell */
888: PetscScalar centroid_y = _prop_coords[j][i].y;
889: PETSC_UNUSED PetscScalar coord_x,coord_y;
892: if (coefficient_structure == 0) { /* isotropic */
893: PetscScalar opts_E,opts_nu;
895: opts_E = 1.0;
896: opts_nu = 0.33;
897: PetscOptionsGetScalar(NULL,NULL,"-iso_E",&opts_E,&flg);
898: PetscOptionsGetScalar(NULL,NULL,"-iso_nu",&opts_nu,&flg);
900: for (p = 0; p < GAUSS_POINTS; p++) {
901: element_props[j][i].E[p] = opts_E;
902: element_props[j][i].nu[p] = opts_nu;
904: element_props[j][i].fx[p] = 0.0;
905: element_props[j][i].fy[p] = 0.0;
906: }
907: } else if (coefficient_structure == 1) { /* step */
908: PetscScalar opts_E0,opts_nu0,opts_xc;
909: PetscScalar opts_E1,opts_nu1;
911: opts_E0 = opts_E1 = 1.0;
912: opts_nu0 = opts_nu1 = 0.333;
913: opts_xc = 0.5;
914: PetscOptionsGetScalar(NULL,NULL,"-step_E0",&opts_E0,&flg);
915: PetscOptionsGetScalar(NULL,NULL,"-step_nu0",&opts_nu0,&flg);
916: PetscOptionsGetScalar(NULL,NULL,"-step_E1",&opts_E1,&flg);
917: PetscOptionsGetScalar(NULL,NULL,"-step_nu1",&opts_nu1,&flg);
918: PetscOptionsGetScalar(NULL,NULL,"-step_xc",&opts_xc,&flg);
920: for (p = 0; p < GAUSS_POINTS; p++) {
921: coord_x = centroid_x;
922: coord_y = centroid_y;
923: if (use_gp_coords) {
924: coord_x = element_props[j][i].gp_coords[2*p];
925: coord_y = element_props[j][i].gp_coords[2*p+1];
926: }
928: element_props[j][i].E[p] = opts_E0;
929: element_props[j][i].nu[p] = opts_nu0;
930: if (PetscRealPart(coord_x) > PetscRealPart(opts_xc)) {
931: element_props[j][i].E[p] = opts_E1;
932: element_props[j][i].nu[p] = opts_nu1;
933: }
935: element_props[j][i].fx[p] = 0.0;
936: element_props[j][i].fy[p] = 0.0;
937: }
938: } else if (coefficient_structure == 2) { /* brick */
939: PetscReal values_E[10];
940: PetscReal values_nu[10];
941: PetscInt nbricks,maxnbricks;
942: PetscInt index,span;
943: PetscInt jj;
945: flg = PETSC_FALSE;
946: maxnbricks = 10;
947: PetscOptionsGetRealArray(NULL,NULL, "-brick_E",values_E,&maxnbricks,&flg);
948: nbricks = maxnbricks;
949: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply a list of E values for each brick");
951: flg = PETSC_FALSE;
952: maxnbricks = 10;
953: PetscOptionsGetRealArray(NULL,NULL, "-brick_nu",values_nu,&maxnbricks,&flg);
954: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply a list of nu values for each brick");
955: if (maxnbricks != nbricks) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"User must supply equal numbers of values for E and nu");
957: span = 1;
958: PetscOptionsGetInt(NULL,NULL,"-brick_span",&span,&flg);
960: /* cycle through the indices so that no two material properties are repeated in lines of x or y */
961: jj = (j/span)%nbricks;
962: index = (jj+i/span)%nbricks;
963: /*printf("j=%d: index = %d \n", j,index); */
965: for (p = 0; p < GAUSS_POINTS; p++) {
966: element_props[j][i].E[p] = values_E[index];
967: element_props[j][i].nu[p] = values_nu[index];
968: }
969: } else if (coefficient_structure == 3) { /* sponge */
970: PetscScalar opts_E0,opts_nu0;
971: PetscScalar opts_E1,opts_nu1;
972: PetscInt opts_t,opts_w;
973: PetscInt ii,jj,ci,cj;
975: opts_E0 = opts_E1 = 1.0;
976: opts_nu0 = opts_nu1 = 0.333;
977: PetscOptionsGetScalar(NULL,NULL,"-sponge_E0",&opts_E0,&flg);
978: PetscOptionsGetScalar(NULL,NULL,"-sponge_nu0",&opts_nu0,&flg);
979: PetscOptionsGetScalar(NULL,NULL,"-sponge_E1",&opts_E1,&flg);
980: PetscOptionsGetScalar(NULL,NULL,"-sponge_nu1",&opts_nu1,&flg);
982: opts_t = opts_w = 1;
983: PetscOptionsGetInt(NULL,NULL,"-sponge_t",&opts_t,&flg);
984: PetscOptionsGetInt(NULL,NULL,"-sponge_w",&opts_w,&flg);
986: ii = (i)/(opts_t+opts_w+opts_t);
987: jj = (j)/(opts_t+opts_w+opts_t);
989: ci = i - ii*(opts_t+opts_w+opts_t);
990: cj = j - jj*(opts_t+opts_w+opts_t);
992: for (p = 0; p < GAUSS_POINTS; p++) {
993: element_props[j][i].E[p] = opts_E0;
994: element_props[j][i].nu[p] = opts_nu0;
995: }
996: if ((ci >= opts_t) && (ci < opts_t+opts_w)) {
997: if ((cj >= opts_t) && (cj < opts_t+opts_w)) {
998: for (p = 0; p < GAUSS_POINTS; p++) {
999: element_props[j][i].E[p] = opts_E1;
1000: element_props[j][i].nu[p] = opts_nu1;
1001: }
1002: }
1003: }
1005: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown coefficient_structure");
1006: }
1007: }
1008: DMDAVecRestoreArray(prop_cda,prop_coords,&_prop_coords);
1010: DMDAVecRestoreArray(vel_cda,vel_coords,&_vel_coords);
1012: DMDAVecRestoreArray(da_prop,l_properties,&element_props);
1013: DMLocalToGlobalBegin(da_prop,l_properties,ADD_VALUES,properties);
1014: DMLocalToGlobalEnd(da_prop,l_properties,ADD_VALUES,properties);
1016: PetscOptionsGetBool(NULL,NULL,"-no_view",&no_view,NULL);
1017: if (!no_view) {
1018: DMDAViewCoefficientsGnuplot2d(da_prop,properties,"Coeffcients for elasticity eqn.","properties");
1019: DMDACoordViewGnuplot2d(elas_da,"mesh");
1020: }
1022: /* Generate a matrix with the correct non-zero pattern of type AIJ. This will work in parallel and serial */
1023: DMSetMatType(elas_da,MATAIJ);
1024: DMCreateMatrix(elas_da,&A);
1025: MatSetBlockSize(A,2);
1026: DMGetCoordinates(elas_da,&vel_coords);
1027: MatNullSpaceCreateRigidBody(vel_coords,&matnull);
1028: MatSetNearNullSpace(A,matnull);
1029: MatNullSpaceDestroy(&matnull);
1030: MatCreateVecs(A,&f,&X);
1032: /* assemble A11 */
1033: MatZeroEntries(A);
1034: VecZeroEntries(f);
1036: AssembleA_Elasticity(A,elas_da,da_prop,properties);
1037: /* build force vector */
1038: AssembleF_Elasticity(f,elas_da,da_prop,properties);
1041: KSPCreate(PETSC_COMM_WORLD,&ksp_E);
1042: KSPSetOptionsPrefix(ksp_E,"elas_"); /* elasticity */
1044: PetscOptionsGetBool(NULL,NULL,"-use_nonsymbc",&use_nonsymbc,&flg);
1045: /* solve */
1046: if (!use_nonsymbc) {
1047: Mat AA;
1048: Vec ff,XX;
1049: IS is;
1050: VecScatter scat;
1052: DMDABCApplySymmetricCompression(elas_da,A,f,&is,&AA,&ff);
1053: VecDuplicate(ff,&XX);
1055: KSPSetOperators(ksp_E,AA,AA);
1056: KSPSetFromOptions(ksp_E);
1058: KSPSolve(ksp_E,ff,XX);
1060: /* push XX back into X */
1061: DMDABCApplyCompression(elas_da,NULL,X);
1063: VecScatterCreate(XX,NULL,X,is,&scat);
1064: VecScatterBegin(scat,XX,X,INSERT_VALUES,SCATTER_FORWARD);
1065: VecScatterEnd(scat,XX,X,INSERT_VALUES,SCATTER_FORWARD);
1066: VecScatterDestroy(&scat);
1068: MatDestroy(&AA);
1069: VecDestroy(&ff);
1070: VecDestroy(&XX);
1071: ISDestroy(&is);
1072: } else {
1073: DMDABCApplyCompression(elas_da,A,f);
1075: KSPSetOperators(ksp_E,A,A);
1076: KSPSetFromOptions(ksp_E);
1078: KSPSolve(ksp_E,f,X);
1079: }
1081: if (!no_view) {DMDAViewGnuplot2d(elas_da,X,"Displacement solution for elasticity eqn.","X");}
1082: KSPDestroy(&ksp_E);
1084: VecDestroy(&X);
1085: VecDestroy(&f);
1086: MatDestroy(&A);
1088: DMDestroy(&elas_da);
1089: DMDestroy(&da_prop);
1091: VecDestroy(&properties);
1092: VecDestroy(&l_properties);
1093: return(0);
1094: }
1098: int main(int argc,char **args)
1099: {
1101: PetscInt mx,my;
1103: PetscInitialize(&argc,&args,(char*)0,help);
1105: mx = my = 10;
1106: PetscOptionsGetInt(NULL,NULL,"-mx",&mx,NULL);
1107: PetscOptionsGetInt(NULL,NULL,"-my",&my,NULL);
1109: solve_elasticity_2d(mx,my);
1111: PetscFinalize();
1112: return 0;
1113: }
1115: /* -------------------------- helpers for boundary conditions -------------------------------- */
1119: static PetscErrorCode BCApply_EAST(DM da,PetscInt d_idx,PetscScalar bc_val,Mat A,Vec b)
1120: {
1121: DM cda;
1122: Vec coords;
1123: PetscInt si,sj,nx,ny,i,j;
1124: PetscInt M,N;
1125: DMDACoor2d **_coords;
1126: const PetscInt *g_idx;
1127: PetscInt *bc_global_ids;
1128: PetscScalar *bc_vals;
1129: PetscInt nbcs;
1130: PetscInt n_dofs;
1131: PetscErrorCode ierr;
1132: ISLocalToGlobalMapping ltogm;
1135: /* enforce bc's */
1136: DMGetLocalToGlobalMapping(da,<ogm);
1137: ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);
1139: DMGetCoordinateDM(da,&cda);
1140: DMGetCoordinatesLocal(da,&coords);
1141: DMDAVecGetArray(cda,coords,&_coords);
1142: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1143: DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);
1145: /* --- */
1147: PetscMalloc1(ny*n_dofs,&bc_global_ids);
1148: PetscMalloc1(ny*n_dofs,&bc_vals);
1150: /* init the entries to -1 so VecSetValues will ignore them */
1151: for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;
1153: i = nx-1;
1154: for (j = 0; j < ny; j++) {
1155: PetscInt local_id;
1156: PETSC_UNUSED PetscScalar coordx,coordy;
1158: local_id = i+j*nx;
1160: bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];
1162: coordx = _coords[j+sj][i+si].x;
1163: coordy = _coords[j+sj][i+si].y;
1165: bc_vals[j] = bc_val;
1166: }
1167: ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1168: nbcs = 0;
1169: if ((si+nx) == (M)) nbcs = ny;
1171: if (b) {
1172: VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1173: VecAssemblyBegin(b);
1174: VecAssemblyEnd(b);
1175: }
1176: if (A) {
1177: MatZeroRows(A,nbcs,bc_global_ids,1.0,0,0);
1178: }
1180: PetscFree(bc_vals);
1181: PetscFree(bc_global_ids);
1183: DMDAVecRestoreArray(cda,coords,&_coords);
1184: return(0);
1185: }
1189: static PetscErrorCode BCApply_WEST(DM da,PetscInt d_idx,PetscScalar bc_val,Mat A,Vec b)
1190: {
1191: DM cda;
1192: Vec coords;
1193: PetscInt si,sj,nx,ny,i,j;
1194: PetscInt M,N;
1195: DMDACoor2d **_coords;
1196: const PetscInt *g_idx;
1197: PetscInt *bc_global_ids;
1198: PetscScalar *bc_vals;
1199: PetscInt nbcs;
1200: PetscInt n_dofs;
1201: PetscErrorCode ierr;
1202: ISLocalToGlobalMapping ltogm;
1205: /* enforce bc's */
1206: DMGetLocalToGlobalMapping(da,<ogm);
1207: ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);
1209: DMGetCoordinateDM(da,&cda);
1210: DMGetCoordinatesLocal(da,&coords);
1211: DMDAVecGetArray(cda,coords,&_coords);
1212: DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1213: DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);
1215: /* --- */
1217: PetscMalloc1(ny*n_dofs,&bc_global_ids);
1218: PetscMalloc1(ny*n_dofs,&bc_vals);
1220: /* init the entries to -1 so VecSetValues will ignore them */
1221: for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;
1223: i = 0;
1224: for (j = 0; j < ny; j++) {
1225: PetscInt local_id;
1226: PETSC_UNUSED PetscScalar coordx,coordy;
1228: local_id = i+j*nx;
1230: bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];
1232: coordx = _coords[j+sj][i+si].x;
1233: coordy = _coords[j+sj][i+si].y;
1235: bc_vals[j] = bc_val;
1236: }
1237: ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1238: nbcs = 0;
1239: if (si == 0) nbcs = ny;
1241: if (b) {
1242: VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1243: VecAssemblyBegin(b);
1244: VecAssemblyEnd(b);
1245: }
1246: if (A) {
1247: MatZeroRows(A,nbcs,bc_global_ids,1.0,0,0);
1248: }
1250: PetscFree(bc_vals);
1251: PetscFree(bc_global_ids);
1253: DMDAVecRestoreArray(cda,coords,&_coords);
1254: return(0);
1255: }
1259: static PetscErrorCode DMDABCApplyCompression(DM elas_da,Mat A,Vec f)
1260: {
1264: BCApply_EAST(elas_da,0,-1.0,A,f);
1265: BCApply_EAST(elas_da,1, 0.0,A,f);
1266: BCApply_WEST(elas_da,0,1.0,A,f);
1267: BCApply_WEST(elas_da,1,0.0,A,f);
1268: return(0);
1269: }
1273: static PetscErrorCode Orthogonalize(PetscInt n,Vec *vecs)
1274: {
1275: PetscInt i,j;
1276: PetscScalar dot;
1280: for (i=0; i<n; i++) {
1281: VecNormalize(vecs[i],PETSC_NULL);
1282: for (j=i+1; j<n; j++) {
1283: VecDot(vecs[i],vecs[j],&dot);
1284: VecAXPY(vecs[j],-dot,vecs[i]);
1285: }
1286: }
1287: return(0);
1288: }
1292: static PetscErrorCode DMDABCApplySymmetricCompression(DM elas_da,Mat A,Vec f,IS *dofs,Mat *AA,Vec *ff)
1293: {
1295: PetscInt start,end,m;
1296: PetscInt *unconstrained;
1297: PetscInt cnt,i;
1298: Vec x;
1299: PetscScalar *_x;
1300: IS is;
1301: VecScatter scat;
1304: /* push bc's into f and A */
1305: VecDuplicate(f,&x);
1306: BCApply_EAST(elas_da,0,-1.0,A,x);
1307: BCApply_EAST(elas_da,1, 0.0,A,x);
1308: BCApply_WEST(elas_da,0,1.0,A,x);
1309: BCApply_WEST(elas_da,1,0.0,A,x);
1311: /* define which dofs are not constrained */
1312: VecGetLocalSize(x,&m);
1313: PetscMalloc1(m,&unconstrained);
1314: VecGetOwnershipRange(x,&start,&end);
1315: VecGetArray(x,&_x);
1316: cnt = 0;
1317: for (i = 0; i < m; i+=2) {
1318: PetscReal val1,val2;
1320: val1 = PetscRealPart(_x[i]);
1321: val2 = PetscRealPart(_x[i+1]);
1322: if (PetscAbs(val1) < 0.1 && PetscAbs(val2) < 0.1) {
1323: unconstrained[cnt] = start + i;
1324: cnt++;
1325: unconstrained[cnt] = start + i + 1;
1326: cnt++;
1327: }
1328: }
1329: VecRestoreArray(x,&_x);
1331: ISCreateGeneral(PETSC_COMM_WORLD,cnt,unconstrained,PETSC_COPY_VALUES,&is);
1332: PetscFree(unconstrained);
1333: ISSetBlockSize(is,2);
1335: /* define correction for dirichlet in the rhs */
1336: MatMult(A,x,f);
1337: VecScale(f,-1.0);
1339: /* get new matrix */
1340: MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,AA);
1341: /* get new vector */
1342: MatCreateVecs(*AA,NULL,ff);
1344: VecScatterCreate(f,is,*ff,NULL,&scat);
1345: VecScatterBegin(scat,f,*ff,INSERT_VALUES,SCATTER_FORWARD);
1346: VecScatterEnd(scat,f,*ff,INSERT_VALUES,SCATTER_FORWARD);
1348: { /* Constrain near-null space */
1349: PetscInt nvecs;
1350: const Vec *vecs;
1351: Vec *uvecs;
1352: PetscBool has_const;
1353: MatNullSpace mnull,unull;
1355: MatGetNearNullSpace(A,&mnull);
1356: MatNullSpaceGetVecs(mnull,&has_const,&nvecs,&vecs);
1357: VecDuplicateVecs(*ff,nvecs,&uvecs);
1358: for (i=0; i<nvecs; i++) {
1359: VecScatterBegin(scat,vecs[i],uvecs[i],INSERT_VALUES,SCATTER_FORWARD);
1360: VecScatterEnd(scat,vecs[i],uvecs[i],INSERT_VALUES,SCATTER_FORWARD);
1361: }
1362: Orthogonalize(nvecs,uvecs);
1363: MatNullSpaceCreate(PetscObjectComm((PetscObject)A),PETSC_FALSE,nvecs,uvecs,&unull);
1364: MatSetNearNullSpace(*AA,unull);
1365: MatNullSpaceDestroy(&unull);
1366: VecDestroyVecs(nvecs,&uvecs);
1367: }
1369: VecScatterDestroy(&scat);
1371: *dofs = is;
1372: VecDestroy(&x);
1373: return(0);
1374: }