10 private :: nodalstress_inv3, nodalstress_inv2, inverse_func
18 type(hecmwst_local_mesh) :: hecMESH
20 real(kind=kreal),
pointer :: tnstrain(:), testrain(:), yield_ratio(:)
21 integer(kind=kint),
pointer :: is_rot(:)
23 integer(kind=kint) :: itype, icel, ic, is, iE, jS, i, j, k, m, ic_type, nn, ni, ID_area
24 integer(kind=kint) :: nodlocal(20), ntemp
25 integer(kind=kint),
allocatable :: nnumber(:)
26 real(kind=kreal) :: estrain(6), estress(6), naturalcoord(3)
27 real(kind=kreal) :: enqm(12)
28 real(kind=kreal) :: ndstrain(20,6), ndstress(20,6), tdstrain(20,6)
29 real(kind=kreal) :: ecoord(3, 20), edisp(60), tt(20), t0(20)
30 real(kind=kreal),
allocatable :: func(:,:), inv_func(:,:)
33 integer(kind=kint) :: isect, ihead, ntot_lyr, nlyr, flag33, cid, truss
34 real(kind=kreal) :: thick, thick_lyr, dtot_lyr
37 allocate( nnumber(hecmesh%n_node) )
38 if( .not.
associated(fstrsolid%is_rot) )
allocate( fstrsolid%is_rot(hecmesh%n_node) )
44 tnstrain => fstrsolid%tnstrain
45 testrain => fstrsolid%testrain
46 is_rot => fstrsolid%is_rot
47 yield_ratio => fstrsolid%yield_ratio
49 if(
associated(tnstrain) ) tnstrain = 0.0d0
52 ntot_lyr = fstrsolid%max_lyr
53 flag33 = fstrsolid%is_33shell
54 truss = fstrsolid%is_33beam
59 do itype = 1, hecmesh%n_elem_type
60 is = hecmesh%elem_type_index(itype-1) + 1
61 ie = hecmesh%elem_type_index(itype )
62 ic_type = hecmesh%elem_type_item(itype)
63 if( ic_type == fe_tet10nc ) ic_type = fe_tet10n
64 if( .not. (hecmw_is_etype_solid(ic_type) .or. ic_type == 781 &
65 & .or. ic_type == 761 .or. ic_type == fe_beam341 ) ) cycle
67 nn = hecmw_get_max_node( ic_type )
68 ni = numofquadpoints( ic_type )
69 allocate( func(ni,nn), inv_func(nn,ni) )
70 if( ic_type == fe_tet10n )
then
71 ic = hecmw_get_max_node( fe_tet4n )
73 call getquadpoint( ic_type, i, naturalcoord )
74 call getshapefunc( fe_tet4n, naturalcoord, func(i,1:ic) )
76 call inverse_func( ic, func, inv_func )
77 else if( ic_type == fe_hex8n )
then
79 call getquadpoint( ic_type, i, naturalcoord )
80 call getshapefunc( ic_type, naturalcoord, func(i,1:nn) )
82 call inverse_func( ni, func, inv_func )
83 else if( ic_type == fe_prism15n )
then
86 if( i==1 .or. i==2 .or. i==3 .or. i==7 .or. i==8 .or. i==9 )
then
88 call getquadpoint( ic_type, i, naturalcoord )
89 call getshapefunc( fe_prism6n, naturalcoord, func(ic,1:6) )
92 call inverse_func( ic, func, inv_func )
94 else if( ic_type == fe_hex20n )
then
97 if( i==1 .or. i==3 .or. i==7 .or. i==9 .or. &
98 i==19 .or. i==21 .or. i==25 .or. i==27 )
then
100 call getquadpoint( ic_type, i, naturalcoord )
101 call getshapefunc( fe_hex8n, naturalcoord, func(ic,1:8) )
104 call inverse_func( ic, func, inv_func )
109 js = hecmesh%elem_node_index(icel-1)
110 id_area = hecmesh%elem_ID(icel*2)
111 isect= hecmesh%section_ID(icel)
112 ihead = hecmesh%section%sect_R_index(isect-1)
113 thick = hecmesh%section%sect_R_item(ihead+1)
123 if( ic_type == 641 )
then
125 nodlocal(j) = hecmesh%elem_node_item(js+j)
126 ecoord(1:3,j) = hecmesh%node(3*nodlocal(j)-2:3*nodlocal(j))
127 edisp(3*j-2:3*j) = fstrsolid%unode(3*nodlocal(j)-2:3*nodlocal(j))
130 if(
associated( fstrsolid%temperature ) )
then
133 nodlocal(j) = hecmesh%elem_node_item(js+j)
134 t0(j) = fstrsolid%last_temp( nodlocal(j) )
135 tt(j) = fstrsolid%temperature( nodlocal(j) )
138 call nodalstress_beam_641( ic_type, nn, ecoord, fstrsolid%elements(icel)%gausses, &
139 & hecmesh%section%sect_R_item(ihead+1:), edisp, &
140 & ndstrain(1:nn,1:6), ndstress(1:nn,1:6), tt(1:nn), t0(1:nn), ntemp )
141 call elementalstress_beam_641( fstrsolid%elements(icel)%gausses, estrain, estress, enqm )
142 fstrsolid%ENQM(icel*12-11:icel*12) = enqm(1:12)
145 elseif( ic_type == 781)
then
147 nodlocal(j ) = hecmesh%elem_node_item(js+j )
148 nodlocal(j+4) = hecmesh%elem_node_item(js+j+4)
149 is_rot(nodlocal(j+4)) = 1
150 ecoord(1:3,j ) = hecmesh%node(3*nodlocal(j )-2:3*nodlocal(j ))
151 ecoord(1:3,j+4) = hecmesh%node(3*nodlocal(j+4)-2:3*nodlocal(j+4))
152 edisp(6*j-5:6*j-3) = fstrsolid%unode(3*nodlocal(j )-2:3*nodlocal(j ))
153 edisp(6*j-2:6*j ) = fstrsolid%unode(3*nodlocal(j+4)-2:3*nodlocal(j+4))
155 ntot_lyr = fstrsolid%elements(icel)%gausses(1)%pMaterial%totallyr
157 call elementstress_shell_mitc( 741, 4, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
158 & ndstrain(1:4,1:6), ndstress(1:4,1:6), thick, 1.0d0, nlyr, ntot_lyr)
161 call elementstress_shell_mitc( 741, 4, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
162 & ndstrain(1:4,1:6), ndstress(1:4,1:6), thick,-1.0d0, nlyr, ntot_lyr)
165 call fstr_getavg_shell(4,fstrsolid,icel,nodlocal,ndstrain(1:4,1:6),ndstress(1:4,1:6),estrain,estress)
167 elseif( ic_type == 761)
then
169 nodlocal(j ) = hecmesh%elem_node_item(js+j )
170 nodlocal(j+3) = hecmesh%elem_node_item(js+j+3)
171 is_rot(nodlocal(j+3)) = 1
172 ecoord(1:3,j ) = hecmesh%node(3*nodlocal(j )-2:3*nodlocal(j ))
173 ecoord(1:3,j+3) = hecmesh%node(3*nodlocal(j+3)-2:3*nodlocal(j+3))
174 edisp(6*j-5:6*j-3) = fstrsolid%unode(3*nodlocal(j )-2:3*nodlocal(j ))
175 edisp(6*j-2:6*j ) = fstrsolid%unode(3*nodlocal(j+3)-2:3*nodlocal(j+3))
177 ntot_lyr = fstrsolid%elements(icel)%gausses(1)%pMaterial%totallyr
179 call elementstress_shell_mitc( 731, 3, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
180 & ndstrain(1:3,1:6), ndstress(1:3,1:6), thick, 1.0d0, nlyr, ntot_lyr)
183 call elementstress_shell_mitc( 731, 3, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
184 & ndstrain(1:3,1:6), ndstress(1:3,1:6), thick,-1.0d0, nlyr, ntot_lyr)
187 call fstr_getavg_shell(3,fstrsolid,icel,nodlocal,ndstrain(1:3,1:6),ndstress(1:3,1:6),estrain,estress)
189 else if( ic_type == 301 )
then
190 call nodalstress_c1( ic_type, nn, fstrsolid%elements(icel)%gausses, &
191 ndstrain(1:nn,1:6), ndstress(1:nn,1:6) )
192 call elementstress_c1( ic_type, fstrsolid%elements(icel)%gausses, estrain, estress )
194 else if( ic_type == fe_tet10n .or. ic_type == fe_hex8n .or. &
195 ic_type == fe_prism15n .or. ic_type == fe_hex20n )
then
196 call nodalstress_inv3( ic_type, ni, fstrsolid%elements(icel)%gausses, &
197 inv_func, ndstrain(1:nn,1:6), ndstress(1:nn,1:6), &
199 call elementstress_c3( ic_type, fstrsolid%elements(icel)%gausses, estrain, estress )
201 else if ( ic_type == 881 .or. ic_type == 891 )
then
204 if( ic_type == 341 .and. fstrsolid%sections(isect)%elemopt341 ==
kel341sesns ) cycle
206 call nodalstress_c3( ic_type, nn, fstrsolid%elements(icel)%gausses, &
207 ndstrain(1:nn,1:6), ndstress(1:nn,1:6) )
210 call elementstress_c3( ic_type, fstrsolid%elements(icel)%gausses, estrain, estress )
216 ic = hecmesh%elem_node_item(js+j)
217 fstrsolid%STRAIN(6*(ic-1)+1:6*(ic-1)+6) = fstrsolid%STRAIN(6*(ic-1)+1:6*(ic-1)+6) + ndstrain(j,1:6)
218 fstrsolid%STRESS(6*(ic-1)+1:6*(ic-1)+6) = fstrsolid%STRESS(6*(ic-1)+1:6*(ic-1)+6) + ndstress(j,1:6)
219 if(
associated(tnstrain) )
then
220 tnstrain(6*(ic-1)+1:6*(ic-1)+6) = tnstrain(6*(ic-1)+1:6*(ic-1)+6) + tdstrain(j,1:6)
222 nnumber(ic) = nnumber(ic) + 1
225 fstrsolid%ESTRAIN(6*(icel-1)+1:6*(icel-1)+6) = fstrsolid%ESTRAIN(6*(icel-1)+1:6*(icel-1)+6) + estrain(1:6)
226 fstrsolid%ESTRESS(6*(icel-1)+1:6*(icel-1)+6) = fstrsolid%ESTRESS(6*(icel-1)+1:6*(icel-1)+6) + estress(1:6)
230 deallocate( func, inv_func )
234 do i = 1, hecmesh%n_node
235 if( nnumber(i) == 0 ) cycle
236 fstrsolid%STRAIN(6*(i-1)+1:6*(i-1)+6) = fstrsolid%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
237 fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6) = fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
238 if(
associated(tnstrain) )
then
239 tnstrain(6*(i-1)+1:6*(i-1)+6) = tnstrain(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
244 & hecmesh, fstrsolid, nnumber, fstrsolid%STRAIN, fstrsolid%STRESS, fstrsolid%ESTRAIN, fstrsolid%ESTRESS )
246 if( flag33 == 1 )
then
247 do nlyr = 1, ntot_lyr
248 do i = 1, hecmesh%n_node
249 if( nnumber(i) == 0 ) cycle
250 fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRAIN(6*(i-1)+1:6*(i-1)+6) = &
251 & fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
252 fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6) = &
253 & fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
254 fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRAIN(6*(i-1)+1:6*(i-1)+6) = &
255 & fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
256 fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6) = &
257 & fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
263 do i = 1, hecmesh%n_node
264 fstrsolid%MISES(i) =
get_mises(fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6))
266 do i = 1, hecmesh%n_elem
267 fstrsolid%EMISES(i) =
get_mises(fstrsolid%ESTRESS(6*(i-1)+1:6*(i-1)+6))
270 if( flag33 == 1 )
then
271 if( fstrsolid%output_ctrl(3)%outinfo%on(27) .or. fstrsolid%output_ctrl(4)%outinfo%on(27) )
then
272 do nlyr = 1, ntot_lyr
273 call make_principal(fstrsolid, hecmesh, fstrsolid%SHELL%LAYER(nlyr)%PLUS)
274 call make_principal(fstrsolid, hecmesh, fstrsolid%SHELL%LAYER(nlyr)%MINUS)
282 deallocate( nnumber )
287 integer(kind=kint),
allocatable,
intent(in) :: irow(:)
288 integer(kind=kint),
allocatable,
intent(in) :: asect(:)
289 integer(kind=kint),
intent(in) :: nid
290 integer(kind=kint),
intent(in) :: sid
292 integer(kind=kint) :: i
295 do i=irow(nid-1)+1,irow(nid)
296 if( asect(i) == sid )
then
305 Nodal_STRAIN, Nodal_STRESS, Elemental_STRAIN, Elemental_STRESS )
306 type(hecmwst_local_mesh),
intent(in) :: hecMESH
308 integer(kind=kint),
allocatable,
intent(inout) :: nnumber(:)
309 real(kind=kreal),
pointer,
intent(inout) :: nodal_strain(:)
310 real(kind=kreal),
pointer,
intent(inout) :: nodal_stress(:)
311 real(kind=kreal),
pointer,
intent(inout) :: elemental_strain(:)
312 real(kind=kreal),
pointer,
intent(inout) :: elemental_stress(:)
314 integer(kind=kint) :: itype, iS, iE, jS, ic_type, icel, i, j, isect
315 integer(kind=kint) :: nsize, nid(2), idx(2), nd
316 integer(kind=kint) :: nnode, nlen
317 type(hecmwst_varray_int),
allocatable :: nodal_sections(:)
318 real(kind=kreal) :: tmpval(6), hydval, nsecdup
319 integer(kind=kint),
allocatable :: irow(:), jcol(:), asect(:)
320 real(kind=kreal),
allocatable :: stress_hyd(:), strain_hyd(:)
321 real(kind=kreal),
allocatable :: stress_dev(:)
322 real(kind=kreal) :: stress_hyd_ndave(6), strain_hyd_ndave(6)
323 real(kind=kreal) :: stress_dev_ndave(6), strain_dev_ndave(6)
324 real(kind=kreal),
allocatable :: n_dup_dev(:), n_dup_hyd(:)
325 real(kind=kreal) :: edstrain(6), edstress(6)
327 nnode = hecmesh%n_node
328 nsize =
size(nodal_strain)
331 call hecmw_varray_int_initialize_all( nodal_sections, nnode, 2 )
332 do itype = 1, hecmesh%n_elem_type
333 ic_type = hecmesh%elem_type_item(itype)
334 if( ic_type /= 341 ) cycle
336 is = hecmesh%elem_type_index(itype-1) + 1
337 ie = hecmesh%elem_type_index(itype )
340 isect= hecmesh%section_ID(icel)
341 if( fstrsolid%sections(isect)%elemopt341 /=
kel341sesns ) cycle
342 js = hecmesh%elem_node_index(icel-1)
344 nd = hecmesh%elem_node_item(js+i)
345 call hecmw_varray_int_add_if_not_exits( nodal_sections(nd), isect )
351 allocate(irow(0:nnode))
354 irow(i) = irow(i-1)+hecmw_varray_int_get_nitem(nodal_sections(i))
358 allocate(asect(nlen))
360 if( irow(i-1) == irow(i) ) cycle
361 call hecmw_varray_int_get_item_all( nodal_sections(i), asect(irow(i-1)+1:irow(i)) )
365 allocate(stress_hyd(6*nlen), strain_hyd(6*nlen))
366 allocate(stress_dev(6*nlen))
367 allocate(n_dup_dev(nlen),n_dup_hyd(nlen))
374 do itype = 1, hecmesh%n_elem_type
375 ic_type = hecmesh%elem_type_item(itype)
376 if( ic_type /= 881 .and. ic_type /= 891 ) cycle
378 is = hecmesh%elem_type_index(itype-1) + 1
379 ie = hecmesh%elem_type_index(itype )
382 js = hecmesh%elem_node_index(icel-1)
383 isect= hecmesh%section_ID(icel)
384 if( ic_type == 881 )
then
385 nid(1) = hecmesh%elem_node_item(js+1)
389 strain_hyd(6*idx(1)-5:6*idx(1)) = fstrsolid%elements(icel)%gausses(1)%strain_out(1:6)
391 stress_hyd(6*idx(1)-5:6*idx(1)) = fstrsolid%elements(icel)%gausses(1)%stress_out(1:6)
393 n_dup_hyd(idx(1)) = n_dup_hyd(idx(1)) + 1.d0
394 else if( ic_type == 891 )
then
395 nid(1:2) = hecmesh%elem_node_item(js+1:js+2)
400 tmpval(1:6) = fstrsolid%elements(icel)%gausses(1)%stress_out(1:6)
401 stress_dev(6*idx(1)-5:6*idx(1)) = stress_dev(6*idx(1)-5:6*idx(1)) + tmpval(1:6)
402 stress_dev(6*idx(2)-5:6*idx(2)) = stress_dev(6*idx(2)-5:6*idx(2)) + tmpval(1:6)
404 n_dup_dev(idx(1)) = n_dup_dev(idx(1)) + 1.d0
405 n_dup_dev(idx(2)) = n_dup_dev(idx(2)) + 1.d0
411 if( irow(i-1) == irow(i) ) cycle
412 do j=irow(i-1)+1,irow(i)
413 if( n_dup_dev(j) < 1.0d-8 ) cycle
414 stress_dev(6*j-5:6*j) = stress_dev(6*j-5:6*j)/n_dup_dev(j)
420 if( irow(i-1) == irow(i) ) cycle
421 strain_hyd_ndave(:) = 0.d0
422 stress_hyd_ndave(:) = 0.d0
423 stress_dev_ndave(:) = 0.d0
424 do j=irow(i-1)+1,irow(i)
425 strain_hyd_ndave(1:6) = strain_hyd_ndave(1:6) + strain_hyd(6*j-5:6*j)
426 stress_hyd_ndave(1:6) = stress_hyd_ndave(1:6) + stress_hyd(6*j-5:6*j)
427 stress_dev_ndave(1:6) = stress_dev_ndave(1:6) + stress_dev(6*j-5:6*j)
429 nsecdup = dble(irow(i)-irow(i-1))
430 strain_hyd_ndave(1:6) = strain_hyd_ndave(1:6)/nsecdup
431 stress_hyd_ndave(1:6) = stress_hyd_ndave(1:6)/nsecdup
432 stress_dev_ndave(1:6) = stress_dev_ndave(1:6)/nsecdup
434 if( nnumber(i) == 0 )
then
435 nodal_strain(6*i-5:6*i) = strain_hyd_ndave(1:6)
436 nodal_stress(6*i-5:6*i) = stress_hyd_ndave(1:6)+stress_dev_ndave(1:6)
438 nodal_strain(6*i-5:6*i) = 0.5d0*(nodal_strain(6*i-5:6*i)+strain_hyd_ndave(1:6))
439 nodal_stress(6*i-5:6*i) = 0.5d0*(nodal_stress(6*i-5:6*i)+stress_hyd_ndave(1:6)+stress_dev_ndave(1:6))
444 do itype = 1, hecmesh%n_elem_type
445 ic_type = hecmesh%elem_type_item(itype)
446 if( ic_type /= 341 ) cycle
448 is = hecmesh%elem_type_index(itype-1) + 1
449 ie = hecmesh%elem_type_index(itype )
452 isect= hecmesh%section_ID(icel)
453 if( fstrsolid%sections(isect)%elemopt341 /=
kel341sesns ) cycle
454 js = hecmesh%elem_node_index(icel-1)
458 nd = hecmesh%elem_node_item(js+i)
459 idx(1) =
search_idx_senes( irow, asect, hecmesh%elem_node_item(js+i), isect )
460 edstrain(1:6) = edstrain(1:6) + strain_hyd(6*idx(1)-5:6*idx(1))
461 edstress(1:6) = edstress(1:6) + stress_hyd(6*idx(1)-5:6*idx(1)) + stress_dev(6*idx(1)-5:6*idx(1))
463 edstrain(1:6) = 0.25d0*edstrain(1:6)
464 edstress(1:6) = 0.25d0*edstress(1:6)
466 elemental_strain(6*(icel-1)+1:6*(icel-1)+6) = elemental_strain(6*(icel-1)+1:6*(icel-1)+6) + edstrain(1:6)
467 elemental_stress(6*(icel-1)+1:6*(icel-1)+6) = elemental_stress(6*(icel-1)+1:6*(icel-1)+6) + edstress(1:6)
469 fstrsolid%elements(icel)%gausses(1)%strain_out(1:6) = elemental_strain(6*(icel-1)+1:6*(icel-1)+6)
470 fstrsolid%elements(icel)%gausses(1)%stress_out(1:6) = elemental_stress(6*(icel-1)+1:6*(icel-1)+6)
474 deallocate(stress_hyd, strain_hyd)
475 deallocate(stress_dev)
476 deallocate(n_dup_dev, n_dup_hyd)
483 integer(kind=kint) :: nodlocal(20)
484 integer(kind=kint) :: nn, i, j, k, m, nlyr, weight, icel, flag
485 real(kind=kreal) :: strain(nn, 6), stress(nn, 6)
492 layer => fstrsolid%SHELL%LAYER(nlyr)%PLUS
493 elseif(flag == -1)
then
494 layer => fstrsolid%SHELL%LAYER(nlyr)%MINUS
497 layer%STRAIN(6*(i-1)+k) = layer%STRAIN(6*(i-1)+k) + strain(j,k)
498 layer%STRAIN(6*(m-1)+k) = layer%STRAIN(6*(m-1)+k) + strain(j,k)
499 layer%STRESS(6*(i-1)+k) = layer%STRESS(6*(i-1)+k) + stress(j,k)
500 layer%STRESS(6*(m-1)+k) = layer%STRESS(6*(m-1)+k) + stress(j,k)
501 layer%ESTRAIN(6*(icel-1)+k) = layer%ESTRAIN(6*(icel-1)+k) + strain(j,k)/nn
502 layer%ESTRESS(6*(icel-1)+k) = layer%ESTRESS(6*(icel-1)+k) + stress(j,k)/nn
507 subroutine fstr_getavg_shell(nn,fstrSOLID,icel,nodLOCAL,strain,stress,estrain,estress)
510 integer(kind=kint) :: nodlocal(20)
511 integer(kind=kint) :: nn, i, j, k, m, nlyr, icel, flag, ntot_lyr
512 real(kind=kreal) :: strain(nn,6), stress(nn,6), estrain(6), estress(6), weight
515 ntot_lyr = fstrsolid%elements(icel)%gausses(1)%pMaterial%totallyr
521 do nlyr = 1, ntot_lyr
522 layer => fstrsolid%SHELL%LAYER(nlyr)
523 weight = fstrsolid%elements(icel)%gausses(1)%pMaterial%shell_var(nlyr)%weight
527 strain(j,k) = strain(j,k) &
528 & + weight*(0.5d0*layer%PLUS%STRAIN(6*(i-1)+k) + 0.5d0*layer%MINUS%STRAIN(6*(i-1)+k))
529 stress(j,k) = stress(j,k) &
530 & + weight*(0.5d0*layer%PLUS%STRESS(6*(i-1)+k) + 0.5d0*layer%MINUS%STRESS(6*(i-1)+k))
532 estrain(j) = estrain(j) &
533 & + weight*(0.5d0*layer%PLUS%ESTRAIN(6*(icel-1)+j) + 0.5d0*layer%MINUS%ESTRAIN(6*(icel-1)+j))
534 estress(j) = estress(j) &
535 & + weight*(0.5d0*layer%PLUS%ESTRESS(6*(icel-1)+j) + 0.5d0*layer%MINUS%ESTRESS(6*(icel-1)+j))
541 subroutine nodalstress_inv3( etype, ni, gausses, func, edstrain, edstress, tdstrain )
544 integer(kind=kint) :: etype, ni
546 real(kind=kreal) :: func(:, :), edstrain(:, :), edstress(:, :), tdstrain(:, :)
547 integer :: i, j, k, ic
553 if( etype == fe_hex8n )
then
557 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
558 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
563 else if( etype == fe_tet10n )
then
567 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
568 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
573 edstrain(5,1:6) = ( edstrain(1,1:6) + edstrain(2,1:6) ) / 2.0
574 edstress(5,1:6) = ( edstress(1,1:6) + edstress(2,1:6) ) / 2.0
575 tdstrain(5,1:6) = ( tdstrain(1,1:6) + tdstrain(2,1:6) ) / 2.0
576 edstrain(6,1:6) = ( edstrain(2,1:6) + edstrain(3,1:6) ) / 2.0
577 edstress(6,1:6) = ( edstress(2,1:6) + edstress(3,1:6) ) / 2.0
578 tdstrain(6,1:6) = ( tdstrain(2,1:6) + tdstrain(3,1:6) ) / 2.0
579 edstrain(7,1:6) = ( edstrain(3,1:6) + edstrain(1,1:6) ) / 2.0
580 edstress(7,1:6) = ( edstress(3,1:6) + edstress(1,1:6) ) / 2.0
581 tdstrain(7,1:6) = ( tdstrain(3,1:6) + tdstrain(1,1:6) ) / 2.0
582 edstrain(8,1:6) = ( edstrain(1,1:6) + edstrain(4,1:6) ) / 2.0
583 edstress(8,1:6) = ( edstress(1,1:6) + edstress(4,1:6) ) / 2.0
584 tdstrain(8,1:6) = ( tdstrain(1,1:6) + tdstrain(4,1:6) ) / 2.0
585 edstrain(9,1:6) = ( edstrain(2,1:6) + edstrain(4,1:6) ) / 2.0
586 edstress(9,1:6) = ( edstress(2,1:6) + edstress(4,1:6) ) / 2.0
587 tdstrain(9,1:6) = ( tdstrain(2,1:6) + tdstrain(4,1:6) ) / 2.0
588 edstrain(10,1:6) = ( edstrain(3,1:6) + edstrain(4,1:6) ) / 2.0
589 edstress(10,1:6) = ( edstress(3,1:6) + edstress(4,1:6) ) / 2.0
590 tdstrain(10,1:6) = ( tdstrain(3,1:6) + tdstrain(4,1:6) ) / 2.0
591 else if( etype == fe_prism15n )
then
594 do j = 1, numofquadpoints(etype)
595 if( j==1 .or. j==2 .or. j==3 .or. j==7 .or. j==8 .or. j==9 )
then
598 edstrain(i,k) = edstrain(i,k) + func(i,ic) * gausses(j)%strain_out(k)
599 edstress(i,k) = edstress(i,k) + func(i,ic) * gausses(j)%stress_out(k)
605 edstrain(7,1:6) = ( edstrain(1,1:6) + edstrain(2,1:6) ) / 2.0
606 edstress(7,1:6) = ( edstress(1,1:6) + edstress(2,1:6) ) / 2.0
607 tdstrain(7,1:6) = ( tdstrain(1,1:6) + tdstrain(2,1:6) ) / 2.0
608 edstrain(8,1:6) = ( edstrain(2,1:6) + edstrain(3,1:6) ) / 2.0
609 edstress(8,1:6) = ( edstress(2,1:6) + edstress(3,1:6) ) / 2.0
610 tdstrain(8,1:6) = ( tdstrain(2,1:6) + tdstrain(3,1:6) ) / 2.0
611 edstrain(9,1:6) = ( edstrain(3,1:6) + edstrain(1,1:6) ) / 2.0
612 edstress(9,1:6) = ( edstress(3,1:6) + edstress(1,1:6) ) / 2.0
613 tdstrain(9,1:6) = ( tdstrain(3,1:6) + tdstrain(1,1:6) ) / 2.0
614 edstrain(10,1:6) = ( edstrain(4,1:6) + edstrain(5,1:6) ) / 2.0
615 edstress(10,1:6) = ( edstress(4,1:6) + edstress(5,1:6) ) / 2.0
616 tdstrain(10,1:6) = ( tdstrain(4,1:6) + tdstrain(5,1:6) ) / 2.0
617 edstrain(11,1:6) = ( edstrain(5,1:6) + edstrain(6,1:6) ) / 2.0
618 edstress(11,1:6) = ( edstress(5,1:6) + edstress(6,1:6) ) / 2.0
619 tdstrain(11,1:6) = ( tdstrain(5,1:6) + tdstrain(6,1:6) ) / 2.0
620 edstrain(12,1:6) = ( edstrain(6,1:6) + edstrain(4,1:6) ) / 2.0
621 edstress(12,1:6) = ( edstress(6,1:6) + edstress(4,1:6) ) / 2.0
622 tdstrain(12,1:6) = ( tdstrain(6,1:6) + tdstrain(4,1:6) ) / 2.0
623 edstrain(13,1:6) = ( edstrain(1,1:6) + edstrain(4,1:6) ) / 2.0
624 edstress(13,1:6) = ( edstress(1,1:6) + edstress(4,1:6) ) / 2.0
625 tdstrain(13,1:6) = ( tdstrain(1,1:6) + tdstrain(4,1:6) ) / 2.0
626 edstrain(14,1:6) = ( edstrain(2,1:6) + edstrain(5,1:6) ) / 2.0
627 edstress(14,1:6) = ( edstress(2,1:6) + edstress(5,1:6) ) / 2.0
628 tdstrain(14,1:6) = ( tdstrain(2,1:6) + tdstrain(5,1:6) ) / 2.0
629 edstrain(15,1:6) = ( edstrain(3,1:6) + edstrain(6,1:6) ) / 2.0
630 edstress(15,1:6) = ( edstress(3,1:6) + edstress(6,1:6) ) / 2.0
631 tdstrain(15,1:6) = ( tdstrain(3,1:6) + tdstrain(6,1:6) ) / 2.0
632 else if( etype == fe_hex20n )
then
635 do j = 1, numofquadpoints(etype)
636 if( j==1 .or. j==3 .or. j==7 .or. j==9 .or. &
637 j==19 .or. j==21 .or. j==25 .or. j==27 )
then
640 edstrain(i,k) = edstrain(i,k) + func(i,ic) * gausses(j)%strain_out(k)
641 edstress(i,k) = edstress(i,k) + func(i,ic) * gausses(j)%stress_out(k)
647 edstrain(9,1:6) = ( edstrain(1,1:6) + edstrain(2,1:6) ) / 2.0
648 edstress(9,1:6) = ( edstress(1,1:6) + edstress(2,1:6) ) / 2.0
649 tdstrain(9,1:6) = ( tdstrain(1,1:6) + tdstrain(2,1:6) ) / 2.0
650 edstrain(10,1:6) = ( edstrain(2,1:6) + edstrain(3,1:6) ) / 2.0
651 edstress(10,1:6) = ( edstress(2,1:6) + edstress(3,1:6) ) / 2.0
652 tdstrain(10,1:6) = ( tdstrain(2,1:6) + tdstrain(3,1:6) ) / 2.0
653 edstrain(11,1:6) = ( edstrain(3,1:6) + edstrain(4,1:6) ) / 2.0
654 edstress(11,1:6) = ( edstress(3,1:6) + edstress(4,1:6) ) / 2.0
655 tdstrain(11,1:6) = ( tdstrain(3,1:6) + tdstrain(4,1:6) ) / 2.0
656 edstrain(12,1:6) = ( edstrain(4,1:6) + edstrain(1,1:6) ) / 2.0
657 edstress(12,1:6) = ( edstress(4,1:6) + edstress(1,1:6) ) / 2.0
658 tdstrain(12,1:6) = ( tdstrain(4,1:6) + tdstrain(1,1:6) ) / 2.0
659 edstrain(13,1:6) = ( edstrain(5,1:6) + edstrain(6,1:6) ) / 2.0
660 edstress(13,1:6) = ( edstress(5,1:6) + edstress(6,1:6) ) / 2.0
661 tdstrain(13,1:6) = ( tdstrain(5,1:6) + tdstrain(6,1:6) ) / 2.0
662 edstrain(14,1:6) = ( edstrain(6,1:6) + edstrain(7,1:6) ) / 2.0
663 edstress(14,1:6) = ( edstress(6,1:6) + edstress(7,1:6) ) / 2.0
664 tdstrain(14,1:6) = ( tdstrain(6,1:6) + tdstrain(7,1:6) ) / 2.0
665 edstrain(15,1:6) = ( edstrain(7,1:6) + edstrain(8,1:6) ) / 2.0
666 edstress(15,1:6) = ( edstress(7,1:6) + edstress(8,1:6) ) / 2.0
667 tdstrain(15,1:6) = ( tdstrain(7,1:6) + tdstrain(8,1:6) ) / 2.0
668 edstrain(16,1:6) = ( edstrain(8,1:6) + edstrain(5,1:6) ) / 2.0
669 edstress(16,1:6) = ( edstress(8,1:6) + edstress(5,1:6) ) / 2.0
670 tdstrain(16,1:6) = ( tdstrain(8,1:6) + tdstrain(5,1:6) ) / 2.0
671 edstrain(17,1:6) = ( edstrain(1,1:6) + edstrain(5,1:6) ) / 2.0
672 edstress(17,1:6) = ( edstress(1,1:6) + edstress(5,1:6) ) / 2.0
673 tdstrain(17,1:6) = ( tdstrain(1,1:6) + tdstrain(5,1:6) ) / 2.0
674 edstrain(18,1:6) = ( edstrain(2,1:6) + edstrain(6,1:6) ) / 2.0
675 edstress(18,1:6) = ( edstress(2,1:6) + edstress(6,1:6) ) / 2.0
676 tdstrain(18,1:6) = ( tdstrain(2,1:6) + tdstrain(6,1:6) ) / 2.0
677 edstrain(19,1:6) = ( edstrain(3,1:6) + edstrain(7,1:6) ) / 2.0
678 edstress(19,1:6) = ( edstress(3,1:6) + edstress(7,1:6) ) / 2.0
679 tdstrain(19,1:6) = ( tdstrain(3,1:6) + tdstrain(7,1:6) ) / 2.0
680 edstrain(20,1:6) = ( edstrain(4,1:6) + edstrain(8,1:6) ) / 2.0
681 edstress(20,1:6) = ( edstress(4,1:6) + edstress(8,1:6) ) / 2.0
682 tdstrain(20,1:6) = ( tdstrain(4,1:6) + tdstrain(8,1:6) ) / 2.0
684 end subroutine nodalstress_inv3
689 real(kind=kreal) :: s11, s22, s33, s12, s23, s13, ps, smises
697 ps = ( s11 + s22 + s33 ) / 3.0d0
698 smises = 0.5d0 * ( (s11-ps)**2 + (s22-ps)**2 + (s33-ps)**2 ) + s12**2 + s23**2 + s13**2
708 type (hecmwST_local_mesh) :: hecMESH
709 type (fstr_solid) :: fstrSOLID
710 real(kind=kreal),
pointer :: tnstrain(:), testrain(:)
712 integer(kind=kint) :: itype, icel, ic, is, iE, jS, i, j, ic_type, nn, ni, ID_area
713 real(kind=kreal) :: estrain(4), estress(4), tstrain(4), naturalcoord(4)
714 real(kind=kreal) :: edstrain(8,4), edstress(8,4), tdstrain(8,4)
715 real(kind=kreal) :: s11, s22, s33, s12, s23, s13, ps, smises
716 real(kind=kreal),
allocatable :: func(:,:), inv_func(:,:)
717 integer(kind=kint),
allocatable :: nnumber(:)
719 tnstrain => fstrsolid%tnstrain
720 testrain => fstrsolid%testrain
723 allocate( nnumber(hecmesh%n_node) )
724 if( .not.
associated(fstrsolid%is_rot) )
allocate( fstrsolid%is_rot(hecmesh%n_node) )
731 do itype = 1, hecmesh%n_elem_type
732 is = hecmesh%elem_type_index(itype-1) + 1
733 ie = hecmesh%elem_type_index(itype )
734 ic_type = hecmesh%elem_type_item(itype)
735 if( .not. hecmw_is_etype_surface(ic_type) ) cycle
737 nn = hecmw_get_max_node( ic_type )
738 ni = numofquadpoints( ic_type )
739 allocate( func(ni,nn), inv_func(nn,ni) )
740 if( ic_type == fe_tri6n )
then
741 ic = hecmw_get_max_node( fe_tri3n )
743 call getquadpoint( ic_type, i, naturalcoord )
744 call getshapefunc( fe_tri3n, naturalcoord, func(i,1:ic) )
746 call inverse_func( ic, func, inv_func )
747 else if( ic_type == fe_quad4n )
then
749 call getquadpoint( ic_type, i, naturalcoord )
750 call getshapefunc( ic_type, naturalcoord, func(i,1:nn) )
752 call inverse_func( ni, func, inv_func )
753 else if( ic_type == fe_quad8n )
then
756 if( i==1 .or. i==3 .or. i==7 .or. i==9 )
then
758 call getquadpoint( ic_type, i, naturalcoord )
759 call getshapefunc( fe_quad4n, naturalcoord, func(ic,1:4) )
762 call inverse_func( ic, func, inv_func )
767 js = hecmesh%elem_node_index(icel-1)
768 id_area = hecmesh%elem_ID(icel*2)
770 if( ic_type == fe_tri6n .or. ic_type == fe_quad4n .or. ic_type == fe_quad8n )
then
771 call nodalstress_inv2( ic_type, ni, fstrsolid%elements(icel)%gausses, &
772 inv_func, edstrain(1:nn,1:4), edstress(1:nn,1:4), &
775 call nodalstress_c2( ic_type, nn, fstrsolid%elements(icel)%gausses, &
776 edstrain(1:nn,1:4), edstress(1:nn,1:4) )
781 ic = hecmesh%elem_node_item(js+j)
782 fstrsolid%STRAIN(3*ic-2) = fstrsolid%STRAIN(3*ic-2) + edstrain(j,1)
783 fstrsolid%STRAIN(3*ic-1) = fstrsolid%STRAIN(3*ic-1) + edstrain(j,2)
784 fstrsolid%STRAIN(3*ic-0) = fstrsolid%STRAIN(3*ic-0) + edstrain(j,3)
785 fstrsolid%STRESS(3*ic-2) = fstrsolid%STRESS(3*ic-2) + edstress(j,1)
786 fstrsolid%STRESS(3*ic-1) = fstrsolid%STRESS(3*ic-1) + edstress(j,2)
787 fstrsolid%STRESS(3*ic-0) = fstrsolid%STRESS(3*ic-0) + edstress(j,3)
789 if(
associated(tnstrain) )
then
790 tnstrain(3*ic-2) = tnstrain(3*ic-2) + tdstrain(j,1)
791 tnstrain(3*ic-1) = tnstrain(3*ic-1) + tdstrain(j,2)
792 tnstrain(3*ic ) = tnstrain(3*ic ) + tdstrain(j,3)
794 nnumber(ic) = nnumber(ic) + 1
798 call elementstress_c2( ic_type, fstrsolid%elements(icel)%gausses, estrain, estress )
801 fstrsolid%ESTRAIN(3*icel-2) = estrain(1)
802 fstrsolid%ESTRAIN(3*icel-1) = estrain(2)
803 fstrsolid%ESTRAIN(3*icel-0) = estrain(3)
804 fstrsolid%ESTRESS(3*icel-2) = estress(1)
805 fstrsolid%ESTRESS(3*icel-1) = estress(2)
806 fstrsolid%ESTRESS(3*icel-0) = estress(3)
816 smises = 0.5d0 * ((s11-s22)**2+(s11)**2+(s22)**2) + 3*s12**2
817 fstrsolid%EMISES(icel) = sqrt( smises )
820 deallocate( func, inv_func )
824 do i = 1, hecmesh%n_node
825 if( nnumber(i) == 0 ) cycle
826 fstrsolid%STRAIN(3*i-2:3*i-0) = fstrsolid%STRAIN(3*i-2:3*i-0) / nnumber(i)
827 fstrsolid%STRESS(3*i-2:3*i-0) = fstrsolid%STRESS(3*i-2:3*i-0) / nnumber(i)
828 if(
associated(tnstrain) ) tnstrain(3*i-2:3*i) = tnstrain(3*i-2:3*i) / nnumber(i)
831 do i = 1, hecmesh%n_node
832 s11 = fstrsolid%STRESS(3*i-2)
833 s22 = fstrsolid%STRESS(3*i-1)
834 s12 = fstrsolid%STRESS(3*i-0)
835 smises = 0.5d0 * ((s11-s22)**2+(s11)**2+(s22)**2) + 3*s12**2
836 fstrsolid%MISES(i) = sqrt( smises )
839 deallocate( nnumber )
843 subroutine nodalstress_inv2( etype, ni, gausses, func, edstrain, edstress, tdstrain )
846 integer(kind=kint) :: etype, ni
848 real(kind=kreal) :: func(:,:), edstrain(:,:), edstress(:,:), tdstrain(:,:)
849 integer :: i, j, k, ic
855 if( etype == fe_quad4n )
then
859 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
860 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
865 else if( etype == fe_tri6n )
then
869 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
870 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
875 edstrain(4,1:4) = ( edstrain(1,1:4) + edstrain(2,1:4) ) / 2.0
876 edstress(4,1:4) = ( edstress(1,1:4) + edstress(2,1:4) ) / 2.0
877 tdstrain(4,1:4) = ( tdstrain(1,1:4) + tdstrain(2,1:4) ) / 2.0
878 edstrain(5,1:4) = ( edstrain(2,1:4) + edstrain(3,1:4) ) / 2.0
879 edstress(5,1:4) = ( edstress(2,1:4) + edstress(3,1:4) ) / 2.0
880 tdstrain(5,1:4) = ( tdstrain(2,1:4) + tdstrain(3,1:4) ) / 2.0
881 edstrain(6,1:4) = ( edstrain(3,1:4) + edstrain(1,1:4) ) / 2.0
882 edstress(6,1:4) = ( edstress(3,1:4) + edstress(1,1:4) ) / 2.0
883 tdstrain(6,1:4) = ( tdstrain(3,1:4) + tdstrain(1,1:4) ) / 2.0
884 else if( etype == fe_quad8n )
then
887 do j = 1, numofquadpoints(etype)
888 if( j==1 .or. j==3 .or. j==7 .or. j==9 )
then
891 edstrain(i,k) = edstrain(i,k) + func(i,ic) * gausses(j)%strain_out(k)
892 edstress(i,k) = edstress(i,k) + func(i,ic) * gausses(j)%stress_out(k)
898 edstrain(5,1:4) = ( edstrain(1,1:4) + edstrain(2,1:4) ) / 2.0
899 edstress(5,1:4) = ( edstress(1,1:4) + edstress(2,1:4) ) / 2.0
900 tdstrain(5,1:4) = ( tdstrain(1,1:4) + tdstrain(2,1:4) ) / 2.0
901 edstrain(6,1:4) = ( edstrain(2,1:4) + edstrain(3,1:4) ) / 2.0
902 edstress(6,1:4) = ( edstress(2,1:4) + edstress(3,1:4) ) / 2.0
903 tdstrain(6,1:4) = ( tdstrain(2,1:4) + tdstrain(3,1:4) ) / 2.0
904 edstrain(7,1:4) = ( edstrain(3,1:4) + edstrain(4,1:4) ) / 2.0
905 edstress(7,1:4) = ( edstress(3,1:4) + edstress(4,1:4) ) / 2.0
906 tdstrain(7,1:4) = ( tdstrain(3,1:4) + tdstrain(4,1:4) ) / 2.0
907 edstrain(8,1:4) = ( edstrain(4,1:4) + edstrain(1,1:4) ) / 2.0
908 edstress(8,1:4) = ( edstress(4,1:4) + edstress(1,1:4) ) / 2.0
909 tdstrain(8,1:4) = ( tdstrain(4,1:4) + tdstrain(1,1:4) ) / 2.0
911 end subroutine nodalstress_inv2
914 subroutine inverse_func( n, a, inv_a )
916 integer(kind=kint) :: n
917 real(kind=kreal) :: a(:,:), inv_a(:,:)
918 integer(kind=kint) :: i, j, k
919 real(kind=kreal) :: buf
934 a(i,j) = a(i,j) * buf
935 inv_a(i,j) = inv_a(i,j) *buf
941 a(j,k) = a(j,k) - a(i,k) * buf
942 inv_a(j,k) = inv_a(j,k) - inv_a(i,k) * buf
947 end subroutine inverse_func
954 type (hecmwST_local_mesh) :: hecMESH
955 type (fstr_solid) :: fstrSOLID
957 integer(kind=kint) :: itype, icel, is, iE, jS, i, j, k, it, ic, ic_type, nn, isect, ihead, ID_area
958 integer(kind=kint) :: nodLOCAL(20), n_layer, ntot_lyr, nlyr, n_totlyr, com_total_layer, shellmatl
959 real(kind=kreal) :: ecoord(3,9), edisp(6,9), estrain(6), estress(6), ndstrain(9,6), ndstress(9,6)
960 real(kind=kreal) :: thick, thick_layer
961 real(kind=kreal) :: s11, s22, s33, s12, s23, s13, t11, t22, t33, t12, t23, t13, ps, smises, tmises
962 integer(kind=kint),
allocatable :: nnumber(:)
967 n_totlyr = fstrsolid%max_lyr
969 allocate( nnumber(hecmesh%n_node) )
970 if( .not.
associated(fstrsolid%is_rot) )
allocate( fstrsolid%is_rot(hecmesh%n_node) )
977 do itype = 1, hecmesh%n_elem_type
978 is = hecmesh%elem_type_index(itype-1) + 1
979 ie = hecmesh%elem_type_index(itype )
980 ic_type = hecmesh%elem_type_item(itype)
981 if( .not. hecmw_is_etype_shell(ic_type) )
then
985 nn = hecmw_get_max_node( ic_type )
988 js = hecmesh%elem_node_index(icel-1)
989 id_area = hecmesh%elem_ID(icel*2)
991 nodlocal(j) = hecmesh%elem_node_item(js+j)
992 ecoord(1,j) = hecmesh%node(3*nodlocal(j)-2)
993 ecoord(2,j) = hecmesh%node(3*nodlocal(j)-1)
994 ecoord(3,j) = hecmesh%node(3*nodlocal(j) )
995 edisp(1,j) = fstrsolid%unode(6*nodlocal(j)-5)
996 edisp(2,j) = fstrsolid%unode(6*nodlocal(j)-4)
997 edisp(3,j) = fstrsolid%unode(6*nodlocal(j)-3)
998 edisp(4,j) = fstrsolid%unode(6*nodlocal(j)-2)
999 edisp(5,j) = fstrsolid%unode(6*nodlocal(j)-1)
1000 edisp(6,j) = fstrsolid%unode(6*nodlocal(j) )
1002 isect = hecmesh%section_ID(icel)
1003 ihead = hecmesh%section%sect_R_index(isect-1)
1004 thick = hecmesh%section%sect_R_item(ihead+1)
1006 if( ic_type == 731 .or. ic_type == 741 .or. ic_type == 743 )
then
1007 ntot_lyr = fstrsolid%elements(icel)%gausses(1)%pMaterial%totallyr
1009 call elementstress_shell_mitc( ic_type, nn, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
1010 & ndstrain(1:nn,1:6), ndstress(1:nn,1:6), thick, 1.0d0, nlyr, ntot_lyr)
1013 layer => fstrsolid%SHELL%LAYER(nlyr)%PLUS
1015 layer%STRAIN(6*(i-1)+k) = layer%STRAIN(6*(i-1)+k) + ndstrain(j,k)
1016 layer%STRESS(6*(i-1)+k) = layer%STRESS(6*(i-1)+k) + ndstress(j,k)
1017 layer%ESTRAIN(6*(icel-1)+k) = layer%ESTRAIN(6*(icel-1)+k) + ndstrain(j,k)/nn
1018 layer%ESTRESS(6*(icel-1)+k) = layer%ESTRESS(6*(icel-1)+k) + ndstress(j,k)/nn
1022 call elementstress_shell_mitc( ic_type, nn, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
1023 & ndstrain(1:nn,1:6), ndstress(1:nn,1:6), thick,-1.0d0, nlyr, ntot_lyr)
1026 layer => fstrsolid%SHELL%LAYER(nlyr)%MINUS
1028 layer%STRAIN(6*(i-1)+k) = layer%STRAIN(6*(i-1)+k) + ndstrain(j,k)
1029 layer%STRESS(6*(i-1)+k) = layer%STRESS(6*(i-1)+k) + ndstress(j,k)
1030 layer%ESTRAIN(6*(icel-1)+k) = layer%ESTRAIN(6*(icel-1)+k) + ndstrain(j,k)/nn
1031 layer%ESTRESS(6*(icel-1)+k) = layer%ESTRESS(6*(icel-1)+k) + ndstress(j,k)/nn
1035 call fstr_getavg_shell(nn,fstrsolid,icel,nodlocal,ndstrain(1:nn,1:6),ndstress(1:nn,1:6),estrain,estress)
1041 ic = hecmesh%elem_node_item(js+j)
1042 fstrsolid%STRAIN(6*(ic-1)+1:6*(ic-1)+6) = fstrsolid%STRAIN(6*(ic-1)+1:6*(ic-1)+6) + ndstrain(j,1:6)
1043 fstrsolid%STRESS(6*(ic-1)+1:6*(ic-1)+6) = fstrsolid%STRESS(6*(ic-1)+1:6*(ic-1)+6) + ndstress(j,1:6)
1047 nnumber(ic) = nnumber(ic) + 1
1050 fstrsolid%ESTRAIN(6*(icel-1)+1:6*(icel-1)+6) = fstrsolid%ESTRAIN(6*(icel-1)+1:6*(icel-1)+6) + estrain(1:6)
1051 fstrsolid%ESTRESS(6*(icel-1)+1:6*(icel-1)+6) = fstrsolid%ESTRESS(6*(icel-1)+1:6*(icel-1)+6) + estress(1:6)
1057 do i = 1, hecmesh%n_node
1058 if( nnumber(i) == 0 ) cycle
1059 fstrsolid%STRAIN(6*(i-1)+1:6*(i-1)+6) = fstrsolid%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1060 fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6) = fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1066 do nlyr = 1, ntot_lyr
1067 do i = 1, hecmesh%n_node
1068 fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRAIN(6*(i-1)+1:6*(i-1)+6) = &
1069 & fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1070 fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6) = &
1071 & fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1072 fstrsolid%SHELL%LAYER(nlyr)%PLUS%MISES(i) = &
1073 &
get_mises(fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6))
1075 fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRAIN(6*(i-1)+1:6*(i-1)+6) = &
1076 & fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1077 fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6) = &
1078 & fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1079 fstrsolid%SHELL%LAYER(nlyr)%MINUS%MISES(i) = &
1080 &
get_mises(fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6))
1085 do i = 1, hecmesh%n_node
1086 fstrsolid%MISES(i) =
get_mises(fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6))
1088 do i = 1, hecmesh%n_elem
1089 fstrsolid%EMISES(i) =
get_mises(fstrsolid%ESTRESS(6*(i-1)+1:6*(i-1)+6))
1091 deallocate( nnumber )
1103 integer(kind=kint) :: i, flag
1104 real(kind=
kreal) :: tmat(3, 3), tvec(3), strain(6)
1106 flag=ieor(flag,flag)
1107 if( fstrsolid%output_ctrl(3)%outinfo%on(19) .or. fstrsolid%output_ctrl(4)%outinfo%on(19) )
then
1108 if ( .not.
associated(res%PSTRESS) )
then
1109 allocate(res%PSTRESS( 3*hecmesh%n_node ))
1111 flag=ior(flag,b
'00000001')
1113 if( fstrsolid%output_ctrl(3)%outinfo%on(23) .or. fstrsolid%output_ctrl(4)%outinfo%on(23) )
then
1114 if ( .not.
associated(res%PSTRESS_VECT) )
then
1115 allocate(res%PSTRESS_VECT( 3*hecmesh%n_node ,3))
1117 flag=ior(flag,b
'00000010')
1119 if( fstrsolid%output_ctrl(3)%outinfo%on(21) .or. fstrsolid%output_ctrl(4)%outinfo%on(21) )
then
1120 if ( .not.
associated(res%PSTRAIN) )
then
1121 allocate(res%PSTRAIN( 3*hecmesh%n_node ))
1123 flag=ior(flag,b
'00000100')
1125 if( fstrsolid%output_ctrl(3)%outinfo%on(25) .or. fstrsolid%output_ctrl(4)%outinfo%on(25) )
then
1126 if ( .not.
associated(res%PSTRAIN_VECT) )
then
1127 allocate(res%PSTRAIN_VECT( 3*hecmesh%n_node ,3))
1129 flag=ior(flag,b
'00001000')
1131 if( fstrsolid%output_ctrl(3)%outinfo%on(20) .or. fstrsolid%output_ctrl(4)%outinfo%on(20) )
then
1132 if ( .not.
associated(res%EPSTRESS) )
then
1133 allocate(res%EPSTRESS( 3*hecmesh%n_elem ))
1135 flag=ior(flag,b
'00010000')
1137 if( fstrsolid%output_ctrl(3)%outinfo%on(24) .or. fstrsolid%output_ctrl(4)%outinfo%on(24) )
then
1138 if ( .not.
associated(res%EPSTRESS_VECT) )
then
1139 allocate(res%EPSTRESS_VECT( 3*hecmesh%n_elem ,3))
1141 flag=ior(flag,b
'00100000')
1143 if( fstrsolid%output_ctrl(3)%outinfo%on(22) .or. fstrsolid%output_ctrl(4)%outinfo%on(22) )
then
1144 if ( .not.
associated(res%EPSTRAIN) )
then
1145 allocate(res%EPSTRAIN( 3*hecmesh%n_elem ))
1147 flag=ior(flag,b
'01000000')
1149 if( fstrsolid%output_ctrl(3)%outinfo%on(26) .or. fstrsolid%output_ctrl(4)%outinfo%on(26) )
then
1150 if ( .not.
associated(res%EPSTRAIN_VECT) )
then
1151 allocate(res%EPSTRAIN_VECT( 3*hecmesh%n_elem ,3))
1153 flag=ior(flag,b
'10000000')
1156 if (iand(flag,b
'00000011') /= 0)
then
1157 do i = 1, hecmesh%n_node
1158 call get_principal(res%STRESS(6*i-5:6*i), tvec, tmat)
1159 if (iand(flag,b
'00000001') /= 0) res%PSTRESS(3*(i-1)+1:3*(i-1)+3)=tvec
1160 if (iand(flag,b
'00000010') /= 0) res%PSTRESS_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat
1163 if (iand(flag,b
'00001100') /= 0)
then
1164 do i = 1, hecmesh%n_node
1165 strain(1:6) = res%STRAIN(6*i-5:6*i)
1166 strain(4:6) = 0.5d0*strain(4:6)
1167 call get_principal(strain, tvec, tmat)
1168 if (iand(flag,b
'00000100') /= 0) res%PSTRAIN(3*(i-1)+1:3*(i-1)+3)=tvec
1169 if (iand(flag,b
'00001000') /= 0) res%PSTRAIN_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat
1173 if (iand(flag,b
'00110000') /= 0)
then
1174 do i = 1, hecmesh%n_elem
1175 call get_principal( res%ESTRESS(6*i-5:6*i), tvec, tmat)
1176 if (iand(flag,b
'00010000') /= 0) res%EPSTRESS(3*(i-1)+1:3*(i-1)+3)=tvec
1177 if (iand(flag,b
'00100000') /= 0) res%EPSTRESS_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat
1180 if (iand(flag,b
'11000000') /= 0)
then
1181 do i = 1, hecmesh%n_elem
1182 strain(1:6) = res%ESTRAIN(6*i-5:6*i)
1183 strain(4:6) = 0.5d0*strain(4:6)
1184 call get_principal(strain, tvec, tmat)
1185 if (iand(flag,b
'01000000') /= 0) res%EPSTRAIN(3*(i-1)+1:3*(i-1)+3)=tvec
1186 if (iand(flag,b
'10000000') /= 0) res%EPSTRAIN_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat