Actual source code: ex10.c

  1: /*$Id: ex10.c,v 1.97 2001/08/24 16:21:45 bsmith Exp $*/

  3: static char help[] = "This example calculates the stiffness matrix for a brick in three\n\
  4: dimensions using 20 node serendipity elements and the equations of linear\n\
  5: elasticity. This also demonstrates use of  block\n\
  6: diagonal data structure.  Input arguments are:\n\
  7:   -m : problem size\n\n";

 9:  #include petscksp.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(PetscReal**);
 19: extern int AddElement(Mat,int,int,PetscReal**,int,int);
 20: extern int paulsetup20(void);
 21: extern int paulintegrate20(PetscReal K[60][60]);

 25: int main(int argc,char **args)
 26: {
 27:   Mat         mat;
 28:   int         ierr,i,its,m = 3,rdim,cdim,rstart,rend,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(PETSC_NULL,"-m",&m,PETSC_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);
 56: 
 57:   /* Compute right-hand-side */
 58:   MatMult(mat,u,b);
 59: 
 60:   /* Solve linear system */
 61:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 62:   KSPSetOperators(ksp,mat,mat,SAME_NONZERO_PATTERN);
 63:   KSPGMRESSetRestart(ksp,2*m);
 64:   KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
 65:   KSPSetType(ksp,KSPCG);
 66:   KSPSetFromOptions(ksp);
 67:   KSPSetRhs(ksp,b);
 68:   KSPSetSolution(ksp,x);
 69:   KSPSolve(ksp);
 70:   KSPGetIterationNumber(ksp,&its);
 71:   /* Check error */
 72:   VecAXPY(&neg1,u,x);
 73:   VecNorm(x,NORM_2,&norm);

 75:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Number of iterations %d\n",norm,its);

 77:   /* Free work space */
 78:   KSPDestroy(ksp);
 79:   VecDestroy(u);
 80:   VecDestroy(x);
 81:   VecDestroy(b);
 82:   MatDestroy(mat);

 84:   PetscFinalize();
 85:   return 0;
 86: }
 87: /* -------------------------------------------------------------------- */
 90: /* 
 91:   GetElasticityMatrix - Forms 3D linear elasticity matrix.
 92:  */
 93: int GetElasticityMatrix(int m,Mat *newmat)
 94: {
 95:   int           i,j,k,i1,i2,j_1,j2,k1,k2,h1,h2,shiftx,shifty,shiftz;
 96:   int           ict,nz,base,r1,r2,N,*rowkeep,nstart,ierr;
 97:   IS            iskeep;
 98:   PetscReal     **K,norm;
 99:   Mat           mat,submat = 0,*submatb;
100:   const MatType type = MATSEQBAIJ;

102:   m /= 2;   /* This is done just to be consistent with the old example */
103:   N = 3*(2*m+1)*(2*m+1)*(2*m+1);
104:   PetscPrintf(PETSC_COMM_SELF,"m = %d, N=%d\n",m,N);
105:   MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,80,PETSC_NULL,&mat);

107:   /* Form stiffness for element */
108:   PetscMalloc(81*sizeof(PetscReal *),&K);
109:   for (i=0; i<81; i++) {
110:     PetscMalloc(81*sizeof(PetscReal),&K[i]);
111:   }
112:   Elastic20Stiff(K);

114:   /* Loop over elements and add contribution to stiffness */
115:   shiftx = 3; shifty = 3*(2*m+1); shiftz = 3*(2*m+1)*(2*m+1);
116:   for (k=0; k<m; k++) {
117:     for (j=0; j<m; j++) {
118:       for (i=0; i<m; i++) {
119:         h1 = 0;
120:         base = 2*k*shiftz + 2*j*shifty + 2*i*shiftx;
121:         for (k1=0; k1<3; k1++) {
122:           for (j_1=0; j_1<3; j_1++) {
123:             for (i1=0; i1<3; i1++) {
124:               h2 = 0;
125:               r1 = base + i1*shiftx + j_1*shifty + k1*shiftz;
126:               for (k2=0; k2<3; k2++) {
127:                 for (j2=0; j2<3; j2++) {
128:                   for (i2=0; i2<3; i2++) {
129:                     r2 = base + i2*shiftx + j2*shifty + k2*shiftz;
130:                     AddElement(mat,r1,r2,K,h1,h2);
131:                     h2 += 3;
132:                   }
133:                 }
134:               }
135:               h1 += 3;
136:             }
137:           }
138:         }
139:       }
140:     }
141:   }

