Actual source code: ex10.c
petsc-3.7.1 2016-05-15
2: static char help[] = "Linear elastiticty with dimensions using 20 node serendipity elements.\n\
3: This also demonstrates use of block\n\
4: diagonal data structure. Input arguments are:\n\
5: -m : problem size\n\n";
7: #include <petscksp.h>
9: /* This code is not intended as an efficient implementation, it is only
10: here to produce an interesting sparse matrix quickly.
12: PLEASE DO NOT BASE ANY OF YOUR CODES ON CODE LIKE THIS, THERE ARE MUCH
13: BETTER WAYS TO DO THIS. */
15: extern PetscErrorCode GetElasticityMatrix(PetscInt,Mat*);
16: extern PetscErrorCode Elastic20Stiff(PetscReal**);
17: extern PetscErrorCode AddElement(Mat,PetscInt,PetscInt,PetscReal**,PetscInt,PetscInt);
18: extern PetscErrorCode paulsetup20(void);
19: extern PetscErrorCode paulintegrate20(PetscReal K[60][60]);
23: int main(int argc,char **args)
24: {
25: Mat mat;
27: PetscInt i,its,m = 3,rdim,cdim,rstart,rend;
28: PetscMPIInt rank,size;
29: PetscScalar v,neg1 = -1.0;
30: Vec u,x,b;
31: KSP ksp;
32: PetscReal norm;
34: PetscInitialize(&argc,&args,(char*)0,help);
35: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
36: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
37: MPI_Comm_size(PETSC_COMM_WORLD,&size);
39: /* Form matrix */
40: GetElasticityMatrix(m,&mat);
42: /* Generate vectors */
43: MatGetSize(mat,&rdim,&cdim);
44: MatGetOwnershipRange(mat,&rstart,&rend);
45: VecCreate(PETSC_COMM_WORLD,&u);
46: VecSetSizes(u,PETSC_DECIDE,rdim);
47: VecSetFromOptions(u);
48: VecDuplicate(u,&b);
49: VecDuplicate(b,&x);
50: for (i=rstart; i<rend; i++) {
51: v = (PetscScalar)(i-rstart + 100*rank);
52: VecSetValues(u,1,&i,&v,INSERT_VALUES);
53: }
54: VecAssemblyBegin(u);
55: VecAssemblyEnd(u);
57: /* Compute right-hand-side */
58: MatMult(mat,u,b);
60: /* Solve linear system */
61: KSPCreate(PETSC_COMM_WORLD,&ksp);
62: KSPSetOperators(ksp,mat,mat);
63: KSPSetFromOptions(ksp);
64: KSPSolve(ksp,b,x);
65: KSPGetIterationNumber(ksp,&its);
66: /* Check error */
67: VecAXPY(x,neg1,u);
68: VecNorm(x,NORM_2,&norm);
70: PetscPrintf(PETSC_COMM_WORLD,"Norm of residual %g Number of iterations %D\n",(double)norm,its);
72: /* Free work space */
73: KSPDestroy(&ksp);
74: VecDestroy(&u);
75: VecDestroy(&x);
76: VecDestroy(&b);
77: MatDestroy(&mat);
79: PetscFinalize();
80: return 0;
81: }
82: /* -------------------------------------------------------------------- */
85: /*
86: GetElasticityMatrix - Forms 3D linear elasticity matrix.
87: */
88: PetscErrorCode GetElasticityMatrix(PetscInt m,Mat *newmat)
89: {
90: PetscInt i,j,k,i1,i2,j_1,j2,k1,k2,h1,h2,shiftx,shifty,shiftz;
91: PetscInt ict,nz,base,r1,r2,N,*rowkeep,nstart;
93: IS iskeep;
94: PetscReal **K,norm;
95: Mat mat,submat = 0,*submatb;
96: MatType type = MATSEQBAIJ;
98: m /= 2; /* This is done just to be consistent with the old example */
99: N = 3*(2*m+1)*(2*m+1)*(2*m+1);
100: PetscPrintf(PETSC_COMM_SELF,"m = %D, N=%D\n",m,N);
101: MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,80,NULL,&mat);
103: /* Form stiffness for element */
104: PetscMalloc1(81,&K);
105: for (i=0; i<81; i++) {
106: PetscMalloc1(81,&K[i]);
107: }
108: Elastic20Stiff(K);
110: /* Loop over elements and add contribution to stiffness */
111: shiftx = 3; shifty = 3*(2*m+1); shiftz = 3*(2*m+1)*(2*m+1);
112: for (k=0; k<m; k++) {
113: for (j=0; j<m; j++) {
114: for (i=0; i<m; i++) {
115: h1 = 0;
116: base = 2*k*shiftz + 2*j*shifty + 2*i*shiftx;
117: for (k1=0; k1<3; k1++) {
118: for (j_1=0; j_1<3; j_1++) {
119: for (i1=0; i1<3; i1++) {
120: h2 = 0;
121: r1 = base + i1*shiftx + j_1*shifty + k1*shiftz;
122: for (k2=0; k2<3; k2++) {
123: for (j2=0; j2<3; j2++) {
124: for (i2=0; i2<3; i2++) {
125: r2 = base + i2*shiftx + j2*shifty + k2*shiftz;
126: AddElement(mat,r1,r2,K,h1,h2);
127: h2 += 3;
128: }
129: }
130: }
131: h1 += 3;
132: }
133: }
134: }
135: }
136: }
137: }
139: for (i=0; i<81; i++) {
140: PetscFree(K[i]);
141: }
142: PetscFree(K);
144: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
145: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
147: /* Exclude any superfluous rows and columns */
148: nstart = 3*(2*m+1)*(2*m+1);
149: ict = 0;
150: PetscMalloc1(N-nstart,&rowkeep);
151: for (i=nstart; i<N; i++) {
152: MatGetRow(mat,i,&nz,0,0);
153: if (nz) rowkeep[ict++] = i;
154: MatRestoreRow(mat,i,&nz,0,0);
155: }
156: ISCreateGeneral(PETSC_COMM_SELF,ict,rowkeep,PETSC_COPY_VALUES,&iskeep);
157: MatGetSubMatrices(mat,1,&iskeep,&iskeep,MAT_INITIAL_MATRIX,&submatb);
158: submat = *submatb;
159: PetscFree(submatb);
160: PetscFree(rowkeep);
161: ISDestroy(&iskeep);
162: MatDestroy(&mat);
164: /* Convert storage formats -- just to demonstrate conversion to various
165: formats (in particular, block diagonal storage). This is NOT the
166: recommended means to solve such a problem. */
167: MatConvert(submat,type,MAT_INITIAL_MATRIX,newmat);
168: MatDestroy(&submat);
170: MatNorm(*newmat,NORM_1,&norm);
171: PetscPrintf(PETSC_COMM_WORLD,"matrix 1 norm = %g\n",(double)norm);
173: return 0;
174: }
175: /* -------------------------------------------------------------------- */
178: PetscErrorCode AddElement(Mat mat,PetscInt r1,PetscInt r2,PetscReal **K,PetscInt h1,PetscInt h2)
179: {
180: PetscScalar val;
181: PetscInt l1,l2,row,col;
184: for (l1=0; l1<3; l1++) {
185: for (l2=0; l2<3; l2++) {
186: /*
187: NOTE you should never do this! Inserting values 1 at a time is
188: just too expensive!
189: */
190: if (K[h1+l1][h2+l2] != 0.0) {
191: row = r1+l1; col = r2+l2; val = K[h1+l1][h2+l2];
192: MatSetValues(mat,1,&row,1,&col,&val,ADD_VALUES);
193: row = r2+l2; col = r1+l1;
194: MatSetValues(mat,1,&row,1,&col,&val,ADD_VALUES);
195: }
196: }
197: }
198: return 0;
199: }
200: /* -------------------------------------------------------------------- */
201: PetscReal N[20][64]; /* Interpolation function. */
202: PetscReal part_N[3][20][64]; /* Partials of interpolation function. */
203: PetscReal rst[3][64]; /* Location of integration pts in (r,s,t) */
204: PetscReal weight[64]; /* Gaussian quadrature weights. */
205: PetscReal xyz[20][3]; /* (x,y,z) coordinates of nodes */
206: PetscReal E,nu; /* Physcial constants. */
207: PetscInt n_int,N_int; /* N_int = n_int^3, number of int. pts. */
208: /* Ordering of the vertices, (r,s,t) coordinates, of the canonical cell. */
209: PetscReal r2[20] = {-1.0,0.0,1.0,-1.0,1.0,-1.0,0.0,1.0,
210: -1.0,1.0,-1.0,1.0,
211: -1.0,0.0,1.0,-1.0,1.0,-1.0,0.0,1.0};
212: PetscReal s2[20] = {-1.0,-1.0, -1.0,0.0,0.0,1.0, 1.0, 1.0,
213: -1.0,-1.0,1.0,1.0,
214: -1.0,-1.0, -1.0,0.0,0.0,1.0, 1.0, 1.0};
215: PetscReal t2[20] = {-1.0,-1.0,-1.0,-1.0,-1.0,-1.0,-1.0,-1.0,
216: 0.0,0.0,0.0,0.0,
217: 1.0,1.0,1.0,1.0,1.0,1.0,1.0,1.0};
218: PetscInt rmap[20] = {0,1,2,3,5,6,7,8,9,11,15,17,18,19,20,21,23,24,25,26};
219: /* -------------------------------------------------------------------- */
222: /*
223: Elastic20Stiff - Forms 20 node elastic stiffness for element.
224: */
225: PetscErrorCode Elastic20Stiff(PetscReal **Ke)
226: {
227: PetscReal K[60][60],x,y,z,dx,dy,dz,v;
228: PetscInt i,j,k,l,Ii,J;
230: paulsetup20();
232: x = -1.0; y = -1.0; z = -1.0; dx = 2.0; dy = 2.0; dz = 2.0;
233: xyz[0][0] = x; xyz[0][1] = y; xyz[0][2] = z;
234: xyz[1][0] = x + dx; xyz[1][1] = y; xyz[1][2] = z;
235: xyz[2][0] = x + 2.*dx; xyz[2][1] = y; xyz[2][2] = z;
236: xyz[3][0] = x; xyz[3][1] = y + dy; xyz[3][2] = z;
237: xyz[4][0] = x + 2.*dx; xyz[4][1] = y + dy; xyz[4][2] = z;
238: xyz[5][0] = x; xyz[5][1] = y + 2.*dy; xyz[5][2] = z;
239: xyz[6][0] = x + dx; xyz[6][1] = y + 2.*dy; xyz[6][2] = z;
240: xyz[7][0] = x + 2.*dx; xyz[7][1] = y + 2.*dy; xyz[7][2] = z;
241: xyz[8][0] = x; xyz[8][1] = y; xyz[8][2] = z + dz;
242: xyz[9][0] = x + 2.*dx; xyz[9][1] = y; xyz[9][2] = z + dz;
243: xyz[10][0] = x; xyz[10][1] = y + 2.*dy; xyz[10][2] = z + dz;
244: xyz[11][0] = x + 2.*dx; xyz[11][1] = y + 2.*dy; xyz[11][2] = z + dz;
245: xyz[12][0] = x; xyz[12][1] = y; xyz[12][2] = z + 2.*dz;
246: xyz[13][0] = x + dx; xyz[13][1] = y; xyz[13][2] = z + 2.*dz;
247: xyz[14][0] = x + 2.*dx; xyz[14][1] = y; xyz[14][2] = z + 2.*dz;
248: xyz[15][0] = x; xyz[15][1] = y + dy; xyz[15][2] = z + 2.*dz;
249: xyz[16][0] = x + 2.*dx; xyz[16][1] = y + dy; xyz[16][2] = z + 2.*dz;
250: xyz[17][0] = x; xyz[17][1] = y + 2.*dy; xyz[17][2] = z + 2.*dz;
251: xyz[18][0] = x + dx; xyz[18][1] = y + 2.*dy; xyz[18][2] = z + 2.*dz;
252: xyz[19][0] = x + 2.*dx; xyz[19][1] = y + 2.*dy; xyz[19][2] = z + 2.*dz;
253: paulintegrate20(K);
255: /* copy the stiffness from K into format used by Ke */
256: for (i=0; i<81; i++) {
257: for (j=0; j<81; j++) {
258: Ke[i][j] = 0.0;
259: }
260: }
261: Ii = 0;
262: for (i=0; i<20; i++) {
263: J = 0;
264: for (j=0; j<20; j++) {
265: for (k=0; k<3; k++) {
266: for (l=0; l<3; l++) {
267: Ke[3*rmap[i]+k][3*rmap[j]+l] = v = K[Ii+k][J+l];
268: }
269: }
270: J += 3;
271: }
272: Ii += 3;
273: }
275: /* force the matrix to be exactly symmetric */
276: for (i=0; i<81; i++) {
277: for (j=0; j<i; j++) {
278: Ke[i][j] = (Ke[i][j] + Ke[j][i])/2.0;
279: }
280: }
281: return 0;
282: }
283: /* -------------------------------------------------------------------- */
286: /*
287: paulsetup20 - Sets up data structure for forming local elastic stiffness.
288: */
289: PetscErrorCode paulsetup20(void)
290: {
291: PetscInt i,j,k,cnt;
292: PetscReal x[4],w[4];
293: PetscReal c;
295: n_int = 3;
296: nu = 0.3;
297: E = 1.0;
299: /* Assign integration points and weights for
300: Gaussian quadrature formulae. */
301: if (n_int == 2) {
302: x[0] = (-0.577350269189626);
303: x[1] = (0.577350269189626);
304: w[0] = 1.0000000;
305: w[1] = 1.0000000;
306: } else if (n_int == 3) {
307: x[0] = (-0.774596669241483);
308: x[1] = 0.0000000;
309: x[2] = 0.774596669241483;
310: w[0] = 0.555555555555555;
311: w[1] = 0.888888888888888;
312: w[2] = 0.555555555555555;
313: } else if (n_int == 4) {
314: x[0] = (-0.861136311594053);
315: x[1] = (-0.339981043584856);
316: x[2] = 0.339981043584856;
317: x[3] = 0.861136311594053;
318: w[0] = 0.347854845137454;
319: w[1] = 0.652145154862546;
320: w[2] = 0.652145154862546;
321: w[3] = 0.347854845137454;
322: } else SETERRQ(PETSC_COMM_SELF,1,"Unknown value for n_int");
324: /* rst[][i] contains the location of the i-th integration point
325: in the canonical (r,s,t) coordinate system. weight[i] contains
326: the Gaussian weighting factor. */
328: cnt = 0;
329: for (i=0; i<n_int; i++) {
330: for (j=0; j<n_int; j++) {
331: for (k=0; k<n_int; k++) {
332: rst[0][cnt] =x[i];
333: rst[1][cnt] =x[j];
334: rst[2][cnt] =x[k];
335: weight[cnt] = w[i]*w[j]*w[k];
336: ++cnt;
337: }
338: }
339: }
340: N_int = cnt;
342: /* N[][j] is the interpolation vector, N[][j] .* xyz[] */
343: /* yields the (x,y,z) locations of the integration point. */
344: /* part_N[][][j] is the partials of the N function */
345: /* w.r.t. (r,s,t). */
347: c = 1.0/8.0;
348: for (j=0; j<N_int; j++) {
349: for (i=0; i<20; i++) {
350: if (i==0 || i==2 || i==5 || i==7 || i==12 || i==14 || i== 17 || i==19) {
351: N[i][j] = c*(1.0 + r2[i]*rst[0][j])*
352: (1.0 + s2[i]*rst[1][j])*(1.0 + t2[i]*rst[2][j])*
353: (-2.0 + r2[i]*rst[0][j] + s2[i]*rst[1][j] + t2[i]*rst[2][j]);
354: part_N[0][i][j] = c*r2[i]*(1 + s2[i]*rst[1][j])*(1 + t2[i]*rst[2][j])*
355: (-1.0 + 2.0*r2[i]*rst[0][j] + s2[i]*rst[1][j] +
356: t2[i]*rst[2][j]);
357: part_N[1][i][j] = c*s2[i]*(1 + r2[i]*rst[0][j])*(1 + t2[i]*rst[2][j])*
358: (-1.0 + r2[i]*rst[0][j] + 2.0*s2[i]*rst[1][j] +
359: t2[i]*rst[2][j]);
360: part_N[2][i][j] = c*t2[i]*(1 + r2[i]*rst[0][j])*(1 + s2[i]*rst[1][j])*
361: (-1.0 + r2[i]*rst[0][j] + s2[i]*rst[1][j] +
362: 2.0*t2[i]*rst[2][j]);
363: } else if (i==1 || i==6 || i==13 || i==18) {
364: N[i][j] = .25*(1.0 - rst[0][j]*rst[0][j])*
365: (1.0 + s2[i]*rst[1][j])*(1.0 + t2[i]*rst[2][j]);
366: part_N[0][i][j] = -.5*rst[0][j]*(1 + s2[i]*rst[1][j])*
367: (1 + t2[i]*rst[2][j]);
368: part_N[1][i][j] = .25*s2[i]*(1 + t2[i]*rst[2][j])*
369: (1.0 - rst[0][j]*rst[0][j]);
370: part_N[2][i][j] = .25*t2[i]*(1.0 - rst[0][j]*rst[0][j])*
371: (1 + s2[i]*rst[1][j]);
372: } else if (i==3 || i==4 || i==15 || i==16) {
373: N[i][j] = .25*(1.0 - rst[1][j]*rst[1][j])*
374: (1.0 + r2[i]*rst[0][j])*(1.0 + t2[i]*rst[2][j]);
375: part_N[0][i][j] = .25*r2[i]*(1 + t2[i]*rst[2][j])*
376: (1.0 - rst[1][j]*rst[1][j]);
377: part_N[1][i][j] = -.5*rst[1][j]*(1 + r2[i]*rst[0][j])*
378: (1 + t2[i]*rst[2][j]);
379: part_N[2][i][j] = .25*t2[i]*(1.0 - rst[1][j]*rst[1][j])*
380: (1 + r2[i]*rst[0][j]);
381: } else if (i==8 || i==9 || i==10 || i==11) {
382: N[i][j] = .25*(1.0 - rst[2][j]*rst[2][j])*
383: (1.0 + r2[i]*rst[0][j])*(1.0 + s2[i]*rst[1][j]);
384: part_N[0][i][j] = .25*r2[i]*(1 + s2[i]*rst[1][j])*
385: (1.0 - rst[2][j]*rst[2][j]);
386: part_N[1][i][j] = .25*s2[i]*(1.0 - rst[2][j]*rst[2][j])*
387: (1 + r2[i]*rst[0][j]);
388: part_N[2][i][j] = -.5*rst[2][j]*(1 + r2[i]*rst[0][j])*
389: (1 + s2[i]*rst[1][j]);
390: }
391: }
392: }
393: return 0;
394: }
395: /* -------------------------------------------------------------------- */
398: /*
399: paulintegrate20 - Does actual numerical integration on 20 node element.
400: */
401: PetscErrorCode paulintegrate20(PetscReal K[60][60])
402: {
403: PetscReal det_jac,jac[3][3],inv_jac[3][3];
404: PetscReal B[6][60],B_temp[6][60],C[6][6];
405: PetscReal temp;
406: PetscInt i,j,k,step;
408: /* Zero out K, since we will accumulate the result here */
409: for (i=0; i<60; i++) {
410: for (j=0; j<60; j++) {
411: K[i][j] = 0.0;
412: }
413: }
415: /* Loop over integration points ... */
416: for (step=0; step<N_int; step++) {
418: /* Compute the Jacobian, its determinant, and inverse. */
419: for (i=0; i<3; i++) {
420: for (j=0; j<3; j++) {
421: jac[i][j] = 0;
422: for (k=0; k<20; k++) {
423: jac[i][j] += part_N[i][k][step]*xyz[k][j];
424: }
425: }
426: }
427: det_jac = jac[0][0]*(jac[1][1]*jac[2][2]-jac[1][2]*jac[2][1])
428: + jac[0][1]*(jac[1][2]*jac[2][0]-jac[1][0]*jac[2][2])
429: + jac[0][2]*(jac[1][0]*jac[2][1]-jac[1][1]*jac[2][0]);
430: inv_jac[0][0] = (jac[1][1]*jac[2][2]-jac[1][2]*jac[2][1])/det_jac;
431: inv_jac[0][1] = (jac[0][2]*jac[2][1]-jac[0][1]*jac[2][2])/det_jac;
432: inv_jac[0][2] = (jac[0][1]*jac[1][2]-jac[1][1]*jac[0][2])/det_jac;
433: inv_jac[1][0] = (jac[1][2]*jac[2][0]-jac[1][0]*jac[2][2])/det_jac;
434: inv_jac[1][1] = (jac[0][0]*jac[2][2]-jac[2][0]*jac[0][2])/det_jac;
435: inv_jac[1][2] = (jac[0][2]*jac[1][0]-jac[0][0]*jac[1][2])/det_jac;
436: inv_jac[2][0] = (jac[1][0]*jac[2][1]-jac[1][1]*jac[2][0])/det_jac;
437: inv_jac[2][1] = (jac[0][1]*jac[2][0]-jac[0][0]*jac[2][1])/det_jac;
438: inv_jac[2][2] = (jac[0][0]*jac[1][1]-jac[1][0]*jac[0][1])/det_jac;
440: /* Compute the B matrix. */
441: for (i=0; i<3; i++) {
442: for (j=0; j<20; j++) {
443: B_temp[i][j] = 0.0;
444: for (k=0; k<3; k++) {
445: B_temp[i][j] += inv_jac[i][k]*part_N[k][j][step];
446: }
447: }
448: }
449: for (i=0; i<6; i++) {
450: for (j=0; j<60; j++) {
451: B[i][j] = 0.0;
452: }
453: }
455: /* Put values in correct places in B. */
456: for (k=0; k<20; k++) {
457: B[0][3*k] = B_temp[0][k];
458: B[1][3*k+1] = B_temp[1][k];
459: B[2][3*k+2] = B_temp[2][k];
460: B[3][3*k] = B_temp[1][k];
461: B[3][3*k+1] = B_temp[0][k];
462: B[4][3*k+1] = B_temp[2][k];
463: B[4][3*k+2] = B_temp[1][k];
464: B[5][3*k] = B_temp[2][k];
465: B[5][3*k+2] = B_temp[0][k];
466: }
468: /* Construct the C matrix, uses the constants "nu" and "E". */
469: for (i=0; i<6; i++) {
470: for (j=0; j<6; j++) {
471: C[i][j] = 0.0;
472: }
473: }
474: temp = (1.0 + nu)*(1.0 - 2.0*nu);
475: temp = E/temp;
476: C[0][0] = temp*(1.0 - nu);
477: C[1][1] = C[0][0];
478: C[2][2] = C[0][0];
479: C[3][3] = temp*(0.5 - nu);
480: C[4][4] = C[3][3];
481: C[5][5] = C[3][3];
482: C[0][1] = temp*nu;
483: C[0][2] = C[0][1];
484: C[1][0] = C[0][1];
485: C[1][2] = C[0][1];
486: C[2][0] = C[0][1];
487: C[2][1] = C[0][1];
489: for (i=0; i<6; i++) {
490: for (j=0; j<60; j++) {
491: B_temp[i][j] = 0.0;
492: for (k=0; k<6; k++) {
493: B_temp[i][j] += C[i][k]*B[k][j];
494: }
495: B_temp[i][j] *= det_jac;
496: }
497: }
499: /* Accumulate B'*C*B*det(J)*weight, as a function of (r,s,t), in K. */
500: for (i=0; i<60; i++) {
501: for (j=0; j<60; j++) {
502: temp = 0.0;
503: for (k=0; k<6; k++) {
504: temp += B[k][i]*B_temp[k][j];
505: }
506: K[i][j] += temp*weight[step];
507: }
508: }
509: } /* end of loop over integration points */
510: return 0;
511: }