11 use hecmw,
only : kint, kreal
22 (etype, nn, ecoord, gausses, stiff, cdsys_id, coords, &
23 time, tincr, u, aux, temperature)
33 integer(kind=kint),
intent(in) :: etype
34 integer(kind=kint),
intent(in) :: nn
35 real(kind=kreal),
intent(in) :: ecoord(3, nn)
37 real(kind=kreal),
intent(out) :: stiff(:, :)
38 integer(kind=kint),
intent(in) :: cdsys_id
39 real(kind=kreal),
intent(inout) :: coords(3, 3)
40 real(kind=kreal),
intent(in) :: time
41 real(kind=kreal),
intent(in) :: tincr
42 real(kind=kreal),
intent(in),
optional :: u(3, nn)
43 real(kind=kreal),
intent(in),
optional :: aux(3, 3)
44 real(kind=kreal),
intent(in) :: temperature(nn)
48 integer(kind=kint) :: flag
49 integer(kind=kint),
parameter :: ndof = 3
50 real(kind=kreal) :: d(6, 6), b(6, ndof*(nn+3)), db(6, ndof*(nn+3))
51 real(kind=kreal) :: gderiv(nn+3, 3), stress(6)
52 real(kind=kreal) :: xj(9, 9), jacobian(3, 3), inverse(3, 3)
53 real(kind=kreal) :: tmpstiff((nn+3)*3, (nn+3)*3), tmpk(nn*3, 9)
54 real(kind=kreal) :: det, wg, elem(3, nn), mat(6, 6)
55 integer(kind=kint) :: i, j, lx
56 integer(kind=kint) :: serr
57 real(kind=kreal) :: naturalcoord(3), unode(3, nn+3)
58 real(kind=kreal) :: gdispderiv(3, 3), coordsys(3, 3)
59 real(kind=kreal) :: b1(6, ndof*(nn+3))
60 real(kind=kreal) :: smat(9, 9)
61 real(kind=kreal) :: bn(9, ndof*(nn+3)), sbn(9, ndof*(nn+3))
62 real(kind=kreal) :: spfunc(nn), temp
66 if(
present(u) .AND.
present(aux) )
then
67 unode(:, 1:nn) = u(:, :)
68 unode(:, nn+1:nn+3) = aux(:, :)
72 flag = gausses(1)%pMaterial%nlgeom_flag
73 if( .not.
present(u) ) flag = infinitesimal
74 elem(:, :) = ecoord(:, :)
75 if( flag == updatelag ) elem(:, :) = ecoord(:, :)+unode(:, 1:nn)
78 naturalcoord(:) = 0.0d0
79 call getjacobian(etype, nn, naturalcoord, elem, det, jacobian, inverse)
80 inverse(:, :)= inverse(:, :)*det
84 tmpstiff(:, :) = 0.0d0
85 b1(1:6, 1:(nn+3)*ndof) = 0.0d0
86 bn(1:9, 1:(nn+3)*ndof) = 0.0d0
91 call getglobalderiv( etype, nn, naturalcoord, elem, det, gderiv(1:nn, 1:3) )
93 if( cdsys_id > 0 )
then
95 if( serr == -1 ) stop
"Fail to setup local coordinate"
97 write(*, *)
"WARNING! Cannot setup local coordinate, it is modified automatically"
102 temp = dot_product( temperature, spfunc )
103 call matlmatrix( gausses(lx), d3, d, time, tincr, coordsys, temp )
105 if( flag == updatelag )
then
106 call geomat_c3( gausses(lx)%stress, mat )
107 d(:, :) = d(:, :)-mat
115 gderiv(nn+1, :) = -2.0d0*naturalcoord(1)*inverse(1, :)/det
116 gderiv(nn+2, :) = -2.0d0*naturalcoord(2)*inverse(2, :)/det
117 gderiv(nn+3, :) = -2.0d0*naturalcoord(3)*inverse(3, :)/det
120 b(1:6, 1:(nn+3)*ndof)=0.0d0
122 b(1, 3*j-2) = gderiv(j, 1)
123 b(2, 3*j-1) = gderiv(j, 2)
124 b(3, 3*j ) = gderiv(j, 3)
125 b(4, 3*j-2) = gderiv(j, 2)
126 b(4, 3*j-1) = gderiv(j, 1)
127 b(5, 3*j-1) = gderiv(j, 3)
128 b(5, 3*j ) = gderiv(j, 2)
129 b(6, 3*j-2) = gderiv(j, 3)
130 b(6, 3*j ) = gderiv(j, 1)
134 if( flag == totallag )
then
136 gdispderiv(1:ndof,1:ndof) = matmul( unode(1:ndof, 1:nn+3), gderiv(1:nn+3, 1:ndof) )
139 b1(1, 3*j-2) = gdispderiv(1, 1)*gderiv(j, 1)
140 b1(1, 3*j-1) = gdispderiv(2, 1)*gderiv(j, 1)
141 b1(1, 3*j ) = gdispderiv(3, 1)*gderiv(j, 1)
142 b1(2, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 2)
143 b1(2, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 2)
144 b1(2, 3*j ) = gdispderiv(3, 2)*gderiv(j, 2)
145 b1(3, 3*j-2) = gdispderiv(1, 3)*gderiv(j, 3)
146 b1(3, 3*j-1) = gdispderiv(2, 3)*gderiv(j, 3)
147 b1(3, 3*j ) = gdispderiv(3, 3)*gderiv(j, 3)
148 b1(4, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 1)+gdispderiv(1, 1)*gderiv(j, 2)
149 b1(4, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 1)+gdispderiv(2, 1)*gderiv(j, 2)
150 b1(4, 3*j ) = gdispderiv(3, 2)*gderiv(j, 1)+gdispderiv(3, 1)*gderiv(j, 2)
151 b1(5, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 3)+gdispderiv(1, 3)*gderiv(j, 2)
152 b1(5, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 3)+gdispderiv(2, 3)*gderiv(j, 2)
153 b1(5, 3*j ) = gdispderiv(3, 2)*gderiv(j, 3)+gdispderiv(3, 3)*gderiv(j, 2)
154 b1(6, 3*j-2) = gdispderiv(1, 3)*gderiv(j, 1)+gdispderiv(1, 1)*gderiv(j, 3)
155 b1(6, 3*j-1) = gdispderiv(2, 3)*gderiv(j, 1)+gdispderiv(2, 1)*gderiv(j, 3)
156 b1(6, 3*j ) = gdispderiv(3, 3)*gderiv(j, 1)+gdispderiv(3, 1)*gderiv(j, 3)
160 do j = 1, (nn+3)*ndof
161 b(:, j) = b(:, j)+b1(:, j)
166 db(1:6, 1:(nn+3)*ndof) = matmul( d, b(1:6, 1:(nn+3)*ndof) )
169 tmpstiff(i, j) = tmpstiff(i, j)+dot_product( b(:, i), db(:, j) )*wg
174 if( flag == totallag .OR. flag == updatelag )
then
175 stress(1:6) = gausses(lx)%stress
177 bn(1, 3*j-2) = gderiv(j, 1)
178 bn(2, 3*j-1) = gderiv(j, 1)
179 bn(3, 3*j ) = gderiv(j, 1)
180 bn(4, 3*j-2) = gderiv(j, 2)
181 bn(5, 3*j-1) = gderiv(j, 2)
182 bn(6, 3*j ) = gderiv(j, 2)
183 bn(7, 3*j-2) = gderiv(j, 3)
184 bn(8, 3*j-1) = gderiv(j, 3)
185 bn(9, 3*j ) = gderiv(j, 3)
189 smat(j , j ) = stress(1)
190 smat(j , j+3) = stress(4)
191 smat(j , j+6) = stress(6)
192 smat(j+3, j ) = stress(4)
193 smat(j+3, j+3) = stress(2)
194 smat(j+3, j+6) = stress(5)
195 smat(j+6, j ) = stress(6)
196 smat(j+6, j+3) = stress(5)
197 smat(j+6, j+6) = stress(3)
199 sbn(1:9, 1:(nn+3)*ndof) = matmul( smat(1:9, 1:9), bn(1:9, 1:(nn+3)*ndof) )
202 tmpstiff(i, j) = tmpstiff(i, j)+dot_product( bn(:, i), sbn(:, j) )*wg
210 xj(1:9, 1:9)= tmpstiff(nn*ndof+1:(nn+3)*ndof, nn*ndof+1:(nn+3)*ndof)
212 tmpk = matmul( tmpstiff( 1:nn*ndof, nn*ndof+1:(nn+3)*ndof ), xj )
213 stiff(1:nn*ndof, 1:nn*ndof) = tmpstiff(1:nn*ndof, 1:nn*ndof)-matmul( tmpk, tmpstiff(nn*ndof+1:(nn+3)*ndof, 1:nn*ndof) )
221 (etype,nn,ecoord, u, du, ddu, cdsys_id, coords, &
222 qf, gausses, iter, time, tincr, aux, ddaux, tt, t0, tn )
233 integer(kind=kint),
intent(in) :: etype
234 integer(kind=kint),
intent(in) :: nn
235 real(kind=kreal),
intent(in) :: ecoord(3, nn)
236 real(kind=kreal),
intent(in) :: u(3, nn)
237 real(kind=kreal),
intent(in) :: du(3, nn)
238 real(kind=kreal),
intent(in) :: ddu(3, nn)
239 integer(kind=kint),
intent(in) :: cdsys_id
240 real(kind=kreal),
intent(inout) :: coords(3, 3)
241 real(kind=kreal),
intent(out) :: qf(nn*3)
243 integer,
intent(in) :: iter
244 real(kind=kreal),
intent(in) :: time
245 real(kind=kreal),
intent(in) :: tincr
246 real(kind=kreal),
intent(inout) :: aux(3, 3)
247 real(kind=kreal),
intent(out) :: ddaux(3, 3)
248 real(kind=kreal),
intent(in) :: tt(nn)
249 real(kind=kreal),
intent(in) :: t0(nn)
250 real(kind=kreal),
intent(in) :: tn(nn)
253 integer(kind=kint) :: flag
254 integer(kind=kint),
parameter :: ndof = 3
255 real(kind=kreal) :: d(6,6), b(6,ndof*(nn+3)), b1(6,ndof*(nn+3)), spfunc(nn), ina(1)
256 real(kind=kreal) :: bn(9,ndof*(nn+3)), db(6, ndof*(nn+3)), stress(6), smat(9, 9), sbn(9, ndof*(nn+3))
257 real(kind=kreal) :: gderiv(nn+3,3), gderiv0(nn+3,3), gdispderiv(3,3), f(3,3), det, det0, wg, ttc, tt0, ttn
258 integer(kind=kint) :: i, j, lx, mtype, serr
259 real(kind=kreal) :: naturalcoord(3), rot(3,3), mat(6,6), epsth(6)
260 real(kind=kreal) :: totaldisp(3,nn+3), elem(3,nn), elem1(3,nn), coordsys(3,3)
261 real(kind=kreal) :: dstrain(6)
262 real(kind=kreal) :: alpo(3)
263 logical :: ierr, matlaniso
264 real(kind=kreal) :: stiff_ad(9, nn*3), stiff_aa(9, 9), qf_a(9)
265 real(kind=kreal) :: xj(9, 9)
266 real(kind=kreal) :: tmpforce(9), tmpdisp(9), tmp_d(nn*3), tmp_a(9)
267 real(kind=kreal) :: jacobian(3, 3), inverse(3, 3), inverse1(3, 3), inverse0(3, 3)
271 totaldisp(:, 1:nn) = u(:, :) + du(:, :) - ddu(:, :)
272 totaldisp(:, nn+1:nn+3) = aux(:, :)
275 flag = gausses(1)%pMaterial%nlgeom_flag
276 elem(:, :) = ecoord(:, :)
277 if( flag ==
updatelag ) elem(:, :) = ecoord(:, :)+totaldisp(:, 1:nn)
280 naturalcoord(:) = 0.0d0
281 call getjacobian(etype, nn, naturalcoord, elem, det, jacobian, inverse)
282 inverse(:, :)= inverse(:, :)*det
286 stiff_ad(:, :) = 0.0d0
287 stiff_aa(:, :) = 0.0d0
288 b1(1:6, 1:(nn+3)*ndof) = 0.0d0
289 bn(1:9, 1:(nn+3)*ndof) = 0.0d0
294 call getglobalderiv( etype, nn, naturalcoord, elem, det, gderiv(1:nn, 1:3) )
296 if( cdsys_id > 0 )
then
297 call set_localcoordsys( coords, g_localcoordsys(cdsys_id), coordsys(:, :), serr )
298 if( serr == -1 ) stop
"Fail to setup local coordinate"
299 if( serr == -2 )
then
300 write(*, *)
"WARNING! Cannot setup local coordinate, it is modified automatically"
305 ttc = dot_product( tt, spfunc )
306 call matlmatrix( gausses(lx),
d3, d, time, tincr, coordsys, ttc )
309 call geomat_c3( gausses(lx)%stress, mat )
310 d(:, :) = d(:, :)-mat
318 gderiv(nn+1, :) = -2.0d0*naturalcoord(1)*inverse(1, :)/det
319 gderiv(nn+2, :) = -2.0d0*naturalcoord(2)*inverse(2, :)/det
320 gderiv(nn+3, :) = -2.0d0*naturalcoord(3)*inverse(3, :)/det
323 b(1:6, 1:(nn+3)*ndof)=0.0d0
325 b(1, 3*j-2) = gderiv(j, 1)
326 b(2, 3*j-1) = gderiv(j, 2)
327 b(3, 3*j ) = gderiv(j, 3)
328 b(4, 3*j-2) = gderiv(j, 2)
329 b(4, 3*j-1) = gderiv(j, 1)
330 b(5, 3*j-1) = gderiv(j, 3)
331 b(5, 3*j ) = gderiv(j, 2)
332 b(6, 3*j-2) = gderiv(j, 3)
333 b(6, 3*j ) = gderiv(j, 1)
339 gdispderiv(1:ndof,1:ndof) = matmul( totaldisp(1:ndof, 1:nn+3), gderiv(1:nn+3, 1:ndof) )
342 b1(1, 3*j-2) = gdispderiv(1, 1)*gderiv(j, 1)
343 b1(1, 3*j-1) = gdispderiv(2, 1)*gderiv(j, 1)
344 b1(1, 3*j ) = gdispderiv(3, 1)*gderiv(j, 1)
345 b1(2, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 2)
346 b1(2, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 2)
347 b1(2, 3*j ) = gdispderiv(3, 2)*gderiv(j, 2)
348 b1(3, 3*j-2) = gdispderiv(1, 3)*gderiv(j, 3)
349 b1(3, 3*j-1) = gdispderiv(2, 3)*gderiv(j, 3)
350 b1(3, 3*j ) = gdispderiv(3, 3)*gderiv(j, 3)
351 b1(4, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 1)+gdispderiv(1, 1)*gderiv(j, 2)
352 b1(4, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 1)+gdispderiv(2, 1)*gderiv(j, 2)
353 b1(4, 3*j ) = gdispderiv(3, 2)*gderiv(j, 1)+gdispderiv(3, 1)*gderiv(j, 2)
354 b1(5, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 3)+gdispderiv(1, 3)*gderiv(j, 2)
355 b1(5, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 3)+gdispderiv(2, 3)*gderiv(j, 2)
356 b1(5, 3*j ) = gdispderiv(3, 2)*gderiv(j, 3)+gdispderiv(3, 3)*gderiv(j, 2)
357 b1(6, 3*j-2) = gdispderiv(1, 3)*gderiv(j, 1)+gdispderiv(1, 1)*gderiv(j, 3)
358 b1(6, 3*j-1) = gdispderiv(2, 3)*gderiv(j, 1)+gdispderiv(2, 1)*gderiv(j, 3)
359 b1(6, 3*j ) = gdispderiv(3, 3)*gderiv(j, 1)+gdispderiv(3, 1)*gderiv(j, 3)
363 do j = 1, (nn+3)*ndof
364 b(:, j) = b(:, j)+b1(:, j)
369 db(1:6, 1:(nn+3)*ndof) = matmul( d, b(1:6, 1:(nn+3)*ndof) )
372 stiff_ad(i, j) = stiff_ad(i, j)+dot_product( b(:, nn*ndof+i), db(:, j) )*wg
378 stiff_aa(i, j) = stiff_aa(i, j)+dot_product( b(:, nn*ndof+i), db(:, nn*ndof+j) )*wg
384 stress(1:6) = gausses(lx)%stress
386 bn(1, 3*j-2) = gderiv(j, 1)
387 bn(2, 3*j-1) = gderiv(j, 1)
388 bn(3, 3*j ) = gderiv(j, 1)
389 bn(4, 3*j-2) = gderiv(j, 2)
390 bn(5, 3*j-1) = gderiv(j, 2)
391 bn(6, 3*j ) = gderiv(j, 2)
392 bn(7, 3*j-2) = gderiv(j, 3)
393 bn(8, 3*j-1) = gderiv(j, 3)
394 bn(9, 3*j ) = gderiv(j, 3)
398 smat(j , j ) = stress(1)
399 smat(j , j+3) = stress(4)
400 smat(j , j+6) = stress(6)
401 smat(j+3, j ) = stress(4)
402 smat(j+3, j+3) = stress(2)
403 smat(j+3, j+6) = stress(5)
404 smat(j+6, j ) = stress(6)
405 smat(j+6, j+3) = stress(5)
406 smat(j+6, j+6) = stress(3)
408 sbn(1:9, 1:(nn+3)*ndof) = matmul( smat(1:9, 1:9), bn(1:9, 1:(nn+3)*ndof) )
412 stiff_ad(i, j) = stiff_ad(i, j)+dot_product( bn(:, nn*ndof+i), sbn(:, j) )*wg
417 stiff_aa(i, j) = stiff_aa(i, j)+dot_product( bn(:, nn*ndof+i), sbn(:, nn*ndof+j) )*wg
425 xj(1:3*ndof, 1:3*ndof)= stiff_aa(1:3*ndof, 1:3*ndof)
429 tmp_d((i-1)*ndof+1:i*ndof) = ddu(1:ndof, i)
431 tmpforce(1:3*ndof) = matmul( stiff_ad(1:3*ndof, 1:nn*ndof), tmp_d(1:nn*ndof) )
433 tmp_a(1:3*ndof) = -matmul( xj(1:3*ndof, 1:3*ndof), tmpforce(1:3*ndof) )
435 ddaux(1:ndof, i) = tmp_a((i-1)*ndof+1:i*ndof)
441 totaldisp(:,1:nn) = u(:,:)+ du(:,:)
442 totaldisp(:,nn+1:nn+3) = aux(:,:) + ddaux(:,:)
449 stiff_ad(:, :) = 0.0d0
450 stiff_aa(:, :) = 0.0d0
455 elem(:,:) = (0.5d0*du(:,:)+u(:,:) ) +ecoord(:,:)
456 elem1(:,:) = (du(:,:)+u(:,:) ) +ecoord(:,:)
458 totaldisp(:,1:nn) = du(:,:)
459 if( iter == 1 ) aux(:,:) = 0.0d0
460 totaldisp(:,nn+1:nn+3) = aux(:,:)
465 call fetch_tabledata(
mc_orthoexp, gausses(1)%pMaterial%dict, alpo(:), ierr, ina )
466 if( .not. ierr ) matlaniso = .true.
469 naturalcoord(:) = 0.0d0
470 call getjacobian(etype, nn, naturalcoord, elem, det, jacobian, inverse)
471 inverse(:, :)= inverse(:, :)*det
473 call getjacobian(etype, nn, naturalcoord, elem1, det, jacobian, inverse1)
474 inverse1(:, :)= inverse1(:, :)*det
475 call getjacobian(etype, nn, naturalcoord, ecoord, det, jacobian, inverse0)
476 inverse0(:, :)= inverse0(:, :)*det
481 mtype = gausses(lx)%pMaterial%mtype
486 if( cdsys_id > 0 )
then
487 call set_localcoordsys( coords, g_localcoordsys(cdsys_id), coordsys(:,:), serr )
488 if( serr == -1 ) stop
"Fail to setup local coordinate"
489 if( serr == -2 )
then
490 write(*, *)
"WARNING! Cannot setup local coordinate, it is modified automatically"
501 ttc = dot_product(tt, spfunc)
502 tt0 = dot_product(t0, spfunc)
503 ttn = dot_product(tn, spfunc)
511 gderiv(nn+1, :) = -2.0d0*naturalcoord(1)*inverse(1, :)/det
512 gderiv(nn+2, :) = -2.0d0*naturalcoord(2)*inverse(2, :)/det
513 gderiv(nn+3, :) = -2.0d0*naturalcoord(3)*inverse(3, :)/det
516 gdispderiv(1:ndof, 1:ndof) = matmul( totaldisp(1:ndof, 1:nn+3), gderiv(1:nn+3, 1:ndof) )
517 dstrain(1) = gdispderiv(1, 1)
518 dstrain(2) = gdispderiv(2, 2)
519 dstrain(3) = gdispderiv(3, 3)
520 dstrain(4) = ( gdispderiv(1, 2)+gdispderiv(2, 1) )
521 dstrain(5) = ( gdispderiv(2, 3)+gdispderiv(3, 2) )
522 dstrain(6) = ( gdispderiv(3, 1)+gdispderiv(1, 3) )
523 dstrain(:) = dstrain(:)-epsth(:)
525 f(1:3,1:3) = 0.d0; f(1,1)=1.d0; f(2,2)=1.d0; f(3,3)=1.d0;
527 gausses(lx)%strain(1:6) = dstrain(1:6)+epsth(:)
528 f(1:3,1:3) = f(1:3,1:3) + gdispderiv(1:3,1:3)
532 dstrain(1) = dstrain(1)+0.5d0*dot_product( gdispderiv(:, 1), gdispderiv(:, 1) )
533 dstrain(2) = dstrain(2)+0.5d0*dot_product( gdispderiv(:, 2), gdispderiv(:, 2) )
534 dstrain(3) = dstrain(3)+0.5d0*dot_product( gdispderiv(:, 3), gdispderiv(:, 3) )
535 dstrain(4) = dstrain(4)+( gdispderiv(1, 1)*gdispderiv(1, 2) &
536 +gdispderiv(2, 1)*gdispderiv(2, 2)+gdispderiv(3, 1)*gdispderiv(3, 2) )
537 dstrain(5) = dstrain(5)+( gdispderiv(1, 2)*gdispderiv(1, 3) &
538 +gdispderiv(2, 2)*gdispderiv(2, 3)+gdispderiv(3, 2)*gdispderiv(3, 3) )
539 dstrain(6) = dstrain(6)+( gdispderiv(1, 1)*gdispderiv(1, 3) &
540 +gdispderiv(2, 1)*gdispderiv(2, 3)+gdispderiv(3, 1)*gdispderiv(3, 3) )
542 gausses(lx)%strain(1:6) = dstrain(1:6)+epsth(:)
548 rot(1, 2)= 0.5d0*(gdispderiv(1, 2)-gdispderiv(2, 1) ); rot(2, 1) = -rot(1, 2)
549 rot(2, 3)= 0.5d0*(gdispderiv(2, 3)-gdispderiv(3, 2) ); rot(3, 2) = -rot(2, 3)
550 rot(1, 3)= 0.5d0*(gdispderiv(1, 3)-gdispderiv(3, 1) ); rot(3, 1) = -rot(1, 3)
552 gausses(lx)%strain(1:6) = gausses(lx)%strain_bak(1:6)+ dstrain(1:6)+epsth(:)
553 call getglobalderiv(etype, nn, naturalcoord, ecoord, det0, gderiv0)
554 gderiv0(nn+1, :) = -2.0d0*naturalcoord(1)*inverse0(1, :)/det0
555 gderiv0(nn+2, :) = -2.0d0*naturalcoord(2)*inverse0(2, :)/det0
556 gderiv0(nn+3, :) = -2.0d0*naturalcoord(3)*inverse0(3, :)/det0
557 f(1:3,1:3) = f(1:3,1:3) + matmul( totaldisp(1:ndof, 1:nn+3), gderiv0(1:nn+3, 1:ndof) )
562 call update_stress3d( flag, gausses(lx), rot, dstrain, f, coordsys, time, tincr, ttc, tt0, ttn )
568 b(1:6, 1:(nn+3)*ndof) = 0.0d0
570 b(1,3*j-2) = gderiv(j, 1)
571 b(2,3*j-1) = gderiv(j, 2)
572 b(3,3*j ) = gderiv(j, 3)
573 b(4,3*j-2) = gderiv(j, 2)
574 b(4,3*j-1) = gderiv(j, 1)
575 b(5,3*j-1) = gderiv(j, 3)
576 b(5,3*j ) = gderiv(j, 2)
577 b(6,3*j-2) = gderiv(j, 3)
578 b(6,3*j ) = gderiv(j, 1)
586 gdispderiv(1:ndof, 1:ndof) = matmul( totaldisp(1:ndof, 1:nn+3), gderiv(1:nn+3, 1:ndof) )
587 b1(1:6, 1:(nn+3)*ndof)=0.0d0
589 b1(1, 3*j-2) = gdispderiv(1, 1)*gderiv(j, 1)
590 b1(1, 3*j-1) = gdispderiv(2, 1)*gderiv(j, 1)
591 b1(1, 3*j ) = gdispderiv(3, 1)*gderiv(j, 1)
592 b1(2, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 2)
593 b1(2, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 2)
594 b1(2, 3*j ) = gdispderiv(3, 2)*gderiv(j, 2)
595 b1(3, 3*j-2) = gdispderiv(1, 3)*gderiv(j, 3)
596 b1(3, 3*j-1) = gdispderiv(2, 3)*gderiv(j, 3)
597 b1(3, 3*j ) = gdispderiv(3, 3)*gderiv(j, 3)
598 b1(4, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 1)+gdispderiv(1, 1)*gderiv(j, 2)
599 b1(4, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 1)+gdispderiv(2, 1)*gderiv(j, 2)
600 b1(4, 3*j ) = gdispderiv(3, 2)*gderiv(j, 1)+gdispderiv(3, 1)*gderiv(j, 2)
601 b1(5, 3*j-2) = gdispderiv(1, 2)*gderiv(j, 3)+gdispderiv(1, 3)*gderiv(j, 2)
602 b1(5, 3*j-1) = gdispderiv(2, 2)*gderiv(j, 3)+gdispderiv(2, 3)*gderiv(j, 2)
603 b1(5, 3*j ) = gdispderiv(3, 2)*gderiv(j, 3)+gdispderiv(3, 3)*gderiv(j, 2)
604 b1(6, 3*j-2) = gdispderiv(1, 3)*gderiv(j, 1)+gdispderiv(1, 1)*gderiv(j, 3)
605 b1(6, 3*j-1) = gdispderiv(2, 3)*gderiv(j, 1)+gdispderiv(2, 1)*gderiv(j, 3)
606 b1(6, 3*j ) = gdispderiv(3, 3)*gderiv(j, 1)+gdispderiv(3, 1)*gderiv(j, 3)
610 b(:,j) = b(:,j)+b1(:,j)
615 call getglobalderiv(etype, nn, naturalcoord, elem1, det, gderiv(1:nn,1:3))
622 gderiv(nn+1, :) = -2.0d0*naturalcoord(1)*inverse1(1, :)/det
623 gderiv(nn+2, :) = -2.0d0*naturalcoord(2)*inverse1(2, :)/det
624 gderiv(nn+3, :) = -2.0d0*naturalcoord(3)*inverse1(3, :)/det
626 b(1:6, 1:(nn+3)*ndof) = 0.0d0
628 b(1, 3*j-2) = gderiv(j, 1)
629 b(2, 3*j-1) = gderiv(j, 2)
630 b(3, 3*j ) = gderiv(j, 3)
631 b(4, 3*j-2) = gderiv(j, 2)
632 b(4, 3*j-1) = gderiv(j, 1)
633 b(5, 3*j-1) = gderiv(j, 3)
634 b(5, 3*j ) = gderiv(j, 2)
635 b(6, 3*j-2) = gderiv(j, 3)
636 b(6, 3*j ) = gderiv(j, 1)
643 db(1:6, 1:(nn+3)*ndof) = matmul( d, b(1:6, 1:(nn+3)*ndof) )
647 stiff_ad(i, j) = stiff_ad(i, j)+dot_product( b(:, nn*ndof+i), db(:, j) )*wg
652 stiff_aa(i, j) = stiff_aa(i, j)+dot_product( b(:, nn*ndof+i), db(:, nn*ndof+j) )*wg
658 = qf(1:nn*ndof)+matmul( gausses(lx)%stress(1:6), b(1:6,1:nn*ndof) )*wg
660 = qf_a(1:3*ndof)+matmul( gausses(lx)%stress(1:6), b(1:6,nn*ndof+1:(nn+3)*ndof) )*wg
663 gausses(lx)%strain_energy = gausses(lx)%strain_energy*wg
667 xj(1:9, 1:9)= stiff_aa(1:3*ndof, 1:3*ndof)
669 tmpdisp(:) = matmul( xj(:,:), qf_a(:) )
671 qf(i) = qf(i) - dot_product( stiff_ad(:,i), tmpdisp(:) )
679 (etype, nn, xx, yy, zz, tt, t0, &
680 gausses, vect, cdsys_id, coords)
690 integer(kind=kint),
parameter :: ndof=3
691 integer(kind=kint),
intent(in) :: etype
692 integer(kind=kint),
intent(in) :: nn
693 real(kind=kreal),
intent(in) :: xx(nn), yy(nn), zz(nn)
694 real(kind=kreal),
intent(in) :: tt(nn), t0(nn)
696 real(kind=kreal),
intent(out) :: vect(nn*ndof)
697 integer(kind=kint),
intent(in) :: cdsys_id
698 real(kind=kreal),
intent(inout) :: coords(3, 3)
702 real(kind=kreal) :: alp, alp0
703 real(kind=kreal) :: d(6, 6), b(6, ndof*(nn+3)), db_a(6, ndof*3)
704 real(kind=kreal) :: det, wg, ecoord(3, nn)
705 integer(kind=kint) :: i, j, ic, serr
706 real(kind=kreal) :: epsth(6), sgm(6), h(nn), alpo(3), alpo0(3), coordsys(3, 3), xj(9, 9)
707 real(kind=kreal) :: naturalcoord(3), gderiv(nn+3, 3)
708 real(kind=kreal) :: jacobian(3, 3),inverse(3, 3)
709 real(kind=kreal) :: stiff_da(nn*3, 3*3), stiff_aa(3*3, 3*3)
710 real(kind=kreal) :: vect_a(3*3)
711 real(kind=kreal) :: icdisp(9)
712 real(kind=kreal) :: tempc, temp0, thermal_eps, tm(6, 6), outa(1), ina(1)
713 logical :: ierr, matlaniso
718 if( cdsys_id > 0 )
then
720 call fetch_tabledata( mc_orthoexp, gausses(1)%pMaterial%dict, alpo(:), ierr, ina )
721 if( .not. ierr ) matlaniso = .true.
729 naturalcoord(:) = 0.0d0
730 call getjacobian(etype, nn, naturalcoord, ecoord, det, jacobian, inverse)
731 inverse(:, :) = inverse(:, :)*det
735 stiff_da(:, :) = 0.0d0
736 stiff_aa(:, :) = 0.0d0
737 b(1:6, 1:(nn+3)*ndof) = 0.0d0
738 vect(1:nn*ndof) = 0.0d0
739 vect_a(1:3*ndof) = 0.0d0
744 call getglobalderiv( etype, nn, naturalcoord, ecoord, det, gderiv(1:nn, 1:3) )
746 if( cdsys_id > 0 )
then
748 if( serr == -1 ) stop
"Fail to setup local coordinate"
749 if( serr == -2 )
then
750 write(*,*)
"WARNING! Cannot setup local coordinate, it is modified automatically"
755 tempc = dot_product( h(1:nn), tt(1:nn) )
756 temp0 = dot_product( h(1:nn), t0(1:nn) )
757 call matlmatrix( gausses(ic), d3, d, 1.d0, 1.0d0, coordsys, tempc )
764 gderiv(nn+1, :) = -2.0d0*naturalcoord(1)*inverse(1, :)/det
765 gderiv(nn+2, :) = -2.0d0*naturalcoord(2)*inverse(2, :)/det
766 gderiv(nn+3, :) = -2.0d0*naturalcoord(3)*inverse(3, :)/det
770 b(1, 3*j-2) = gderiv(j, 1)
771 b(2, 3*j-1) = gderiv(j, 2)
772 b(3, 3*j ) = gderiv(j, 3)
773 b(4, 3*j-2) = gderiv(j, 2)
774 b(4, 3*j-1) = gderiv(j, 1)
775 b(5, 3*j-1) = gderiv(j, 3)
776 b(5, 3*j ) = gderiv(j, 2)
777 b(6, 3*j-2) = gderiv(j, 3)
778 b(6, 3*j ) = gderiv(j, 1)
781 db_a(1:6, 1:3*ndof) = matmul( d, b(1:6, nn*ndof+1:(nn+3)*ndof) )
784 stiff_da(i, j) = stiff_da(i, j) + dot_product( b(1:6, i), db_a(1:6, j) ) * wg
789 stiff_aa(i, j) = stiff_aa(i, j) + dot_product( b(1:6, nn*ndof+i), db_a(1:6, j) ) * wg
795 call fetch_tabledata( mc_orthoexp, gausses(ic)%pMaterial%dict, alpo(:), ierr, ina )
796 if( ierr ) stop
"Fails in fetching orthotropic expansion coefficient!"
798 call fetch_tabledata( mc_themoexp, gausses(ic)%pMaterial%dict, outa(:), ierr, ina )
799 if( ierr ) outa(1) = gausses(ic)%pMaterial%variables(m_exapnsion)
804 call fetch_tabledata( mc_orthoexp, gausses(ic)%pMaterial%dict, alpo0(:), ierr, ina )
805 if( ierr ) stop
"Fails in fetching orthotropic expansion coefficient!"
807 call fetch_tabledata( mc_themoexp, gausses(ic)%pMaterial%dict, outa(:), ierr, ina )
808 if( ierr ) outa(1) = gausses(ic)%pMaterial%variables(m_exapnsion)
821 epsth(:) = matmul( epsth(:), tm )
822 epsth(4:6) = epsth(4:6)*2.0d0
825 epsth(1:3) = thermal_eps
832 sgm(:) = matmul( d(:, :), epsth(:) )
837 vect(1:nn*ndof) = vect(1:nn*ndof)+matmul( sgm(1:6), b(1:6, 1:nn*ndof) )*wg
838 vect_a(1:3*ndof) = vect_a(1:3*ndof)+matmul( sgm(1:6), b(1:6, nn*ndof+1:(nn+3)*ndof) )*wg
843 xj(1:9,1:9)= stiff_aa(1:9, 1:9)
846 icdisp(1:9) = matmul( xj(:, :), vect_a(1:3*ndof) )
848 vect(1:nn*ndof) = vect(1:nn*ndof) - matmul( stiff_da(1:nn*ndof, 1:3*ndof), icdisp(1:3*ndof) )