143:   for (i=0; i<81; i++) {
144:     PetscFree(K[i]);
145:   }
146:   PetscFree(K);

148:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
149:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

151:   /* Exclude any superfluous rows and columns */
152:   nstart = 3*(2*m+1)*(2*m+1);
153:   ict = 0;
154:   PetscMalloc((N-nstart)*sizeof(int),&rowkeep);
155:   for (i=nstart; i<N; i++) {
156:     MatGetRow(mat,i,&nz,0,0);
157:     if (nz) rowkeep[ict++] = i;
158:     MatRestoreRow(mat,i,&nz,0,0);
159:   }
160:   ISCreateGeneral(PETSC_COMM_SELF,ict,rowkeep,&iskeep);
161:   MatGetSubMatrices(mat,1,&iskeep,&iskeep,MAT_INITIAL_MATRIX,&submatb);
162:   submat = *submatb;
163:   PetscFree(submatb);
164:   PetscFree(rowkeep);
165:   ISDestroy(iskeep);
166:   MatDestroy(mat);

168:   /* Convert storage formats -- just to demonstrate conversion to various
169:      formats (in particular, block diagonal storage).  This is NOT the
170:      recommended means to solve such a problem.  */
171:   MatConvert(submat,type,newmat);
172:   MatDestroy(submat);

174:   PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
175:   MatView(*newmat,PETSC_VIEWER_STDOUT_WORLD);
176:   MatNorm(*newmat,NORM_1,&norm);
177:   PetscPrintf(PETSC_COMM_WORLD,"matrix 1 norm = %g\n",norm);

