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