Actual source code: ex10.c

petsc-3.4.2 2013-07-02
  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,"-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,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:   KSPSolve(ksp,b,x);
 68:   KSPGetIterationNumber(ksp,&its);
 69:   /* Check error */
 70:   VecAXPY(x,neg1,u);
 71:   VecNorm(x,NORM_2,&norm);

 73:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %G Number of iterations %D\n",norm,its);

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

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

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

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

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

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

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

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

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

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

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

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: PetscInt  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: PetscInt  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: PetscErrorCode Elastic20Stiff(PetscReal **Ke)
231: {
232:   PetscReal K[60][60],x,y,z,dx,dy,dz,m,v;
233:   PetscInt  i,j,k,l,Ii,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:   Ii = 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[Ii+k][J+l];

275:           m = PetscMax(m,PetscAbsReal(v));
276:         }
277:       }
278:       J += 3;
279:     }
280:     Ii += 3;
281:   }
282:   /* zero out the extremely small values */
283:   m = (1.e-8)*m;
284:   for (i=0; i<81; i++) {
285:     for (j=0; j<81; j++) {
286:       if (PetscAbsReal(Ke[i][j]) < m) Ke[i][j] = 0.0;
287:     }
288:   }
289:   /* force the matrix to be exactly symmetric */
290:   for (i=0; i<81; i++) {
291:     for (j=0; j<i; j++) {
292:       Ke[i][j] = (Ke[i][j] + Ke[j][i])/2.0;
293:     }
294:   }
295:   return 0;
296: }
297: /* -------------------------------------------------------------------- */
300: /*
301:   paulsetup20 - Sets up data structure for forming local elastic stiffness.
302:  */
303: PetscErrorCode paulsetup20(void)
304: {
305:   PetscInt  i,j,k,cnt;
306:   PetscReal x[4],w[4];
307:   PetscReal c;

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

313:   /* Assign integration points and weights for
314:        Gaussian quadrature formulae. */
315:   if (n_int == 2) {
316:     x[0] = (-0.577350269189626);
317:     x[1] = (0.577350269189626);
318:     w[0] = 1.0000000;
319:     w[1] = 1.0000000;
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:   } else if (n_int == 4) {
328:     x[0] = (-0.861136311594053);
329:     x[1] = (-0.339981043584856);
330:     x[2] = 0.339981043584856;
331:     x[3] = 0.861136311594053;
332:     w[0] = 0.347854845137454;
333:     w[1] = 0.652145154862546;
334:     w[2] = 0.652145154862546;
335:     w[3] = 0.347854845137454;
336:   } else SETERRQ(PETSC_COMM_SELF,1,"Unknown value for n_int");

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

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

356:   /* N[][j] is the interpolation vector, N[][j] .* xyz[] */
357:   /* yields the (x,y,z)  locations of the integration point. */
358:   /*  part_N[][][j] is the partials of the N function */
359:   /*  w.r.t. (r,s,t). */

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

422:   /* Zero out K, since we will accumulate the result here */
423:   for (i=0; i<60; i++) {
424:     for (j=0; j<60; j++) {
425:       K[i][j] = 0.0;
426:     }
427:   }

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

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

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

469:     /* Put values in correct places in B. */
470:     for (k=0; k<20; k++) {
471:       B[0][3*k]   = B_temp[0][k];
472:       B[1][3*k+1] = B_temp[1][k];
473:       B[2][3*k+2] = B_temp[2][k];
474:       B[3][3*k]   = B_temp[1][k];
475:       B[3][3*k+1] = B_temp[0][k];
476:       B[4][3*k+1] = B_temp[2][k];
477:       B[4][3*k+2] = B_temp[1][k];
478:       B[5][3*k]   = B_temp[2][k];
479:       B[5][3*k+2] = B_temp[0][k];
480:     }

482:     /* Construct the C matrix, uses the constants "nu" and "E". */
483:     for (i=0; i<6; i++) {
484:       for (j=0; j<6; j++) {
485:         C[i][j] = 0.0;
486:       }
487:     }
488:     temp = (1.0 + nu)*(1.0 - 2.0*nu);
489:     temp = E/temp;
490:     C[0][0] = temp*(1.0 - nu);
491:     C[1][1] = C[0][0];
492:     C[2][2] = C[0][0];
493:     C[3][3] = temp*(0.5 - nu);
494:     C[4][4] = C[3][3];
495:     C[5][5] = C[3][3];
496:     C[0][1] = temp*nu;
497:     C[0][2] = C[0][1];
498:     C[1][0] = C[0][1];
499:     C[1][2] = C[0][1];
500:     C[2][0] = C[0][1];
501:     C[2][1] = C[0][1];

503:     for (i=0; i<6; i++) {
504:       for (j=0; j<60; j++) {
505:         B_temp[i][j] = 0.0;
506:         for (k=0; k<6; k++) {
507:           B_temp[i][j] += C[i][k]*B[k][j];
508:         }
509:         B_temp[i][j] *= det_jac;
510:       }
511:     }

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