179:   return 0;
180: }
181: /* -------------------------------------------------------------------- */
184: int AddElement(Mat mat,int r1,int r2,PetscReal **K,int h1,int h2)
185: {
186:   PetscScalar val;
187:   int    l1,l2,row,col,ierr;

189:   for (l1=0; l1<3; l1++) {
190:     for (l2=0; l2<3; l2++) {
191: /*
192:    NOTE you should never do this! Inserting values 1 at a time is 
193:    just too expensive!
194: */
195:       if (K[h1+l1][h2+l2] != 0.0) {
196:         row = r1+l1; col = r2+l2; val = K[h1+l1][h2+l2];
197:         MatSetValues(mat,1,&row,1,&col,&val,ADD_VALUES);
198:         row = r2+l2; col = r1+l1;
199:         MatSetValues(mat,1,&row,1,&col,&val,ADD_VALUES);
200:       }
201:     }
202:   }
203:   return 0;
204: }
205: /* -------------------------------------------------------------------- */
206: PetscReal        N[20][64];           /* Interpolation function. */
207: PetscReal        part_N[3][20][64]; /* Partials of interpolation function. */
208: PetscReal        rst[3][64];           /* Location of integration pts in (r,s,t) */
209: PetscReal        weight[64];           /* Gaussian quadrature weights. */
210: PetscReal        xyz[20][3];           /* (x,y,z) coordinates of nodes  */
211: PetscReal        E,nu;                   /* Physcial constants. */
212: int        n_int,N_int;           /* N_int = n_int^3, number of int. pts. */
213: /* Ordering of the vertices, (r,s,t) coordinates, of the canonical cell. */
214: PetscReal        r2[20] = {-1.0,0.0,1.0,-1.0,1.0,-1.0,0.0,1.0,
215:                  -1.0,1.0,-1.0,1.0,
216:                  -1.0,0.0,1.0,-1.0,1.0,-1.0,0.0,1.0};
217: PetscReal        s2[20] = {-1.0,-1.0, -1.0,0.0,0.0,1.0, 1.0, 1.0,
218:                  -1.0,-1.0,1.0,1.0,
219:                  -1.0,-1.0, -1.0,0.0,0.0,1.0, 1.0, 1.0};
220: PetscReal        t2[20] =  {-1.0,-1.0,-1.0,-1.0,-1.0,-1.0,-1.0,-1.0,
221:                  0.0,0.0,0.0,0.0,
222:                  1.0,1.0,1.0,1.0,1.0,1.0,1.0,1.0};
223: int     rmap[20] = {0,1,2,3,5,6,7,8,9,11,15,17,18,19,20,21,23,24,25,26};
224: /* -------------------------------------------------------------------- */
227: /* 
228:   Elastic20Stiff - Forms 20 node elastic stiffness for element.
229:  */
230: int Elastic20Stiff(PetscReal **Ke)
231: {
232:   PetscReal K[60][60],x,y,z,dx,dy,dz,m,v;
233:   int       i,j,k,l,I,J;

235:   paulsetup20();

237:   x = -1.0;  y = -1.0; z = -1.0; dx = 2.0; dy = 2.0; dz = 2.0;
238:   xyz[0][0] = x;          xyz[0][1] = y;          xyz[0][2] = z;
239:   xyz[1][0] = x + dx;     xyz[1][1] = y;          xyz[1][2] = z;
240:   xyz[2][0] = x + 2.*dx;  xyz[2][1] = y;          xyz[2][2] = z;
241:   xyz[3][0] = x;          xyz[3][1] = y + dy;     xyz[3][2] = z;
242:   xyz[4][0] = x + 2.*dx;  xyz[4][1] = y + dy;     xyz[4][2] = z;
243:   xyz[5][0] = x;          xyz[5][1] = y + 2.*dy;  xyz[5][2] = z;
244:   xyz[6][0] = x + dx;     xyz[6][1] = y + 2.*dy;  xyz[6][2] = z;
245:   xyz[7][0] = x + 2.*dx;  xyz[7][1] = y + 2.*dy;  xyz[7][2] = z;
246:   xyz[8][0] = x;          xyz[8][1] = y;          xyz[8][2] = z + dz;
247:   xyz[9][0] = x + 2.*dx;  xyz[9][1] = y;          xyz[9][2] = z + dz;
248:   xyz[10][0] = x;         xyz[10][1] = y + 2.*dy; xyz[10][2] = z + dz;
249:   xyz[11][0] = x + 2.*dx; xyz[11][1] = y + 2.*dy; xyz[11][2] = z + dz;
250:   xyz[12][0] = x;         xyz[12][1] = y;         xyz[12][2] = z + 2.*dz;
251:   xyz[13][0] = x + dx;    xyz[13][1] = y;         xyz[13][2] = z + 2.*dz;
252:   xyz[14][0] = x + 2.*dx; xyz[14][1] = y;         xyz[14][2] = z + 2.*dz;
253:   xyz[15][0] = x;         xyz[15][1] = y + dy;    xyz[15][2] = z + 2.*dz;
254:   xyz[16][0] = x + 2.*dx; xyz[16][1] = y + dy;    xyz[16][2] = z + 2.*dz;
255:   xyz[17][0] = x;         xyz[17][1] = y + 2.*dy; xyz[17][2] = z + 2.*dz;
256:   xyz[18][0] = x + dx;    xyz[18][1] = y + 2.*dy; xyz[18][2] = z + 2.*dz;
257:   xyz[19][0] = x + 2.*dx; xyz[19][1] = y + 2.*dy; xyz[19][2] = z + 2.*dz;
258:   paulintegrate20(K);

260:   /* copy the stiffness from K into format used by Ke */
261:   for (i=0; i<81; i++) {
262:     for (j=0; j<81; j++) {
263:           Ke[i][j] = 0.0;
264:     }
265:   }
266:   I = 0;
267:   m = 0.0;
268:   for (i=0; i<20; i++) {
269:     J = 0;
270:     for (j=0; j<20; j++) {
271:       for (k=0; k<3; k++) {
272:         for (l=0; l<3; l++) {
273:           Ke[3*rmap[i]+k][3*rmap[j]+l] = v = K[I+k][J+l];
274:           m = PetscMax(m,PetscAbsReal(v));
275:         }
276:       }
277:       J += 3;
278:     }
279:     I += 3;
280:   }
281:   /* zero out the extremely small values */
282:   m = (1.e-8)*m;
283:   for (i=0; i<81; i++) {
284:     for (j=0; j<81; j++) {
285:       if (PetscAbsReal(Ke[i][j]) < m)  Ke[i][j] = 0.0;
286:     }
287:   }
288:   /* force the matrix to be exactly symmetric */
289:   for (i=0; i<81; i++) {
290:     for (j=0; j<i; j++) {
291:       Ke[i][j] = (Ke[i][j] + Ke[j][i])/2.0;
292:     }
293:   }
294:   return 0;
295: }
296: /* -------------------------------------------------------------------- */
299: /* 
300:   paulsetup20 - Sets up data structure for forming local elastic stiffness.
301:  */
302: int paulsetup20(void)
303: {
304:   int       i,j,k,cnt;
305:   PetscReal x[4],w[4];
306:   PetscReal c;

308:   n_int = 3;
309:   nu = 0.3;
310:   E = 1.0;

312:   /* Assign integration points and weights for
313:        Gaussian quadrature formulae. */
314:   if(n_int == 2)  {
315:                 x[0] = (-0.577350269189626);
316:                 x[1] = (0.577350269189626);
317:                 w[0] = 1.0000000;
318:                 w[1] = 1.0000000;
319:   }
320:   else if(n_int == 3) {
321:                 x[0] = (-0.774596669241483);
322:                 x[1] = 0.0000000;
323:                 x[2] = 0.774596669241483;
324:                 w[0] = 0.555555555555555;
325:                 w[1] = 0.888888888888888;
326:                 w[2] = 0.555555555555555;
327:   }
328:   else if(n_int == 4) {
329:                 x[0] = (-0.861136311594053);
330:                 x[1] = (-0.339981043584856);
331:                 x[2] = 0.339981043584856;
332:                 x[3] = 0.861136311594053;
333:                 w[0] = 0.347854845137454;
334:                 w[1] = 0.652145154862546;
335:                 w[2] = 0.652145154862546;
336:                 w[3] = 0.347854845137454;
337:   }
338:   else {
339:     SETERRQ(1,"Unknown value for n_int");
340:   }

342:   /* rst[][i] contains the location of the i-th integration point
343:       in the canonical (r,s,t) coordinate system.  weight[i] contains
344:       the Gaussian weighting factor. */

346:   cnt = 0;
347:   for (i=0; i<n_int;i++) {
348:     for (j=0; j<n_int;j++) {
349:       for (k=0; k<n_int;k++) {
350:         rst[0][cnt]=x[i];
351:         rst[1][cnt]=x[j];
352:         rst[2][cnt]=x[k];
353:         weight[cnt] = w[i]*w[j]*w[k];
354:         ++cnt;
355:       }
356:     }
357:   }
358:   N_int = cnt;

360:   /* N[][j] is the interpolation vector, N[][j] .* xyz[] */
361:   /* yields the (x,y,z)  locations of the integration point. */
362:   /*  part_N[][][j] is the partials of the N function */
363:   /*  w.r.t. (r,s,t). */

365:   c = 1.0/8.0;
366:   for (j=0; j<N_int; j++) {
367:     for (i=0; i<20; i++) {
368:       if (i==0 || i==2 || i==5 || i==7 || i==12 || i==14 || i== 17 || i==19){
369:         N[i][j] = c*(1.0 + r2[i]*rst[0][j])*
370:                 (1.0 + s2[i]*rst[1][j])*(1.0 + t2[i]*rst[2][j])*
371:                 (-2.0 + r2[i]*rst[0][j] + s2[i]*rst[1][j] + t2[i]*rst[2][j]);
372:         part_N[0][i][j] = c*r2[i]*(1 + s2[i]*rst[1][j])*(1 + t2[i]*rst[2][j])*
373:                  (-1.0 + 2.0*r2[i]*rst[0][j] + s2[i]*rst[1][j] +
374:                  t2[i]*rst[2][j]);
375:         part_N[1][i][j] = c*s2[i]*(1 + r2[i]*rst[0][j])*(1 + t2[i]*rst[2][j])*
376:                  (-1.0 + r2[i]*rst[0][j] + 2.0*s2[i]*rst[1][j] +
377:                  t2[i]*rst[2][j]);
378:         part_N[2][i][j] = c*t2[i]*(1 + r2[i]*rst[0][j])*(1 + s2[i]*rst[1][j])*
379:                  (-1.0 + r2[i]*rst[0][j] + s2[i]*rst[1][j] +
380:                  2.0*t2[i]*rst[2][j]);
381:       }
382:       else if (i==1 || i==6 || i==13 || i==18) {
383:         N[i][j] = .25*(1.0 - rst[0][j]*rst[0][j])*
384:                 (1.0 + s2[i]*rst[1][j])*(1.0 + t2[i]*rst[2][j]);
385:         part_N[0][i][j] = -.5*rst[0][j]*(1 + s2[i]*rst[1][j])*
386:                          (1 + t2[i]*rst[2][j]);
387:         part_N[1][i][j] = .25*s2[i]*(1 + t2[i]*rst[2][j])*
388:                           (1.0 - rst[0][j]*rst[0][j]);
389:         part_N[2][i][j] = .25*t2[i]*(1.0 - rst[0][j]*rst[0][j])*
390:                           (1 + s2[i]*rst[1][j]);
391:       }
392:       else if (i==3 || i==4 || i==15 || i==16) {
393:         N[i][j] = .25*(1.0 - rst[1][j]*rst[1][j])*
394:                 (1.0 + r2[i]*rst[0][j])*(1.0 + t2[i]*rst[2][j]);
395:         part_N[0][i][j] = .25*r2[i]*(1 + t2[i]*rst[2][j])*
396:                            (1.0 - rst[1][j]*rst[1][j]);
397:         part_N[1][i][j] = -.5*rst[1][j]*(1 + r2[i]*rst[0][j])*
398:                            (1 + t2[i]*rst[2][j]);
399:         part_N[2][i][j] = .25*t2[i]*(1.0 - rst[1][j]*rst[1][j])*
400:                           (1 + r2[i]*rst[0][j]);
401:       }
402:       else if (i==8 || i==9 || i==10 || i==11) {
403:         N[i][j] = .25*(1.0 - rst[2][j]*rst[2][j])*
404:                 (1.0 + r2[i]*rst[0][j])*(1.0 + s2[i]*rst[1][j]);
405:         part_N[0][i][j] = .25*r2[i]*(1 + s2[i]*rst[1][j])*
406:                           (1.0 - rst[2][j]*rst[2][j]);
407:         part_N[1][i][j] = .25*s2[i]*(1.0 - rst[2][j]*rst[2][j])*
408:                           (1 + r2[i]*rst[0][j]);
409:         part_N[2][i][j] = -.5*rst[2][j]*(1 + r2[i]*rst[0][j])*
410:                            (1 + s2[i]*rst[1][j]);
411:       }
412:     }
413:   }
414:   return 0;
415: }
416: /* -------------------------------------------------------------------- */
419: /* 
420:    paulintegrate20 - Does actual numerical integration on 20 node element.
421:  */
422: int paulintegrate20(PetscReal K[60][60])
423: {
424:   PetscReal  det_jac,jac[3][3],inv_jac[3][3];
425:   PetscReal  B[6][60],B_temp[6][60],C[6][6];
426:   PetscReal  temp;
427:   int        i,j,k,step;

429:   /* Zero out K, since we will accumulate the result here */
430:   for (i=0; i<60; i++) {
431:     for (j=0; j<60; j++) {
432:       K[i][j] = 0.0;
433:     }
434:   }

436:   /* Loop over integration points ... */
437:   for (step=0; step<N_int; step++) {

439:     /* Compute the Jacobian, its determinant, and inverse. */
440:     for (i=0; i<3; i++) {
441:       for (j=0; j<3; j++) {
442:         jac[i][j] = 0;
443:         for (k=0; k<20; k++) {
444:           jac[i][j] += part_N[i][k][step]*xyz[k][j];
445:         }
446:       }
447:     }
448:     det_jac = jac[0][0]*(jac[1][1]*jac[2][2]-jac[1][2]*jac[2][1])
449:               + jac[0][1]*(jac[1][2]*jac[2][0]-jac[1][0]*jac[2][2])
450:               + jac[0][2]*(jac[1][0]*jac[2][1]-jac[1][1]*jac[2][0]);
451:     inv_jac[0][0] = (jac[1][1]*jac[2][2]-jac[1][2]*jac[2][1])/det_jac;
452:     inv_jac[0][1] = (jac[0][2]*jac[2][1]-jac[0][1]*jac[2][2])/det_jac;
453:     inv_jac[0][2] = (jac[0][1]*jac[1][2]-jac[1][1]*jac[0][2])/det_jac;
454:     inv_jac[1][0] = (jac[1][2]*jac[2][0]-jac[1][0]*jac[2][2])/det_jac;
455:     inv_jac[1][1] = (jac[0][0]*jac[2][2]-jac[2][0]*jac[0][2])/det_jac;
456:     inv_jac[1][2] = (jac[0][2]*jac[1][0]-jac[0][0]*jac[1][2])/det_jac;
457:     inv_jac[2][0] = (jac[1][0]*jac[2][1]-jac[1][1]*jac[2][0])/det_jac;
458:     inv_jac[2][1] = (jac[0][1]*jac[2][0]-jac[0][0]*jac[2][1])/det_jac;
459:     inv_jac[2][2] = (jac[0][0]*jac[1][1]-jac[1][0]*jac[0][1])/det_jac;

461:     /* Compute the B matrix. */
462:     for (i=0; i<3; i++) {
463:       for (j=0; j<20; j++) {
464:         B_temp[i][j] = 0.0;
465:         for (k=0; k<3; k++) {
466:           B_temp[i][j] += inv_jac[i][k]*part_N[k][j][step];
467:         }
468:       }
469:     }
470:     for (i=0; i<6; i++) {
471:       for (j=0; j<60; j++) {
472:         B[i][j] = 0.0;
473:       }
474:     }

476:     /* Put values in correct places in B. */
477:     for (k=0; k<20; k++) {
478:       B[0][3*k]   = B_temp[0][k];
479:       B[1][3*k+1] = B_temp[1][k];
480:       B[2][3*k+2] = B_temp[2][k];
481:       B[3][3*k]   = B_temp[1][k];
482:       B[3][3*k+1] = B_temp[0][k];
483:       B[4][3*k+1] = B_temp[2][k];
484:       B[4][3*k+2] = B_temp[1][k];
485:       B[5][3*k]   = B_temp[2][k];
486:       B[5][3*k+2] = B_temp[0][k];
487:     }
488: 
489:     /* Construct the C matrix, uses the constants "nu" and "E". */
490:     for (i=0; i<6; i++) {
491:       for (j=0; j<6; j++) {
492:         C[i][j] = 0.0;
493:       }
494:     }
495:     temp = (1.0 + nu)*(1.0 - 2.0*nu);
496:     temp = E/temp;
497:     C[0][0] = temp*(1.0 - nu);
498:     C[1][1] = C[0][0];
499:     C[2][2] = C[0][0];
500:     C[3][3] = temp*(0.5 - nu);
501:     C[4][4] = C[3][3];
502:     C[5][5] = C[3][3];
503:     C[0][1] = temp*nu;
504:     C[0][2] = C[0][1];
505:     C[1][0] = C[0][1];
506:     C[1][2] = C[0][1];
507:     C[2][0] = C[0][1];
508:     C[2][1] = C[0][1];
509: 
510:     for (i=0; i<6; i++) {
511:       for (j=0; j<60; j++) {
512:         B_temp[i][j] = 0.0;
513:         for (k=0; k<6; k++) {
514:           B_temp[i][j] += C[i][k]*B[k][j];
515:         }
516:         B_temp[i][j] *= det_jac;
517:       }
518:     }

520:     /* Accumulate B'*C*B*det(J)*weight, as a function of (r,s,t), in K. */
521:     for (i=0; i<60; i++) {
522:       for (j=0; j<60; j++) {
523:         temp = 0.0;
524:         for (k=0; k<6; k++) {
525:           temp += B[k][i]*B_temp[k][j];
526:         }
527:         K[i][j] += temp*weight[step];
528:       }
529:     }
530:   }  /* end of loop over integration points */
531:   return 0;
532: }