FrontISTR  5.8.0
Large-scale structural analysis program with finit element method
fstr_Update.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
7  use m_fstr
8  implicit none
9 
10  private :: update_abort
11 
12 contains
13 
14  !=====================================================================*
15  ! UPDATE_C3
25  subroutine fstr_updatenewton ( hecMESH, hecMAT, fstrSOLID, time, tincr,iter, strainEnergy)
26  !=====================================================================*
27  use m_static_lib
28  use m_elemact
29 
30  type (hecmwST_matrix) :: hecMAT
31  type (hecmwST_local_mesh) :: hecMESH
32  type (fstr_solid) :: fstrSOLID
33  real(kind=kreal),intent(in) :: time
34  real(kind=kreal),intent(in) :: tincr
35  integer, intent(in) :: iter
36 
37  integer(kind=kint) :: nodLOCAL(fstrSOLID%max_ncon)
38  real(kind=kreal) :: ecoord(3, fstrsolid%max_ncon)
39  real(kind=kreal) :: thick, thick0(6)
40  integer(kind=kint) :: ndof, itype, is, iE, ic_type, nn, icel, iiS, i, j
41 
42  real(kind=kreal) :: total_disp(6, fstrsolid%max_ncon), du(6, fstrsolid%max_ncon), ddu(6, fstrsolid%max_ncon)
43  real(kind=kreal) :: tt(fstrsolid%max_ncon), tt0(fstrsolid%max_ncon), ttn(fstrsolid%max_ncon)
44  real(kind=kreal) :: qf(fstrsolid%max_ncon*6), coords(3, 3)
45  integer :: isect, ihead, cdsys_ID
46  integer :: ndim, initt
47 
48  real(kind=kreal), optional :: strainenergy
49  real(kind=kreal) :: tmp
50  real(kind=kreal) :: ddaux(3,3)
51 
52  ndof = hecmat%NDOF
53  fstrsolid%QFORCE=0.0d0
54 
55  tt0 = 0.d0
56  ttn = 0.d0
57  tt = 0.d0
58 
59  ! if initial temperature exists
60  initt = 0
61  if( associated(g_initialcnd) ) then
62  do j=1,size(g_initialcnd)
63  if( g_initialcnd(j)%cond_name=="temperature" ) then
64  initt=j
65  exit
66  endif
67  end do
68  endif
69 
70  ! --------------------------------------------------------------------
71  ! updated
72  ! 1. stress and strain : ep^(k) = ep^(k-1)+dep^(k)
73  ! sgm^(k) = sgm^(k-1)+dsgm^(k)
74  ! 2. Internal Force : Q^(k-1) ( u^(k-1) )
75  ! --------------------------------------------------------------------
76  ! ----------------------------------------------------------------------------------
77  ! calculate the Strain and Stress and Internal Force ( Equivalent Nodal Force )
78  ! ----------------------------------------------------------------------------------
79 
80  do itype = 1, hecmesh%n_elem_type
81  is = hecmesh%elem_type_index(itype-1)+1
82  ie = hecmesh%elem_type_index(itype )
83  ic_type= hecmesh%elem_type_item(itype)
84  if (hecmw_is_etype_link(ic_type)) cycle
85  if (hecmw_is_etype_patch(ic_type)) cycle
86 
87  !element loop
88  !$omp parallel default(none), &
89  !$omp& private(icel,iiS,j,nn,nodLOCAL,i,ecoord,ddu,du,total_disp, &
90  !$omp& cdsys_ID,coords,thick,qf,isect,ihead,tmp,ndim,ddaux,thick0), &
91  !$omp& shared(iS,iE,hecMESH,fstrSOLID,ndof,hecMAT,ic_type,fstrPR, &
92  !$omp& strainEnergy,iter,time,tincr,initt,g_InitialCnd), &
93  !$omp& firstprivate(tt0,ttn,tt)
94  !$omp do
95  do icel = is, ie
96 
97  ! ----- nodal coordinate, displacement and temperature
98  iis = hecmesh%elem_node_index(icel-1)
99  nn = hecmesh%elem_node_index(icel)-iis
100  !if( nn>150 ) stop "elemental nodes > 150!"
101 
102  thick = 0.d0
103  do j = 1, size(fstrsolid%elements(icel)%gausses)
104  thick0(1:6) = fstrsolid%elements(icel)%gausses(j)%stress(1:6)
105  thick = thick + dsqrt(dot_product(thick0(1:6),thick0(1:6)))
106  enddo
107  if( thick < 1.d-10 ) then
108  do j = 1, size(fstrsolid%elements(icel)%gausses)
109  if( associated(fstrsolid%elements(icel)%gausses(j)%fstatus) ) &
110  & fstrsolid%elements(icel)%gausses(j)%fstatus = 0.d0
111  enddo
112  end if
113 
114  do j = 1, nn
115  nodlocal(j) = hecmesh%elem_node_item (iis+j)
116  do i = 1, 3
117  ecoord(i,j) = hecmesh%node(3*nodlocal(j)+i-3)
118  enddo
119  do i = 1, ndof
120  ddu(i,j) = hecmat%X(ndof*nodlocal(j)+i-ndof)
121  du(i,j) = fstrsolid%dunode(ndof*nodlocal(j)+i-ndof)
122  total_disp(i,j) = fstrsolid%unode(ndof*nodlocal(j)+i-ndof)
123  enddo
124 
125  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
126  if( iselastoplastic(fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype) .or. &
127  fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype == norton ) then
128  tt0(j)=fstrsolid%last_temp( nodlocal(j) )
129  else
130  tt0(j) = 0.d0
131  if( hecmesh%hecmw_flag_initcon == 1 ) tt0(j) = hecmesh%node_init_val_item(nodlocal(j))
132  if( initt>0 ) tt0(j) = g_initialcnd(initt)%realval(nodlocal(j))
133  endif
134  ttn(j) = fstrsolid%last_temp( nodlocal(j) )
135  tt(j) = fstrsolid%temperature( nodlocal(j) )
136  endif
137  enddo
138 
139  isect = hecmesh%section_ID(icel)
140  ihead = hecmesh%section%sect_R_index(isect-1)
141  thick = hecmesh%section%sect_R_item(ihead+1)
142  cdsys_id = hecmesh%section%sect_orien_ID(isect)
143  if( cdsys_id > 0 ) call get_coordsys(cdsys_id, hecmesh, fstrsolid, coords)
144 
145  ! ===== calculate the Internal Force
146  if( ic_type == 241 .or. ic_type == 242 .or. ic_type == 231 .or. ic_type == 232 .or. ic_type == 2322 ) then
147  call update_c2( ic_type,nn,ecoord(1:3,1:nn),fstrsolid%elements(icel)%gausses(:), &
148  thick,fstrsolid%elements(icel)%iset, &
149  total_disp(1:2,1:nn), ddu(1:2,1:nn), qf(1:nn*ndof), &
150  tt(1:nn), tt0(1:nn), ttn(1:nn) )
151 
152  else if( ic_type == 301 ) then
153  call update_c1( ic_type,nn,ecoord(:,1:nn), thick, total_disp(1:3,1:nn), du(1:3,1:nn), &
154  qf(1:nn*ndof),fstrsolid%elements(icel)%gausses(:) )
155 
156  else if( ic_type == 361 ) then
157  if( fstrsolid%sections(isect)%elemopt361 == kel361fi ) then ! full integration element
158  call update_c3( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
159  qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn) )
160  else if( fstrsolid%sections(isect)%elemopt361 == kel361bbar ) then ! B-bar element
161  call update_c3d8bbar( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
162  qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn) )
163  else if( fstrsolid%sections(isect)%elemopt361 == kel361ic ) then ! incompatible element
164  call update_c3d8ic( ic_type,nn,ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), ddu(1:3,1:nn), cdsys_id, coords,&
165  qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, &
166  fstrsolid%elements(icel)%aux, ddaux(1:3,1:3), tt(1:nn), tt0(1:nn), ttn(1:nn) )
167  fstrsolid%elements(icel)%aux(1:3,1:3) = fstrsolid%elements(icel)%aux(1:3,1:3) + ddaux(1:3,1:3)
168  else if( fstrsolid%sections(isect)%elemopt361 == kel361fbar ) then ! F-bar element
169  call update_c3d8fbar( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
170  qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn) )
171  endif
172 
173  else if (ic_type == 341 .or. ic_type == 351 .or. ic_type == 342 .or. ic_type == 352 .or. ic_type == 362 ) then
174  if( ic_type==341 .and. fstrsolid%sections(isect)%elemopt341 == kel341sesns ) cycle ! skip smoothed fem
175  call update_c3( ic_type, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
176  qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), iter, time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn) )
177 
178  else if( ic_type == 511) then
179  call update_connector( ic_type,nn,ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), &
180  qf(1:nn*ndof),fstrsolid%elements(icel)%gausses(:), tincr )
181 
182  else if( ic_type == 611) then
183  if( fstrpr%nlgeom ) call update_abort( ic_type, 2 )
184  CALL updatest_beam(ic_type, nn, ecoord, total_disp(1:6,1:nn), du(1:6,1:nn), &
185  & hecmesh%section%sect_R_item(ihead+1:), fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof))
186 
187  else if( ic_type == 641 ) then
188  if( fstrpr%nlgeom ) call update_abort( ic_type, 2 )
189  call updatest_beam_641(ic_type, nn, ecoord, total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
190  & fstrsolid%elements(icel)%gausses(:), hecmesh%section%sect_R_item(ihead+1:), qf(1:nn*ndof))
191 
192  else if( ( ic_type == 741 ) .or. ( ic_type == 743 ) .or. ( ic_type == 731 ) ) then
193  if( fstrpr%nlgeom ) call update_abort( ic_type, 2 )
194  call updatest_shell_mitc(ic_type, nn, ndof, ecoord(1:3, 1:nn), total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
195  & fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof), thick, 0)
196 
197  else if( ic_type == 761 ) then !for shell-solid mixed analysis
198  if( fstrpr%nlgeom ) call update_abort( ic_type, 2 )
199  call updatest_shell_mitc33(731, 3, 6, ecoord(1:3, 1:3), total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
200  & fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof), thick, 2)
201 
202  else if( ic_type == 781 ) then !for shell-solid mixed analysis
203  if( fstrpr%nlgeom ) call update_abort( ic_type, 2 )
204  call updatest_shell_mitc33(741, 4, 6, ecoord(1:3, 1:4), total_disp(1:ndof,1:nn), du(1:ndof,1:nn), &
205  & fstrsolid%elements(icel)%gausses(:), qf(1:nn*ndof), thick, 1)
206 
207  else if ( ic_type == 3414 ) then
208  if(fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype /= incomp_newtonian) &
209  & call update_abort( ic_type, 3, fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype )
210  call update_c3_vp &
211  ( ic_type, nn, ecoord(:,1:nn), total_disp(1:4,1:nn), du(1:4,1:nn), &
212  fstrsolid%elements(icel)%gausses(:) )
213  qf = 0.0d0
214 
215  else if ( ic_type == 881 .or. ic_type == 891 ) then !for selective es/ns smoothed fem
216  call update_c3_sesns( ic_type, nn, nodlocal, ecoord(:,1:nn), total_disp(1:3,1:nn), du(1:3,1:nn), cdsys_id, coords, &
217  qf(1:nn*ndof), fstrsolid%elements(icel)%gausses(:), time, tincr, tt(1:nn), tt0(1:nn), ttn(1:nn) )
218 
219  else
220  write(*, *) '###ERROR### : Element type not supported for nonlinear static analysis'
221  write(*, *) ' ic_type = ', ic_type
222  call hecmw_abort(hecmw_comm_get_comm())
223 
224  endif
225 
226  ! elemact element
227  if( fstrsolid%elements(icel)%elemact_flag == kelact_inactive ) then
228  call update_dummy( ndof, nn, ecoord(:,1:nn), total_disp(1:3,1:nn), &
229  & du(1:3,1:nn), qf(1:nn*ndof), fstrsolid%elements(icel) )
230  !qf(:) = fstrSOLID%elements(icel)%elemact_coeff*qf(:)
231  end if
232 
233  ! ----- calculate the global internal force ( Q(u_{n+1}^{k-1}) )
234  do j = 1, nn
235  do i = 1, ndof
236  !$omp atomic
237  fstrsolid%QFORCE(ndof*(nodlocal(j)-1)+i) = fstrsolid%QFORCE(ndof*(nodlocal(j)-1)+i)+qf(ndof*(j-1)+i)
238  enddo
239  enddo
240 
241  ! ----- calculate strain energy
242  if(present(strainenergy))then
243  ndim = getspacedimension( fstrsolid%elements(icel)%etype )
244  do j = 1, nn
245  do i = 1, ndim
246  tmp = 0.5d0*( fstrsolid%elements(icel)%equiForces(ndim*(j-1)+i)+qf(ndim*(j-1)+i) )*ddu(i,j)
247  !$omp atomic
248  strainenergy = strainenergy+tmp
249  fstrsolid%elements(icel)%equiForces(ndim*(j-1)+i) = qf(ndim*(j-1)+i)
250  enddo
251  enddo
252  endif
253 
254  enddo ! icel
255  !$omp end do
256  !$omp end parallel
257  enddo ! itype
258 
259  !C
260  !C Update for fstrSOLID%QFORCE
261  !C
262  call hecmw_update_r(hecmesh,fstrsolid%QFORCE,hecmesh%n_node, ndof)
263  end subroutine fstr_updatenewton
264 
265 
267  subroutine fstr_updatestate( hecMESH, fstrSOLID, tincr)
268  use m_fstr
269  use m_static_lib
270  use m_elastoplastic
271  use mcreep
272  use mviscoelastic
273  type(hecmwst_local_mesh) :: hecmesh
274  type(fstr_solid) :: fstrSOLID
275  real(kind=kreal) :: tincr
276  integer(kind=kint) :: itype, is, iE, ic_type, icel, ngauss, i
277 
278  if( associated( fstrsolid%temperature ) ) then
279  do i = 1, hecmesh%n_node
280  fstrsolid%last_temp(i) = fstrsolid%temperature(i)
281  end do
282  endif
283 
284  do itype = 1, hecmesh%n_elem_type
285  is = hecmesh%elem_type_index(itype-1) + 1
286  ie = hecmesh%elem_type_index(itype )
287  ic_type= hecmesh%elem_type_item(itype)
288  if( ic_type == 301 ) ic_type = 111
289  if( hecmw_is_etype_link(ic_type) ) cycle
290  if( hecmw_is_etype_patch(ic_type) ) cycle
291 
292  ngauss = numofquadpoints( ic_type )
293  do icel = is, ie
294  if( iselastoplastic( fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype ) ) then
295  do i = 1, ngauss
296  call updateepstate( fstrsolid%elements(icel)%gausses(i) )
297  enddo
298  elseif( fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype == norton ) then
299  if( tincr>0.d0 ) then
300  do i = 1, ngauss
301  call updateviscostate( fstrsolid%elements(icel)%gausses(i) )
302  enddo
303  endif
304  elseif( isviscoelastic( fstrsolid%elements(icel)%gausses(1)%pMaterial%mtype ) ) then
305  if( tincr > 0.d0 ) then
306  do i = 1, ngauss
307  call updateviscoelasticstate( fstrsolid%elements(icel)%gausses(i) )
308  enddo
309  endif
310  endif
311 
312  do i = 1, ngauss
313  fstrsolid%elements(icel)%gausses(i)%strain_bak = fstrsolid%elements(icel)%gausses(i)%strain
314  fstrsolid%elements(icel)%gausses(i)%stress_bak = fstrsolid%elements(icel)%gausses(i)%stress
315  fstrsolid%elements(icel)%gausses(i)%strain_energy_bak = fstrsolid%elements(icel)%gausses(i)%strain_energy
316  enddo
317  enddo
318  enddo
319 
320  do i = 1, hecmesh%n_node
321  fstrsolid%QFORCE_bak(i) = fstrsolid%QFORCE(i)
322  end do
323 
324  end subroutine fstr_updatestate
325 
326  subroutine update_abort( ic_type, flag, mtype )
327  integer(kind=kint), intent(in) :: ic_type
328  integer(kind=kint), intent(in) :: flag
329  integer(kind=kint), intent(in), optional :: mtype
330 
331  if( flag == 1 ) then
332  write(*,*) '###ERROR### : Element type not supported for static analysis'
333  else if( flag == 2 ) then
334  write(*,*) '###ERROR### : Element type not supported for nonlinear static analysis'
335  else if( flag == 3 ) then
336  write(*,*) '###ERROR### : This element is not supported for this material'
337  endif
338  write(*,*) ' ic_type = ', ic_type
339  if( present(mtype) ) write(*,*) ' mtype = ', mtype
340  call hecmw_abort(hecmw_comm_get_comm())
341  end subroutine
342 
343 end module m_fstr_update
This module provide functions for elastoplastic calculation.
subroutine, public updateepstate(gauss)
Clear elatoplastic state.
This module defined elemact data and function.
subroutine update_dummy(ndof, nn, ecoord, u, du, qf, element)
integer, parameter kelact_inactive
This module provides function to calculate to do updates.
Definition: fstr_Update.f90:6
subroutine fstr_updatestate(hecMESH, fstrSOLID, tincr)
Update elastiplastic status.
subroutine fstr_updatenewton(hecMESH, hecMAT, fstrSOLID, time, tincr, iter, strainEnergy)
Update displacement, stress, strain and internal forces.
Definition: fstr_Update.f90:26
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
integer(kind=kint), parameter kel361bbar
Definition: m_fstr.f90:79
subroutine get_coordsys(cdsys_ID, hecMESH, fstrSOLID, coords)
This subroutine fetch coords defined by local coordinate system.
Definition: m_fstr.f90:1093
integer(kind=kint), parameter kel341sesns
Definition: m_fstr.f90:76
integer(kind=kint), parameter kel361fi
Definition: m_fstr.f90:78
integer(kind=kint), parameter kel361ic
Definition: m_fstr.f90:80
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
Definition: m_fstr.f90:210
integer(kind=kint), parameter kel361fbar
Definition: m_fstr.f90:81
type(tinitialcondition), dimension(:), pointer, save g_initialcnd
Definition: m_fstr.f90:152
This modules just summarizes all modules used in static analysis.
Definition: static_LIB.f90:6
This module provides functions for creep calculation.
Definition: creep.f90:6
subroutine updateviscostate(gauss)
Update viscoplastic state.
Definition: creep.f90:214
This module provides functions for viscoelastic calculation.
Definition: Viscoelastic.f90:6
subroutine updateviscoelasticstate(gauss)
Update viscoplastic state.