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))
271 do i = 1, hecmesh%n_elem
272 if (.not.
associated(fstrsolid%elements(i)%gausses)) cycle
273 fstrsolid%EPLSTRAIN(i) =
get_pl_estrain(fstrsolid%elements(i)%gausses)
276 if( flag33 == 1 )
then
277 if( fstrsolid%output_ctrl(3)%outinfo%on(27) .or. fstrsolid%output_ctrl(4)%outinfo%on(27) )
then
278 do nlyr = 1, ntot_lyr
279 call make_principal(fstrsolid, hecmesh, fstrsolid%SHELL%LAYER(nlyr)%PLUS)
280 call make_principal(fstrsolid, hecmesh, fstrsolid%SHELL%LAYER(nlyr)%MINUS)
288 deallocate( nnumber )
293 integer(kind=kint),
allocatable,
intent(in) :: irow(:)
294 integer(kind=kint),
allocatable,
intent(in) :: asect(:)
295 integer(kind=kint),
intent(in) :: nid
296 integer(kind=kint),
intent(in) :: sid
298 integer(kind=kint) :: i
301 do i=irow(nid-1)+1,irow(nid)
302 if( asect(i) == sid )
then
311 Nodal_STRAIN, Nodal_STRESS, Elemental_STRAIN, Elemental_STRESS )
312 type(hecmwst_local_mesh),
intent(in) :: hecMESH
314 integer(kind=kint),
allocatable,
intent(inout) :: nnumber(:)
315 real(kind=kreal),
pointer,
intent(inout) :: nodal_strain(:)
316 real(kind=kreal),
pointer,
intent(inout) :: nodal_stress(:)
317 real(kind=kreal),
pointer,
intent(inout) :: elemental_strain(:)
318 real(kind=kreal),
pointer,
intent(inout) :: elemental_stress(:)
320 integer(kind=kint) :: itype, iS, iE, jS, ic_type, icel, i, j, isect
321 integer(kind=kint) :: nsize, nid(2), idx(2), nd
322 integer(kind=kint) :: nnode, nlen
323 type(hecmwst_varray_int),
allocatable :: nodal_sections(:)
324 real(kind=kreal) :: tmpval(6), hydval, nsecdup
325 integer(kind=kint),
allocatable :: irow(:), jcol(:), asect(:)
326 real(kind=kreal),
allocatable :: stress_hyd(:), strain_hyd(:)
327 real(kind=kreal),
allocatable :: stress_dev(:)
328 real(kind=kreal),
allocatable :: plstrain_dev(:)
329 real(kind=kreal) :: stress_hyd_ndave(6), strain_hyd_ndave(6)
330 real(kind=kreal) :: stress_dev_ndave(6), strain_dev_ndave(6)
331 real(kind=kreal),
allocatable :: n_dup_dev(:), n_dup_hyd(:)
332 real(kind=kreal) :: edstrain(6), edstress(6)
333 real(kind=kreal) :: edplstrain
335 nnode = hecmesh%n_node
336 nsize =
size(nodal_strain)
339 call hecmw_varray_int_initialize_all( nodal_sections, nnode, 2 )
340 do itype = 1, hecmesh%n_elem_type
341 ic_type = hecmesh%elem_type_item(itype)
342 if( ic_type /= 341 ) cycle
344 is = hecmesh%elem_type_index(itype-1) + 1
345 ie = hecmesh%elem_type_index(itype )
348 isect= hecmesh%section_ID(icel)
349 if( fstrsolid%sections(isect)%elemopt341 /=
kel341sesns ) cycle
350 js = hecmesh%elem_node_index(icel-1)
352 nd = hecmesh%elem_node_item(js+i)
353 call hecmw_varray_int_add_if_not_exits( nodal_sections(nd), isect )
359 allocate(irow(0:nnode))
362 irow(i) = irow(i-1)+hecmw_varray_int_get_nitem(nodal_sections(i))
366 allocate(asect(nlen))
368 if( irow(i-1) == irow(i) ) cycle
369 call hecmw_varray_int_get_item_all( nodal_sections(i), asect(irow(i-1)+1:irow(i)) )
373 allocate(stress_hyd(6*nlen), strain_hyd(6*nlen))
374 allocate(stress_dev(6*nlen))
375 allocate(plstrain_dev(nlen))
376 allocate(n_dup_dev(nlen),n_dup_hyd(nlen))
381 plstrain_dev(:) = 0.d0
384 do itype = 1, hecmesh%n_elem_type
385 ic_type = hecmesh%elem_type_item(itype)
386 if( ic_type /= 881 .and. ic_type /= 891 ) cycle
388 is = hecmesh%elem_type_index(itype-1) + 1
389 ie = hecmesh%elem_type_index(itype )
392 js = hecmesh%elem_node_index(icel-1)
393 isect= hecmesh%section_ID(icel)
394 if( ic_type == 881 )
then
395 nid(1) = hecmesh%elem_node_item(js+1)
399 strain_hyd(6*idx(1)-5:6*idx(1)) = fstrsolid%elements(icel)%gausses(1)%strain_out(1:6)
401 stress_hyd(6*idx(1)-5:6*idx(1)) = fstrsolid%elements(icel)%gausses(1)%stress_out(1:6)
403 n_dup_hyd(idx(1)) = n_dup_hyd(idx(1)) + 1.d0
404 else if( ic_type == 891 )
then
405 nid(1:2) = hecmesh%elem_node_item(js+1:js+2)
410 tmpval(1:6) = fstrsolid%elements(icel)%gausses(1)%stress_out(1:6)
411 stress_dev(6*idx(1)-5:6*idx(1)) = stress_dev(6*idx(1)-5:6*idx(1)) + tmpval(1:6)
412 stress_dev(6*idx(2)-5:6*idx(2)) = stress_dev(6*idx(2)-5:6*idx(2)) + tmpval(1:6)
414 plstrain_dev(idx(1)) = plstrain_dev(idx(1)) + fstrsolid%elements(icel)%gausses(1)%plstrain
415 plstrain_dev(idx(2)) = plstrain_dev(idx(2)) + fstrsolid%elements(icel)%gausses(1)%plstrain
417 n_dup_dev(idx(1)) = n_dup_dev(idx(1)) + 1.d0
418 n_dup_dev(idx(2)) = n_dup_dev(idx(2)) + 1.d0
424 if( irow(i-1) == irow(i) ) cycle
425 do j=irow(i-1)+1,irow(i)
426 if( n_dup_dev(j) < 1.0d-8 ) cycle
427 stress_dev(6*j-5:6*j) = stress_dev(6*j-5:6*j)/n_dup_dev(j)
428 plstrain_dev(j) = plstrain_dev(j)/n_dup_dev(j)
434 if( irow(i-1) == irow(i) ) cycle
435 strain_hyd_ndave(:) = 0.d0
436 stress_hyd_ndave(:) = 0.d0
437 stress_dev_ndave(:) = 0.d0
438 do j=irow(i-1)+1,irow(i)
439 strain_hyd_ndave(1:6) = strain_hyd_ndave(1:6) + strain_hyd(6*j-5:6*j)
440 stress_hyd_ndave(1:6) = stress_hyd_ndave(1:6) + stress_hyd(6*j-5:6*j)
441 stress_dev_ndave(1:6) = stress_dev_ndave(1:6) + stress_dev(6*j-5:6*j)
443 nsecdup = dble(irow(i)-irow(i-1))
444 strain_hyd_ndave(1:6) = strain_hyd_ndave(1:6)/nsecdup
445 stress_hyd_ndave(1:6) = stress_hyd_ndave(1:6)/nsecdup
446 stress_dev_ndave(1:6) = stress_dev_ndave(1:6)/nsecdup
448 if( nnumber(i) == 0 )
then
449 nodal_strain(6*i-5:6*i) = strain_hyd_ndave(1:6)
450 nodal_stress(6*i-5:6*i) = stress_hyd_ndave(1:6)+stress_dev_ndave(1:6)
452 nodal_strain(6*i-5:6*i) = 0.5d0*(nodal_strain(6*i-5:6*i)+strain_hyd_ndave(1:6))
453 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))
458 do itype = 1, hecmesh%n_elem_type
459 ic_type = hecmesh%elem_type_item(itype)
460 if( ic_type /= 341 ) cycle
462 is = hecmesh%elem_type_index(itype-1) + 1
463 ie = hecmesh%elem_type_index(itype )
466 isect= hecmesh%section_ID(icel)
467 if( fstrsolid%sections(isect)%elemopt341 /=
kel341sesns ) cycle
468 js = hecmesh%elem_node_index(icel-1)
473 nd = hecmesh%elem_node_item(js+i)
474 idx(1) =
search_idx_senes( irow, asect, hecmesh%elem_node_item(js+i), isect )
475 edstrain(1:6) = edstrain(1:6) + strain_hyd(6*idx(1)-5:6*idx(1))
476 edstress(1:6) = edstress(1:6) + stress_hyd(6*idx(1)-5:6*idx(1)) + stress_dev(6*idx(1)-5:6*idx(1))
477 edplstrain = edplstrain + plstrain_dev(idx(1))
479 edstrain(1:6) = 0.25d0*edstrain(1:6)
480 edstress(1:6) = 0.25d0*edstress(1:6)
481 edplstrain = 0.25d0*edplstrain
483 elemental_strain(6*(icel-1)+1:6*(icel-1)+6) = elemental_strain(6*(icel-1)+1:6*(icel-1)+6) + edstrain(1:6)
484 elemental_stress(6*(icel-1)+1:6*(icel-1)+6) = elemental_stress(6*(icel-1)+1:6*(icel-1)+6) + edstress(1:6)
486 fstrsolid%elements(icel)%gausses(1)%strain_out(1:6) = elemental_strain(6*(icel-1)+1:6*(icel-1)+6)
487 fstrsolid%elements(icel)%gausses(1)%stress_out(1:6) = elemental_stress(6*(icel-1)+1:6*(icel-1)+6)
488 fstrsolid%elements(icel)%gausses(1)%plstrain = edplstrain
492 deallocate(stress_hyd, strain_hyd)
493 deallocate(stress_dev, plstrain_dev)
494 deallocate(n_dup_dev, n_dup_hyd)
501 integer(kind=kint) :: nodlocal(20)
502 integer(kind=kint) :: nn, i, j, k, m, nlyr, weight, icel, flag
503 real(kind=kreal) :: strain(nn, 6), stress(nn, 6)
510 layer => fstrsolid%SHELL%LAYER(nlyr)%PLUS
511 elseif(flag == -1)
then
512 layer => fstrsolid%SHELL%LAYER(nlyr)%MINUS
515 layer%STRAIN(6*(i-1)+k) = layer%STRAIN(6*(i-1)+k) + strain(j,k)
516 layer%STRAIN(6*(m-1)+k) = layer%STRAIN(6*(m-1)+k) + strain(j,k)
517 layer%STRESS(6*(i-1)+k) = layer%STRESS(6*(i-1)+k) + stress(j,k)
518 layer%STRESS(6*(m-1)+k) = layer%STRESS(6*(m-1)+k) + stress(j,k)
519 layer%ESTRAIN(6*(icel-1)+k) = layer%ESTRAIN(6*(icel-1)+k) + strain(j,k)/nn
520 layer%ESTRESS(6*(icel-1)+k) = layer%ESTRESS(6*(icel-1)+k) + stress(j,k)/nn
528 integer(kind=kint) :: nodlocal(20)
529 integer(kind=kint) :: nn, i, j, k, m, nlyr, icel, flag, ntot_lyr
530 real(kind=kreal) :: strain(nn,6), stress(nn,6), estrain(6), estress(6), weight
533 ntot_lyr = fstrsolid%elements(icel)%gausses(1)%pMaterial%totallyr
539 do nlyr = 1, ntot_lyr
540 layer => fstrsolid%SHELL%LAYER(nlyr)
541 weight = fstrsolid%elements(icel)%gausses(1)%pMaterial%shell_var(nlyr)%weight
545 strain(j,k) = strain(j,k) &
546 & + weight*(0.5d0*layer%PLUS%STRAIN(6*(i-1)+k) + 0.5d0*layer%MINUS%STRAIN(6*(i-1)+k))
547 stress(j,k) = stress(j,k) &
548 & + weight*(0.5d0*layer%PLUS%STRESS(6*(i-1)+k) + 0.5d0*layer%MINUS%STRESS(6*(i-1)+k))
550 estrain(j) = estrain(j) &
551 & + weight*(0.5d0*layer%PLUS%ESTRAIN(6*(icel-1)+j) + 0.5d0*layer%MINUS%ESTRAIN(6*(icel-1)+j))
552 estress(j) = estress(j) &
553 & + weight*(0.5d0*layer%PLUS%ESTRESS(6*(icel-1)+j) + 0.5d0*layer%MINUS%ESTRESS(6*(icel-1)+j))
559 subroutine nodalstress_inv3( etype, ni, gausses, func, edstrain, edstress, tdstrain )
562 integer(kind=kint) :: etype, ni
564 real(kind=kreal) :: func(:, :), edstrain(:, :), edstress(:, :), tdstrain(:, :)
565 integer :: i, j, k, ic
571 if( etype == fe_hex8n )
then
575 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
576 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
581 else if( etype == fe_tet10n )
then
585 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
586 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
591 edstrain(5,1:6) = ( edstrain(1,1:6) + edstrain(2,1:6) ) / 2.0
592 edstress(5,1:6) = ( edstress(1,1:6) + edstress(2,1:6) ) / 2.0
593 tdstrain(5,1:6) = ( tdstrain(1,1:6) + tdstrain(2,1:6) ) / 2.0
594 edstrain(6,1:6) = ( edstrain(2,1:6) + edstrain(3,1:6) ) / 2.0
595 edstress(6,1:6) = ( edstress(2,1:6) + edstress(3,1:6) ) / 2.0
596 tdstrain(6,1:6) = ( tdstrain(2,1:6) + tdstrain(3,1:6) ) / 2.0
597 edstrain(7,1:6) = ( edstrain(3,1:6) + edstrain(1,1:6) ) / 2.0
598 edstress(7,1:6) = ( edstress(3,1:6) + edstress(1,1:6) ) / 2.0
599 tdstrain(7,1:6) = ( tdstrain(3,1:6) + tdstrain(1,1:6) ) / 2.0
600 edstrain(8,1:6) = ( edstrain(1,1:6) + edstrain(4,1:6) ) / 2.0
601 edstress(8,1:6) = ( edstress(1,1:6) + edstress(4,1:6) ) / 2.0
602 tdstrain(8,1:6) = ( tdstrain(1,1:6) + tdstrain(4,1:6) ) / 2.0
603 edstrain(9,1:6) = ( edstrain(2,1:6) + edstrain(4,1:6) ) / 2.0
604 edstress(9,1:6) = ( edstress(2,1:6) + edstress(4,1:6) ) / 2.0
605 tdstrain(9,1:6) = ( tdstrain(2,1:6) + tdstrain(4,1:6) ) / 2.0
606 edstrain(10,1:6) = ( edstrain(3,1:6) + edstrain(4,1:6) ) / 2.0
607 edstress(10,1:6) = ( edstress(3,1:6) + edstress(4,1:6) ) / 2.0
608 tdstrain(10,1:6) = ( tdstrain(3,1:6) + tdstrain(4,1:6) ) / 2.0
609 else if( etype == fe_prism15n )
then
612 do j = 1, numofquadpoints(etype)
613 if( j==1 .or. j==2 .or. j==3 .or. j==7 .or. j==8 .or. j==9 )
then
616 edstrain(i,k) = edstrain(i,k) + func(i,ic) * gausses(j)%strain_out(k)
617 edstress(i,k) = edstress(i,k) + func(i,ic) * gausses(j)%stress_out(k)
623 edstrain(7,1:6) = ( edstrain(1,1:6) + edstrain(2,1:6) ) / 2.0
624 edstress(7,1:6) = ( edstress(1,1:6) + edstress(2,1:6) ) / 2.0
625 tdstrain(7,1:6) = ( tdstrain(1,1:6) + tdstrain(2,1:6) ) / 2.0
626 edstrain(8,1:6) = ( edstrain(2,1:6) + edstrain(3,1:6) ) / 2.0
627 edstress(8,1:6) = ( edstress(2,1:6) + edstress(3,1:6) ) / 2.0
628 tdstrain(8,1:6) = ( tdstrain(2,1:6) + tdstrain(3,1:6) ) / 2.0
629 edstrain(9,1:6) = ( edstrain(3,1:6) + edstrain(1,1:6) ) / 2.0
630 edstress(9,1:6) = ( edstress(3,1:6) + edstress(1,1:6) ) / 2.0
631 tdstrain(9,1:6) = ( tdstrain(3,1:6) + tdstrain(1,1:6) ) / 2.0
632 edstrain(10,1:6) = ( edstrain(4,1:6) + edstrain(5,1:6) ) / 2.0
633 edstress(10,1:6) = ( edstress(4,1:6) + edstress(5,1:6) ) / 2.0
634 tdstrain(10,1:6) = ( tdstrain(4,1:6) + tdstrain(5,1:6) ) / 2.0
635 edstrain(11,1:6) = ( edstrain(5,1:6) + edstrain(6,1:6) ) / 2.0
636 edstress(11,1:6) = ( edstress(5,1:6) + edstress(6,1:6) ) / 2.0
637 tdstrain(11,1:6) = ( tdstrain(5,1:6) + tdstrain(6,1:6) ) / 2.0
638 edstrain(12,1:6) = ( edstrain(6,1:6) + edstrain(4,1:6) ) / 2.0
639 edstress(12,1:6) = ( edstress(6,1:6) + edstress(4,1:6) ) / 2.0
640 tdstrain(12,1:6) = ( tdstrain(6,1:6) + tdstrain(4,1:6) ) / 2.0
641 edstrain(13,1:6) = ( edstrain(1,1:6) + edstrain(4,1:6) ) / 2.0
642 edstress(13,1:6) = ( edstress(1,1:6) + edstress(4,1:6) ) / 2.0
643 tdstrain(13,1:6) = ( tdstrain(1,1:6) + tdstrain(4,1:6) ) / 2.0
644 edstrain(14,1:6) = ( edstrain(2,1:6) + edstrain(5,1:6) ) / 2.0
645 edstress(14,1:6) = ( edstress(2,1:6) + edstress(5,1:6) ) / 2.0
646 tdstrain(14,1:6) = ( tdstrain(2,1:6) + tdstrain(5,1:6) ) / 2.0
647 edstrain(15,1:6) = ( edstrain(3,1:6) + edstrain(6,1:6) ) / 2.0
648 edstress(15,1:6) = ( edstress(3,1:6) + edstress(6,1:6) ) / 2.0
649 tdstrain(15,1:6) = ( tdstrain(3,1:6) + tdstrain(6,1:6) ) / 2.0
650 else if( etype == fe_hex20n )
then
653 do j = 1, numofquadpoints(etype)
654 if( j==1 .or. j==3 .or. j==7 .or. j==9 .or. &
655 j==19 .or. j==21 .or. j==25 .or. j==27 )
then
658 edstrain(i,k) = edstrain(i,k) + func(i,ic) * gausses(j)%strain_out(k)
659 edstress(i,k) = edstress(i,k) + func(i,ic) * gausses(j)%stress_out(k)
665 edstrain(9,1:6) = ( edstrain(1,1:6) + edstrain(2,1:6) ) / 2.0
666 edstress(9,1:6) = ( edstress(1,1:6) + edstress(2,1:6) ) / 2.0
667 tdstrain(9,1:6) = ( tdstrain(1,1:6) + tdstrain(2,1:6) ) / 2.0
668 edstrain(10,1:6) = ( edstrain(2,1:6) + edstrain(3,1:6) ) / 2.0
669 edstress(10,1:6) = ( edstress(2,1:6) + edstress(3,1:6) ) / 2.0
670 tdstrain(10,1:6) = ( tdstrain(2,1:6) + tdstrain(3,1:6) ) / 2.0
671 edstrain(11,1:6) = ( edstrain(3,1:6) + edstrain(4,1:6) ) / 2.0
672 edstress(11,1:6) = ( edstress(3,1:6) + edstress(4,1:6) ) / 2.0
673 tdstrain(11,1:6) = ( tdstrain(3,1:6) + tdstrain(4,1:6) ) / 2.0
674 edstrain(12,1:6) = ( edstrain(4,1:6) + edstrain(1,1:6) ) / 2.0
675 edstress(12,1:6) = ( edstress(4,1:6) + edstress(1,1:6) ) / 2.0
676 tdstrain(12,1:6) = ( tdstrain(4,1:6) + tdstrain(1,1:6) ) / 2.0
677 edstrain(13,1:6) = ( edstrain(5,1:6) + edstrain(6,1:6) ) / 2.0
678 edstress(13,1:6) = ( edstress(5,1:6) + edstress(6,1:6) ) / 2.0
679 tdstrain(13,1:6) = ( tdstrain(5,1:6) + tdstrain(6,1:6) ) / 2.0
680 edstrain(14,1:6) = ( edstrain(6,1:6) + edstrain(7,1:6) ) / 2.0
681 edstress(14,1:6) = ( edstress(6,1:6) + edstress(7,1:6) ) / 2.0
682 tdstrain(14,1:6) = ( tdstrain(6,1:6) + tdstrain(7,1:6) ) / 2.0
683 edstrain(15,1:6) = ( edstrain(7,1:6) + edstrain(8,1:6) ) / 2.0
684 edstress(15,1:6) = ( edstress(7,1:6) + edstress(8,1:6) ) / 2.0
685 tdstrain(15,1:6) = ( tdstrain(7,1:6) + tdstrain(8,1:6) ) / 2.0
686 edstrain(16,1:6) = ( edstrain(8,1:6) + edstrain(5,1:6) ) / 2.0
687 edstress(16,1:6) = ( edstress(8,1:6) + edstress(5,1:6) ) / 2.0
688 tdstrain(16,1:6) = ( tdstrain(8,1:6) + tdstrain(5,1:6) ) / 2.0
689 edstrain(17,1:6) = ( edstrain(1,1:6) + edstrain(5,1:6) ) / 2.0
690 edstress(17,1:6) = ( edstress(1,1:6) + edstress(5,1:6) ) / 2.0
691 tdstrain(17,1:6) = ( tdstrain(1,1:6) + tdstrain(5,1:6) ) / 2.0
692 edstrain(18,1:6) = ( edstrain(2,1:6) + edstrain(6,1:6) ) / 2.0
693 edstress(18,1:6) = ( edstress(2,1:6) + edstress(6,1:6) ) / 2.0
694 tdstrain(18,1:6) = ( tdstrain(2,1:6) + tdstrain(6,1:6) ) / 2.0
695 edstrain(19,1:6) = ( edstrain(3,1:6) + edstrain(7,1:6) ) / 2.0
696 edstress(19,1:6) = ( edstress(3,1:6) + edstress(7,1:6) ) / 2.0
697 tdstrain(19,1:6) = ( tdstrain(3,1:6) + tdstrain(7,1:6) ) / 2.0
698 edstrain(20,1:6) = ( edstrain(4,1:6) + edstrain(8,1:6) ) / 2.0
699 edstress(20,1:6) = ( edstress(4,1:6) + edstress(8,1:6) ) / 2.0
700 tdstrain(20,1:6) = ( tdstrain(4,1:6) + tdstrain(8,1:6) ) / 2.0
702 end subroutine nodalstress_inv3
707 real(kind=kreal) :: s11, s22, s33, s12, s23, s13, ps, smises
715 ps = ( s11 + s22 + s33 ) / 3.0d0
716 smises = 0.5d0 * ( (s11-ps)**2 + (s22-ps)**2 + (s33-ps)**2 ) + s12**2 + s23**2 + s13**2
724 type(tgaussstatus) :: gausses(:)
725 integer(kind=kint) :: i
728 do i = 1,
size(gausses)
740 type (hecmwst_local_mesh) :: hecMESH
741 type (fstr_solid) :: fstrSOLID
742 real(kind=kreal),
pointer :: tnstrain(:), testrain(:)
744 integer(kind=kint) :: itype, icel, ic, is, iE, jS, i, j, ic_type, nn, ni, ID_area
745 real(kind=kreal) :: estrain(4), estress(4), tstrain(4), naturalcoord(4)
746 real(kind=kreal) :: edstrain(8,4), edstress(8,4), tdstrain(8,4)
747 real(kind=kreal) :: s11, s22, s33, s12, s23, s13, ps, smises
748 real(kind=kreal),
allocatable :: func(:,:), inv_func(:,:)
749 integer(kind=kint),
allocatable :: nnumber(:)
751 tnstrain => fstrsolid%tnstrain
752 testrain => fstrsolid%testrain
755 allocate( nnumber(hecmesh%n_node) )
756 if( .not.
associated(fstrsolid%is_rot) )
allocate( fstrsolid%is_rot(hecmesh%n_node) )
763 do itype = 1, hecmesh%n_elem_type
764 is = hecmesh%elem_type_index(itype-1) + 1
765 ie = hecmesh%elem_type_index(itype )
766 ic_type = hecmesh%elem_type_item(itype)
767 if( .not. hecmw_is_etype_surface(ic_type) ) cycle
769 nn = hecmw_get_max_node( ic_type )
770 ni = numofquadpoints( ic_type )
771 allocate( func(ni,nn), inv_func(nn,ni) )
772 if( ic_type == fe_tri6n )
then
773 ic = hecmw_get_max_node( fe_tri3n )
775 call getquadpoint( ic_type, i, naturalcoord )
776 call getshapefunc( fe_tri3n, naturalcoord, func(i,1:ic) )
778 call inverse_func( ic, func, inv_func )
779 else if( ic_type == fe_quad4n )
then
781 call getquadpoint( ic_type, i, naturalcoord )
782 call getshapefunc( ic_type, naturalcoord, func(i,1:nn) )
784 call inverse_func( ni, func, inv_func )
785 else if( ic_type == fe_quad8n )
then
788 if( i==1 .or. i==3 .or. i==7 .or. i==9 )
then
790 call getquadpoint( ic_type, i, naturalcoord )
791 call getshapefunc( fe_quad4n, naturalcoord, func(ic,1:4) )
794 call inverse_func( ic, func, inv_func )
799 js = hecmesh%elem_node_index(icel-1)
800 id_area = hecmesh%elem_ID(icel*2)
802 if( ic_type == fe_tri6n .or. ic_type == fe_quad4n .or. ic_type == fe_quad8n )
then
803 call nodalstress_inv2( ic_type, ni, fstrsolid%elements(icel)%gausses, &
804 inv_func, edstrain(1:nn,1:4), edstress(1:nn,1:4), &
807 call nodalstress_c2( ic_type, nn, fstrsolid%elements(icel)%gausses, &
808 edstrain(1:nn,1:4), edstress(1:nn,1:4) )
813 ic = hecmesh%elem_node_item(js+j)
814 fstrsolid%STRAIN(3*ic-2) = fstrsolid%STRAIN(3*ic-2) + edstrain(j,1)
815 fstrsolid%STRAIN(3*ic-1) = fstrsolid%STRAIN(3*ic-1) + edstrain(j,2)
816 fstrsolid%STRAIN(3*ic-0) = fstrsolid%STRAIN(3*ic-0) + edstrain(j,3)
817 fstrsolid%STRESS(3*ic-2) = fstrsolid%STRESS(3*ic-2) + edstress(j,1)
818 fstrsolid%STRESS(3*ic-1) = fstrsolid%STRESS(3*ic-1) + edstress(j,2)
819 fstrsolid%STRESS(3*ic-0) = fstrsolid%STRESS(3*ic-0) + edstress(j,3)
821 if(
associated(tnstrain) )
then
822 tnstrain(3*ic-2) = tnstrain(3*ic-2) + tdstrain(j,1)
823 tnstrain(3*ic-1) = tnstrain(3*ic-1) + tdstrain(j,2)
824 tnstrain(3*ic ) = tnstrain(3*ic ) + tdstrain(j,3)
826 nnumber(ic) = nnumber(ic) + 1
830 call elementstress_c2( ic_type, fstrsolid%elements(icel)%gausses, estrain, estress )
833 fstrsolid%ESTRAIN(3*icel-2) = estrain(1)
834 fstrsolid%ESTRAIN(3*icel-1) = estrain(2)
835 fstrsolid%ESTRAIN(3*icel-0) = estrain(3)
836 fstrsolid%ESTRESS(3*icel-2) = estress(1)
837 fstrsolid%ESTRESS(3*icel-1) = estress(2)
838 fstrsolid%ESTRESS(3*icel-0) = estress(3)
848 smises = 0.5d0 * ((s11-s22)**2+(s11)**2+(s22)**2) + 3*s12**2
849 fstrsolid%EMISES(icel) = sqrt( smises )
852 deallocate( func, inv_func )
856 do i = 1, hecmesh%n_node
857 if( nnumber(i) == 0 ) cycle
858 fstrsolid%STRAIN(3*i-2:3*i-0) = fstrsolid%STRAIN(3*i-2:3*i-0) / nnumber(i)
859 fstrsolid%STRESS(3*i-2:3*i-0) = fstrsolid%STRESS(3*i-2:3*i-0) / nnumber(i)
860 if(
associated(tnstrain) ) tnstrain(3*i-2:3*i) = tnstrain(3*i-2:3*i) / nnumber(i)
863 do i = 1, hecmesh%n_node
864 s11 = fstrsolid%STRESS(3*i-2)
865 s22 = fstrsolid%STRESS(3*i-1)
866 s12 = fstrsolid%STRESS(3*i-0)
867 smises = 0.5d0 * ((s11-s22)**2+(s11)**2+(s22)**2) + 3*s12**2
868 fstrsolid%MISES(i) = sqrt( smises )
871 deallocate( nnumber )
875 subroutine nodalstress_inv2( etype, ni, gausses, func, edstrain, edstress, tdstrain )
878 integer(kind=kint) :: etype, ni
880 real(kind=kreal) :: func(:,:), edstrain(:,:), edstress(:,:), tdstrain(:,:)
881 integer :: i, j, k, ic
887 if( etype == fe_quad4n )
then
891 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
892 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
897 else if( etype == fe_tri6n )
then
901 edstrain(i,k) = edstrain(i,k) + func(i,j) * gausses(j)%strain_out(k)
902 edstress(i,k) = edstress(i,k) + func(i,j) * gausses(j)%stress_out(k)
907 edstrain(4,1:4) = ( edstrain(1,1:4) + edstrain(2,1:4) ) / 2.0
908 edstress(4,1:4) = ( edstress(1,1:4) + edstress(2,1:4) ) / 2.0
909 tdstrain(4,1:4) = ( tdstrain(1,1:4) + tdstrain(2,1:4) ) / 2.0
910 edstrain(5,1:4) = ( edstrain(2,1:4) + edstrain(3,1:4) ) / 2.0
911 edstress(5,1:4) = ( edstress(2,1:4) + edstress(3,1:4) ) / 2.0
912 tdstrain(5,1:4) = ( tdstrain(2,1:4) + tdstrain(3,1:4) ) / 2.0
913 edstrain(6,1:4) = ( edstrain(3,1:4) + edstrain(1,1:4) ) / 2.0
914 edstress(6,1:4) = ( edstress(3,1:4) + edstress(1,1:4) ) / 2.0
915 tdstrain(6,1:4) = ( tdstrain(3,1:4) + tdstrain(1,1:4) ) / 2.0
916 else if( etype == fe_quad8n )
then
919 do j = 1, numofquadpoints(etype)
920 if( j==1 .or. j==3 .or. j==7 .or. j==9 )
then
923 edstrain(i,k) = edstrain(i,k) + func(i,ic) * gausses(j)%strain_out(k)
924 edstress(i,k) = edstress(i,k) + func(i,ic) * gausses(j)%stress_out(k)
930 edstrain(5,1:4) = ( edstrain(1,1:4) + edstrain(2,1:4) ) / 2.0
931 edstress(5,1:4) = ( edstress(1,1:4) + edstress(2,1:4) ) / 2.0
932 tdstrain(5,1:4) = ( tdstrain(1,1:4) + tdstrain(2,1:4) ) / 2.0
933 edstrain(6,1:4) = ( edstrain(2,1:4) + edstrain(3,1:4) ) / 2.0
934 edstress(6,1:4) = ( edstress(2,1:4) + edstress(3,1:4) ) / 2.0
935 tdstrain(6,1:4) = ( tdstrain(2,1:4) + tdstrain(3,1:4) ) / 2.0
936 edstrain(7,1:4) = ( edstrain(3,1:4) + edstrain(4,1:4) ) / 2.0
937 edstress(7,1:4) = ( edstress(3,1:4) + edstress(4,1:4) ) / 2.0
938 tdstrain(7,1:4) = ( tdstrain(3,1:4) + tdstrain(4,1:4) ) / 2.0
939 edstrain(8,1:4) = ( edstrain(4,1:4) + edstrain(1,1:4) ) / 2.0
940 edstress(8,1:4) = ( edstress(4,1:4) + edstress(1,1:4) ) / 2.0
941 tdstrain(8,1:4) = ( tdstrain(4,1:4) + tdstrain(1,1:4) ) / 2.0
943 end subroutine nodalstress_inv2
946 subroutine inverse_func( n, a, inv_a )
948 integer(kind=kint) :: n
949 real(kind=kreal) :: a(:,:), inv_a(:,:)
950 integer(kind=kint) :: i, j, k
951 real(kind=kreal) :: buf
966 a(i,j) = a(i,j) * buf
967 inv_a(i,j) = inv_a(i,j) *buf
973 a(j,k) = a(j,k) - a(i,k) * buf
974 inv_a(j,k) = inv_a(j,k) - inv_a(i,k) * buf
979 end subroutine inverse_func
986 type (hecmwST_local_mesh) :: hecMESH
987 type (fstr_solid) :: fstrSOLID
989 integer(kind=kint) :: itype, icel, is, iE, jS, i, j, k, it, ic, ic_type, nn, isect, ihead, ID_area
990 integer(kind=kint) :: nodLOCAL(20), n_layer, ntot_lyr, nlyr, n_totlyr, com_total_layer, shellmatl
991 real(kind=kreal) :: ecoord(3,9), edisp(6,9), estrain(6), estress(6), ndstrain(9,6), ndstress(9,6)
992 real(kind=kreal) :: thick, thick_layer
993 real(kind=kreal) :: s11, s22, s33, s12, s23, s13, t11, t22, t33, t12, t23, t13, ps, smises, tmises
994 integer(kind=kint),
allocatable :: nnumber(:)
999 n_totlyr = fstrsolid%max_lyr
1001 allocate( nnumber(hecmesh%n_node) )
1002 if( .not.
associated(fstrsolid%is_rot) )
allocate( fstrsolid%is_rot(hecmesh%n_node) )
1004 fstrsolid%is_rot = 0
1009 do itype = 1, hecmesh%n_elem_type
1010 is = hecmesh%elem_type_index(itype-1) + 1
1011 ie = hecmesh%elem_type_index(itype )
1012 ic_type = hecmesh%elem_type_item(itype)
1013 if( .not. hecmw_is_etype_shell(ic_type) )
then
1017 nn = hecmw_get_max_node( ic_type )
1020 js = hecmesh%elem_node_index(icel-1)
1021 id_area = hecmesh%elem_ID(icel*2)
1023 nodlocal(j) = hecmesh%elem_node_item(js+j)
1024 ecoord(1,j) = hecmesh%node(3*nodlocal(j)-2)
1025 ecoord(2,j) = hecmesh%node(3*nodlocal(j)-1)
1026 ecoord(3,j) = hecmesh%node(3*nodlocal(j) )
1027 edisp(1,j) = fstrsolid%unode(6*nodlocal(j)-5)
1028 edisp(2,j) = fstrsolid%unode(6*nodlocal(j)-4)
1029 edisp(3,j) = fstrsolid%unode(6*nodlocal(j)-3)
1030 edisp(4,j) = fstrsolid%unode(6*nodlocal(j)-2)
1031 edisp(5,j) = fstrsolid%unode(6*nodlocal(j)-1)
1032 edisp(6,j) = fstrsolid%unode(6*nodlocal(j) )
1034 isect = hecmesh%section_ID(icel)
1035 ihead = hecmesh%section%sect_R_index(isect-1)
1036 thick = hecmesh%section%sect_R_item(ihead+1)
1038 if( ic_type == 731 .or. ic_type == 741 .or. ic_type == 743 )
then
1039 ntot_lyr = fstrsolid%elements(icel)%gausses(1)%pMaterial%totallyr
1041 call elementstress_shell_mitc( ic_type, nn, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
1042 & ndstrain(1:nn,1:6), ndstress(1:nn,1:6), thick, 1.0d0, nlyr, ntot_lyr)
1045 layer => fstrsolid%SHELL%LAYER(nlyr)%PLUS
1047 layer%STRAIN(6*(i-1)+k) = layer%STRAIN(6*(i-1)+k) + ndstrain(j,k)
1048 layer%STRESS(6*(i-1)+k) = layer%STRESS(6*(i-1)+k) + ndstress(j,k)
1049 layer%ESTRAIN(6*(icel-1)+k) = layer%ESTRAIN(6*(icel-1)+k) + ndstrain(j,k)/nn
1050 layer%ESTRESS(6*(icel-1)+k) = layer%ESTRESS(6*(icel-1)+k) + ndstress(j,k)/nn
1054 call elementstress_shell_mitc( ic_type, nn, 6, ecoord, fstrsolid%elements(icel)%gausses, edisp, &
1055 & ndstrain(1:nn,1:6), ndstress(1:nn,1:6), thick,-1.0d0, nlyr, ntot_lyr)
1058 layer => fstrsolid%SHELL%LAYER(nlyr)%MINUS
1060 layer%STRAIN(6*(i-1)+k) = layer%STRAIN(6*(i-1)+k) + ndstrain(j,k)
1061 layer%STRESS(6*(i-1)+k) = layer%STRESS(6*(i-1)+k) + ndstress(j,k)
1062 layer%ESTRAIN(6*(icel-1)+k) = layer%ESTRAIN(6*(icel-1)+k) + ndstrain(j,k)/nn
1063 layer%ESTRESS(6*(icel-1)+k) = layer%ESTRESS(6*(icel-1)+k) + ndstress(j,k)/nn
1067 call fstr_getavg_shell(nn,fstrsolid,icel,nodlocal,ndstrain(1:nn,1:6),ndstress(1:nn,1:6),estrain,estress)
1073 ic = hecmesh%elem_node_item(js+j)
1074 fstrsolid%STRAIN(6*(ic-1)+1:6*(ic-1)+6) = fstrsolid%STRAIN(6*(ic-1)+1:6*(ic-1)+6) + ndstrain(j,1:6)
1075 fstrsolid%STRESS(6*(ic-1)+1:6*(ic-1)+6) = fstrsolid%STRESS(6*(ic-1)+1:6*(ic-1)+6) + ndstress(j,1:6)
1079 nnumber(ic) = nnumber(ic) + 1
1082 fstrsolid%ESTRAIN(6*(icel-1)+1:6*(icel-1)+6) = fstrsolid%ESTRAIN(6*(icel-1)+1:6*(icel-1)+6) + estrain(1:6)
1083 fstrsolid%ESTRESS(6*(icel-1)+1:6*(icel-1)+6) = fstrsolid%ESTRESS(6*(icel-1)+1:6*(icel-1)+6) + estress(1:6)
1089 do i = 1, hecmesh%n_node
1090 if( nnumber(i) == 0 ) cycle
1091 fstrsolid%STRAIN(6*(i-1)+1:6*(i-1)+6) = fstrsolid%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1092 fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6) = fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1098 do nlyr = 1, ntot_lyr
1099 do i = 1, hecmesh%n_node
1100 fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRAIN(6*(i-1)+1:6*(i-1)+6) = &
1101 & fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1102 fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6) = &
1103 & fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1104 fstrsolid%SHELL%LAYER(nlyr)%PLUS%MISES(i) = &
1105 &
get_mises(fstrsolid%SHELL%LAYER(nlyr)%PLUS%STRESS(6*(i-1)+1:6*(i-1)+6))
1107 fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRAIN(6*(i-1)+1:6*(i-1)+6) = &
1108 & fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRAIN(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1109 fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6) = &
1110 & fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6) / nnumber(i)
1111 fstrsolid%SHELL%LAYER(nlyr)%MINUS%MISES(i) = &
1112 &
get_mises(fstrsolid%SHELL%LAYER(nlyr)%MINUS%STRESS(6*(i-1)+1:6*(i-1)+6))
1117 do i = 1, hecmesh%n_node
1118 fstrsolid%MISES(i) =
get_mises(fstrsolid%STRESS(6*(i-1)+1:6*(i-1)+6))
1120 do i = 1, hecmesh%n_elem
1121 fstrsolid%EMISES(i) =
get_mises(fstrsolid%ESTRESS(6*(i-1)+1:6*(i-1)+6))
1125 do i = 1, hecmesh%n_elem
1126 if (.not.
associated(fstrsolid%elements(i)%gausses)) cycle
1127 fstrsolid%EPLSTRAIN(i) =
get_pl_estrain(fstrsolid%elements(i)%gausses)
1130 deallocate( nnumber )
1142 integer(kind=kint) :: i, flag
1143 real(kind=
kreal) :: tmat(3, 3), tvec(3), strain(6)
1145 flag=ieor(flag,flag)
1146 if( fstrsolid%output_ctrl(3)%outinfo%on(19) .or. fstrsolid%output_ctrl(4)%outinfo%on(19) )
then
1147 if ( .not.
associated(res%PSTRESS) )
then
1148 allocate(res%PSTRESS( 3*hecmesh%n_node ))
1150 flag=ior(flag,b
'00000001')
1152 if( fstrsolid%output_ctrl(3)%outinfo%on(23) .or. fstrsolid%output_ctrl(4)%outinfo%on(23) )
then
1153 if ( .not.
associated(res%PSTRESS_VECT) )
then
1154 allocate(res%PSTRESS_VECT( 3*hecmesh%n_node ,3))
1156 flag=ior(flag,b
'00000010')
1158 if( fstrsolid%output_ctrl(3)%outinfo%on(21) .or. fstrsolid%output_ctrl(4)%outinfo%on(21) )
then
1159 if ( .not.
associated(res%PSTRAIN) )
then
1160 allocate(res%PSTRAIN( 3*hecmesh%n_node ))
1162 flag=ior(flag,b
'00000100')
1164 if( fstrsolid%output_ctrl(3)%outinfo%on(25) .or. fstrsolid%output_ctrl(4)%outinfo%on(25) )
then
1165 if ( .not.
associated(res%PSTRAIN_VECT) )
then
1166 allocate(res%PSTRAIN_VECT( 3*hecmesh%n_node ,3))
1168 flag=ior(flag,b
'00001000')
1170 if( fstrsolid%output_ctrl(3)%outinfo%on(20) .or. fstrsolid%output_ctrl(4)%outinfo%on(20) )
then
1171 if ( .not.
associated(res%EPSTRESS) )
then
1172 allocate(res%EPSTRESS( 3*hecmesh%n_elem ))
1174 flag=ior(flag,b
'00010000')
1176 if( fstrsolid%output_ctrl(3)%outinfo%on(24) .or. fstrsolid%output_ctrl(4)%outinfo%on(24) )
then
1177 if ( .not.
associated(res%EPSTRESS_VECT) )
then
1178 allocate(res%EPSTRESS_VECT( 3*hecmesh%n_elem ,3))
1180 flag=ior(flag,b
'00100000')
1182 if( fstrsolid%output_ctrl(3)%outinfo%on(22) .or. fstrsolid%output_ctrl(4)%outinfo%on(22) )
then
1183 if ( .not.
associated(res%EPSTRAIN) )
then
1184 allocate(res%EPSTRAIN( 3*hecmesh%n_elem ))
1186 flag=ior(flag,b
'01000000')
1188 if( fstrsolid%output_ctrl(3)%outinfo%on(26) .or. fstrsolid%output_ctrl(4)%outinfo%on(26) )
then
1189 if ( .not.
associated(res%EPSTRAIN_VECT) )
then
1190 allocate(res%EPSTRAIN_VECT( 3*hecmesh%n_elem ,3))
1192 flag=ior(flag,b
'10000000')
1195 if (iand(flag,b
'00000011') /= 0)
then
1196 do i = 1, hecmesh%n_node
1197 call get_principal(res%STRESS(6*i-5:6*i), tvec, tmat)
1198 if (iand(flag,b
'00000001') /= 0) res%PSTRESS(3*(i-1)+1:3*(i-1)+3)=tvec
1199 if (iand(flag,b
'00000010') /= 0) res%PSTRESS_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat
1202 if (iand(flag,b
'00001100') /= 0)
then
1203 do i = 1, hecmesh%n_node
1204 strain(1:6) = res%STRAIN(6*i-5:6*i)
1205 strain(4:6) = 0.5d0*strain(4:6)
1206 call get_principal(strain, tvec, tmat)
1207 if (iand(flag,b
'00000100') /= 0) res%PSTRAIN(3*(i-1)+1:3*(i-1)+3)=tvec
1208 if (iand(flag,b
'00001000') /= 0) res%PSTRAIN_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat
1212 if (iand(flag,b
'00110000') /= 0)
then
1213 do i = 1, hecmesh%n_elem
1214 call get_principal( res%ESTRESS(6*i-5:6*i), tvec, tmat)
1215 if (iand(flag,b
'00010000') /= 0) res%EPSTRESS(3*(i-1)+1:3*(i-1)+3)=tvec
1216 if (iand(flag,b
'00100000') /= 0) res%EPSTRESS_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat
1219 if (iand(flag,b
'11000000') /= 0)
then
1220 do i = 1, hecmesh%n_elem
1221 strain(1:6) = res%ESTRAIN(6*i-5:6*i)
1222 strain(4:6) = 0.5d0*strain(4:6)
1223 call get_principal(strain, tvec, tmat)
1224 if (iand(flag,b
'01000000') /= 0) res%EPSTRAIN(3*(i-1)+1:3*(i-1)+3)=tvec
1225 if (iand(flag,b
'10000000') /= 0) res%EPSTRAIN_VECT(3*(i-1)+1:3*(i-1)+3,1:3)=tmat
integer(kind=4), parameter kreal
This module provides functions to calculation nodal stress.
subroutine fstr_stress_add_shelllyr(nn, fstrSOLID, icel, nodLOCAL, nlyr, strain, stress, flag)
real(kind=kreal) function get_pl_estrain(gausses)
subroutine fstr_nodalstress3d(hecMESH, fstrSOLID)
Calculate NODAL STRESS of solid elements.
subroutine fstr_nodalstress6d(hecMESH, fstrSOLID)
Calculate NODAL STRESS of shell elements.
real(kind=kreal) function get_mises(s)
subroutine fstr_nodalstress3d_c3d4_sesns(hecMESH, fstrSOLID, nnumber, Nodal_STRAIN, Nodal_STRESS, Elemental_STRAIN, Elemental_STRESS)
subroutine fstr_nodalstress2d(hecMESH, fstrSOLID)
Calculate NODAL STRESS of plane elements.
subroutine make_principal(fstrSOLID, hecMESH, RES)
integer(kind=kint) function search_idx_senes(irow, asect, nid, sid)
subroutine fstr_getavg_shell(nn, fstrSOLID, icel, nodLOCAL, strain, stress, estrain, estress)
This module defines common data and basic structures for analysis.
integer(kind=kint), parameter kel341sesns
subroutine fstr_solid_phys_clear(fstrSOLID)
This module manages step information.
This modules just summarizes all modules used in static analysis.
This modules defines a structure to record history dependent parameter in static analysis.
Data for STATIC ANSLYSIS (fstrSOLID)
All data should be recorded in every quadrature points.