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 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(PetscReal**);
 19: extern int AddElement(Mat,int,int,PetscReal**,int,int);
 20: extern int paulsetup20(void);
 21: extern int paulintegrate20(PetscReal 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:   PetscScalar v,neg1 = -1.0;
 28:   Vec         u,x,b;
 29:   SLES        sles;
 30:   KSP         ksp;
 31:   PetscReal   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,&u);
 45:   VecSetSizes(u,PETSC_DECIDE,rdim);
 46:   VecSetFromOptions(u);
 47:   VecDuplicate(u,&b);
 48:   VecDuplicate(b,&x);
 49:   for (i=rstart; i<rend; i++) {
 50:     v = (PetscScalar)(i-rstart + 100*rank);
 51:     VecSetValues(u,1,&i,&v,INSERT_VALUES);
 52:   }
 53:   VecAssemblyBegin(u);
 54:   VecAssemblyEnd(u);
 55: 
 56:   /* Compute right-hand-side */
 57:   MatMult(mat,u,b);
 58: 
 59:   /* Solve linear system */
 60:   SLESCreate(PETSC_COMM_WORLD,&sles);
 61:   SLESSetOperators(sles,mat,mat,SAME_NONZERO_PATTERN);
 62:   SLESGetKSP(sles,&ksp);
 63:   KSPGMRESSetRestart(ksp,2*m);
 64:   KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
 65:   KSPSetType(ksp,KSPCG);
 66:   SLESSetFromOptions(sles);
 67:   SLESSolve(sles,b,x,&its);
 68: 
 69:   /* Check error */
 70:   VecAXPY(&neg1,u,x);
 71:   VecNorm(x,NORM_2,&norm);

 73:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Number of iterations %dn",norm,its);

 75:   /* Free work space */
 76:   SLESDestroy(sles);
 77:   VecDestroy(u);
 78:   VecDestroy(x);
 79:   VecDestroy(b);
 80:   MatDestroy(mat);

 82:   PetscFinalize();
 83:   return 0;
 84: }
 85: /* -------------------------------------------------------------------- */
 86: /* 
 87:   GetElasticityMatrix - Forms 3D linear elasticity matrix.
 88:  */
 89: int GetElasticityMatrix(int m,Mat *newmat)
 90: {
 91:   int        i,j,k,i1,i2,j_1,j2,k1,k2,h1,h2,shiftx,shifty,shiftz;
 92:   int        ict,nz,base,r1,r2,N,*rowkeep,nstart,ierr;
 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=%dn",m,N);
101:   MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,80,PETSC_NULL,&mat);

103:   /* Form stiffness for element */
104:   PetscMalloc(81*sizeof(PetscReal *),&K);
105:   for (i=0; i<81; i++) {
106:     PetscMalloc(81*sizeof(PetscReal),&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:   PetscMalloc((N-nstart)*sizeof(int),&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,&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,newmat);
168:   MatDestroy(submat);

170:   PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
171:   MatView(*newmat,PETSC_VIEWER_STDOUT_WORLD);
172:   MatNorm(*newmat,NORM_1,&norm);
173:   PetscPrintf(PETSC_COMM_WORLD,"matrix 1 norm = %gn",norm);

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

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

227:   paulsetup20();

229:   x = -1.0;  y = -1.0; z = -1.0; dx = 2.0; dy = 2.0; dz = 2.0;
230:   xyz[0][0] = x;          xyz[0][1] = y;          xyz[0][2] = z;
231:   xyz[1][0] = x + dx;     xyz[1][1] = y;          xyz[1][2] = z;
232:   xyz[2][0] = x + 2.*dx;  xyz[2][1] = y;          xyz[2][2] = z;
233:   xyz[3][0] = x;          xyz[3][1] = y + dy;     xyz[3][2] = z;
234:   xyz[4][0] = x + 2.*dx;  xyz[4][1] = y + dy;     xyz[4][2] = z;
235:   xyz[5][0] = x;          xyz[5][1] = y + 2.*dy;  xyz[5][2] = z;
236:   xyz[6][0] = x + dx;     xyz[6][1] = y + 2.*dy;  xyz[6][2] = z;
237:   xyz[7][0] = x + 2.*dx;  xyz[7][1] = y + 2.*dy;  xyz[7][2] = z;
238:   xyz[8][0] = x;          xyz[8][1] = y;          xyz[8][2] = z + dz;
239:   xyz[9][0] = x + 2.*dx;  xyz[9][1] = y;          xyz[9][2] = z + dz;
240:   xyz[10][0] = x;         xyz[10][1] = y + 2.*dy; xyz[10][2] = z + dz;
241:   xyz[11][0] = x + 2.*dx; xyz[11][1] = y + 2.*dy; xyz[11][2] = z + dz;
242:   xyz[12][0] = x;         xyz[12][1] = y;         xyz[12][2] = z + 2.*dz;
243:   xyz[13][0] = x + dx;    xyz[13][1] = y;         xyz[13][2] = z + 2.*dz;
244:   xyz[14][0] = x + 2.*dx; xyz[14][1] = y;         xyz[14][2] = z + 2.*dz;
245:   xyz[15][0] = x;         xyz[15][1] = y + dy;    xyz[15][2] = z + 2.*dz;
246:   xyz[16][0] = x + 2.*dx; xyz[16][1] = y + dy;    xyz[16][2] = z + 2.*dz;
247:   xyz[17][0] = x;         xyz[17][1] = y + 2.*dy; xyz[17][2] = z + 2.*dz;
248:   xyz[18][0] = x + dx;    xyz[18][1] = y + 2.*dy; xyz[18][2] = z + 2.*dz;
249:   xyz[19][0] = x + 2.*dx; xyz[19][1] = y + 2.*dy; xyz[19][2] = z + 2.*dz;
250:   paulintegrate20(K);

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

298:   n_int = 3;
299:   nu = 0.3;
300:   E = 1.0;

302:   /* Assign integration points and weights for
303:        Gaussian quadrature formulae. */
304:   if(n_int == 2)  {
305:                 x[0] = (-0.577350269189626);
306:                 x[1] = (0.577350269189626);
307:                 w[0] = 1.0000000;
308:                 w[1] = 1.0000000;
309:   }
310:   else if(n_int == 3) {
311:                 x[0] = (-0.774596669241483);
312:                 x[1] = 0.0000000;
313:                 x[2] = 0.774596669241483;
314:                 w[0] = 0.555555555555555;
315:                 w[1] = 0.888888888888888;
316:                 w[2] = 0.555555555555555;
317:   }
318:   else if(n_int == 4) {
319:                 x[0] = (-0.861136311594053);
320:                 x[1] = (-0.339981043584856);
321:                 x[2] = 0.339981043584856;
322:                 x[3] = 0.861136311594053;
323:                 w[0] = 0.347854845137454;
324:                 w[1] = 0.652145154862546;
325:                 w[2] = 0.652145154862546;
326:                 w[3] = 0.347854845137454;
327:   }
328:   else {
329:     SETERRQ(1,"Unknown value for n_int");
330:   }

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

336:   cnt = 0;
337:   for (i=0; i<n_int;i++) {
338:     for (j=0; j<n_int;j++) {
339:       for (k=0; k<n_int;k++) {
340:         rst[0][cnt]=x[i];
341:         rst[1][cnt]=x[j];
342:         rst[2][cnt]=x[k];
343:         weight[cnt] = w[i]*w[j]*w[k];
344:         ++cnt;
345:       }
346:     }
347:   }
348:   N_int = cnt;

350:   /* N[][j] is the interpolation vector, N[][j] .* xyz[] */
351:   /* yields the (x,y,z)  locations of the integration point. */
352:   /*  part_N[][][j] is the partials of the N function */
353:   /*  w.r.t. (r,s,t). */

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

417:   /* Zero out K, since we will accumulate the result here */
418:   for (i=0; i<60; i++) {
419:     for (j=0; j<60; j++) {
420:       K[i][j] = 0.0;
421:     }
422:   }

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

427:     /* Compute the Jacobian, its determinant, and inverse. */
428:     for (i=0; i<3; i++) {
429:       for (j=0; j<3; j++) {
430:         jac[i][j] = 0;
431:         for (k=0; k<20; k++) {
432:           jac[i][j] += part_N[i][k][step]*xyz[k][j];
433:         }
434:       }
435:     }
436:     det_jac = jac[0][0]*(jac[1][1]*jac[2][2]-jac[1][2]*jac[2][1])
437:               + jac[0][1]*(jac[1][2]*jac[2][0]-jac[1][0]*jac[2][2])
438:               + jac[0][2]*(jac[1][0]*jac[2][1]-jac[1][1]*jac[2][0]);
439:     inv_jac[0][0] = (jac[1][1]*jac[2][2]-jac[1][2]*jac[2][1])/det_jac;
440:     inv_jac[0][1] = (jac[0][2]*jac[2][1]-jac[0][1]*jac[2][2])/det_jac;
441:     inv_jac[0][2] = (jac[0][1]*jac[1][2]-jac[1][1]*jac[0][2])/det_jac;
442:     inv_jac[1][0] = (jac[1][2]*jac[2][0]-jac[1][0]*jac[2][2])/det_jac;
443:     inv_jac[1][1] = (jac[0][0]*jac[2][2]-jac[2][0]*jac[0][2])/det_jac;
444:     inv_jac[1][2] = (jac[0][2]*jac[1][0]-jac[0][0]*jac[1][2])/det_jac;
445:     inv_jac[2][0] = (jac[1][0]*jac[2][1]-jac[1][1]*jac[2][0])/det_jac;
446:     inv_jac[2][1] = (jac[0][1]*jac[2][0]-jac[0][0]*jac[2][1])/det_jac;
447:     inv_jac[2][2] = (jac[0][0]*jac[1][1]-jac[1][0]*jac[0][1])/det_jac;

449:     /* Compute the B matrix. */
450:     for (i=0; i<3; i++) {
451:       for (j=0; j<20; j++) {
452:         B_temp[i][j] = 0.0;
453:         for (k=0; k<3; k++) {
454:           B_temp[i][j] += inv_jac[i][k]*part_N[k][j][step];
455:         }
456:       }
457:     }
458:     for (i=0; i<6; i++) {
459:       for (j=0; j<60; j++) {
460:         B[i][j] = 0.0;
461:       }
462:     }

464:     /* Put values in correct places in B. */
465:     for (k=0; k<20; k++) {
466:       B[0][3*k]   = B_temp[0][k];
467:       B[1][3*k+1] = B_temp[1][k];
468:       B[2][3*k+2] = B_temp[2][k];
469:       B[3][3*k]   = B_temp[1][k];
470:       B[3][3*k+1] = B_temp[0][k];
471:       B[4][3*k+1] = B_temp[2][k];
472:       B[4][3*k+2] = B_temp[1][k];
473:       B[5][3*k]   = B_temp[2][k];
474:       B[5][3*k+2] = B_temp[0][k];
475:     }
476: 
477:     /* Construct the C matrix, uses the constants "nu" and "E". */
478:     for (i=0; i<6; i++) {
479:       for (j=0; j<6; j++) {
480:         C[i][j] = 0.0;
481:       }
482:     }
483:     temp = (1.0 + nu)*(1.0 - 2.0*nu);
484:     temp = E/temp;
485:     C[0][0] = temp*(1.0 - nu);
486:     C[1][1] = C[0][0];
487:     C[2][2] = C[0][0];
488:     C[3][3] = temp*(0.5 - nu);
489:     C[4][4] = C[3][3];
490:     C[5][5] = C[3][3];
491:     C[0][1] = temp*nu;
492:     C[0][2] = C[0][1];
493:     C[1][0] = C[0][1];
494:     C[1][2] = C[0][1];
495:     C[2][0] = C[0][1];
496:     C[2][1] = C[0][1];
497: 
498:     for (i=0; i<6; i++) {
499:       for (j=0; j<60; j++) {
500:         B_temp[i][j] = 0.0;
501:         for (k=0; k<6; k++) {
502:           B_temp[i][j] += C[i][k]*B[k][j];
503:         }
504:         B_temp[i][j] *= det_jac;
505:       }
506:     }

508:     /* Accumulate B'*C*B*det(J)*weight, as a function of (r,s,t), in K. */
509:     for (i=0; i<60; i++) {
510:       for (j=0; j<60; j++) {
511:         temp = 0.0;
512:         for (k=0; k<6; k++) {
513:           temp += B[k][i]*B_temp[k][j];
514:         }
515:         K[i][j] += temp*weight[step];
516:       }
517:     }
518:   }  /* end of loop over integration points */
519:   return 0;
520: }