FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_dynamic_nlimplicit.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 !-------------------------------------------------------------------------------
6 
8  use m_fstr
16  use m_fstr_update
17  use m_fstr_restart
19  use m_fstr_residual
20  use mcontact
23 
24  !-------- for couple -------
26  use m_fstr_rcap_io
27 
28 contains
29 
31 
32  subroutine fstr_solve_dynamic_nlimplicit(cstep, hecMESH,hecMAT,fstrSOLID,fstrEIG, &
33  fstrDYNAMIC,fstrRESULT,fstrPARAM,fstrCPL, restrt_step_num )
34  implicit none
35  !C-- global variable
36  integer, intent(in) :: cstep
37  type(hecmwst_local_mesh) :: hecMESH
38  type(hecmwst_matrix) :: hecMAT
39  type(fstr_eigen) :: fstrEIG
40  type(fstr_solid) :: fstrSOLID
41  type(hecmwst_result_data) :: fstrRESULT
42  type(fstr_param) :: fstrPARAM
43  type(fstr_dynamic) :: fstrDYNAMIC
44  type(hecmwst_matrix_lagrange) :: hecLagMAT
45  type(fstr_couple) :: fstrCPL !for COUPLE
46 
47  !C-- local variable
48  type(hecmwst_local_mesh), pointer :: hecMESHmpc
49  type(hecmwst_matrix), pointer :: hecMATmpc
50  type(hecmwst_matrix), pointer :: hecMAT0
51  integer(kind=kint) :: nnod, ndof, numnp, nn
52  integer(kind=kint) :: i, j, ids, ide, ims, ime, kk, idm, imm
53  integer(kind=kint) :: iter
54  integer(kind=kint) :: iiii5, iexit
55  integer(kind=kint) :: revocap_flag
56  integer(kind=kint) :: kkk0, kkk1
57  integer(kind=kint) :: restrt_step_num
58  integer(kind=kint) :: n_node_global
59  integer(kind=kint) :: ierr
60 
61  real(kind=kreal) :: a1, a2, a3, b1, b2, b3, c1, c2
62  real(kind=kreal) :: bsize, res, resb
63  real(kind=kreal) :: time_1, time_2
64  real(kind=kreal), parameter :: pi = 3.14159265358979323846d0
65  real(kind=kreal), allocatable :: coord(:)
66 
67  logical :: is_cycle
68 
69  iexit = 0
70  resb = 0.0d0
71 
72  call hecmw_mpc_mat_init(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
73  nullify(hecmat0)
74 
75  ! sum of n_node among all subdomains (to be used to calc res)
76  n_node_global = hecmesh%nn_internal
77  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
78 
79  hecmat%NDOF=hecmesh%n_dof
80  nnod=hecmesh%n_node
81  ndof=hecmat%NDOF
82  nn =ndof*ndof
83 
84  allocate(coord(hecmesh%n_node*ndof))
85 
86  !!-- initial value
87  time_1 = hecmw_wtime()
88 
89  !C-- check parameters
90  if(dabs(fstrdynamic%beta) < 1.0e-20) then
91  if( hecmesh%my_rank == 0 ) then
92  write(imsg,*) 'stop due to Newmark-beta = 0'
93  endif
94  call hecmw_abort( hecmw_comm_get_comm())
95  endif
96 
97  !C-- matrix [M] lumped mass matrix
98  if(fstrdynamic%idx_mas == 1) then
99  call setmass(fstrsolid,hecmesh,hecmat,fstreig)
100 
101  !C-- consistent mass matrix
102  else if(fstrdynamic%idx_mas == 2) then
103  if( hecmesh%my_rank .eq. 0 ) then
104  write(imsg,*) 'stop: consistent mass matrix is not yet available !'
105  endif
106  call hecmw_abort( hecmw_comm_get_comm())
107  endif
108 
109  hecmat%Iarray(98) = 1 !Assembly complete
110  hecmat%Iarray(97) = 1 !Need numerical factorization
111 
112  !C-- time step loop
113  a1 = 0.5d0/fstrdynamic%beta - 1.0d0
114  a2 = 1.0d0/(fstrdynamic%beta*fstrdynamic%t_delta)
115  a3 = 1.0d0/(fstrdynamic%beta*fstrdynamic%t_delta*fstrdynamic%t_delta)
116  b1 = (0.5d0*fstrdynamic%gamma/fstrdynamic%beta - 1.0d0 )*fstrdynamic%t_delta
117  b2 = fstrdynamic%gamma/fstrdynamic%beta - 1.0d0
118  b3 = fstrdynamic%gamma/(fstrdynamic%beta*fstrdynamic%t_delta)
119  c1 = 1.0d0 + fstrdynamic%ray_k*b3
120  c2 = a3 + fstrdynamic%ray_m*b3
121 
122  !C-- output of initial state
123  if( restrt_step_num == 1 ) then
124  call fstr_dynamic_output(hecmesh, fstrsolid, fstrdynamic, fstrparam)
125  call dynamic_output_monit(hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
126  endif
127 
128  fstrdynamic%VEC3(:) =0.0d0
129  hecmat%X(:) =0.0d0
130 
131  !! step = 1,2,....,fstrDYNAMIC%n_step
132  do i = restrt_step_num, fstrdynamic%n_step
133  if(ndof == 4 .and. hecmesh%my_rank==0) write(*,'(a,i5)')"iter: ",i
134 
135  fstrdynamic%i_step = i
136  fstrdynamic%t_curr = fstrdynamic%t_delta * i
137 
138  if(hecmesh%my_rank==0) then
139  !write(*,'('' time step='',i10,'' time='',1pe13.4e3)') i,fstrDYNAMIC%t_curr
140  write(ista,'('' time step='',i10,'' time='',1pe13.4e3)') i,fstrdynamic%t_curr
141  endif
142 
143  do j = 1 ,ndof*nnod
144  fstrdynamic%VEC1(j) = a1*fstrdynamic%ACC(j,1) + a2*fstrdynamic%VEL(j,1)
145  fstrdynamic%VEC2(j) = b1*fstrdynamic%ACC(j,1) + b2*fstrdynamic%VEL(j,1)
146  enddo
147 
148  !C ********************************************************************************
149  !C for couple analysis
150  do
151  fstrsolid%dunode(:) =0.d0
152  ! call fstr_UpdateEPState( hecMESH, fstrSOLID )
153  call fstr_solve_dynamic_nlimplicit_couple_init(fstrparam, fstrcpl)
154 
155  do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
156  if (fstrparam%nlgeom) then
157  call fstr_stiffmatrix( hecmesh, hecmat, fstrsolid, fstrdynamic%t_curr, fstrdynamic%t_delta )
158  else
159  if (.not. associated(hecmat0)) then
160  call fstr_stiffmatrix( hecmesh, hecmat, fstrsolid, fstrdynamic%t_curr, fstrdynamic%t_delta )
161  allocate(hecmat0)
162  call hecmw_mat_init(hecmat0)
163  call hecmw_mat_copy_profile(hecmat, hecmat0)
164  call hecmw_mat_copy_val(hecmat, hecmat0)
165  else
166  call hecmw_mat_copy_val(hecmat0, hecmat)
167  endif
168  endif
169 
170  if( fstrdynamic%ray_k/=0.d0 .or. fstrdynamic%ray_m/=0.d0 ) then
171  do j = 1 ,ndof*nnod
172  hecmat%X(j) = fstrdynamic%VEC2(j) - b3*fstrsolid%dunode(j)
173  enddo
174  endif
175  if( fstrdynamic%ray_k/=0.d0 ) then
176  if( hecmesh%n_dof == 3 ) then
177  call hecmw_matvec (hecmesh, hecmat, hecmat%X, fstrdynamic%VEC3)
178  else if( hecmesh%n_dof == 2 ) then
179  call hecmw_matvec (hecmesh, hecmat, hecmat%X, fstrdynamic%VEC3)
180  else if( hecmesh%n_dof == 6 ) then
181  call matvec(fstrdynamic%VEC3, hecmat%X, hecmat, ndof, hecmat%D, hecmat%AU, hecmat%AL)
182  endif
183  endif
184 
185  !C-- mechanical boundary condition
186  call dynamic_mat_ass_load (hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, iter)
187  do j=1, hecmesh%n_node* hecmesh%n_dof
188  hecmat%B(j) = hecmat%B(j)- fstrsolid%QFORCE(j) + fstreig%mass(j)*( fstrdynamic%VEC1(j)-a3*fstrsolid%dunode(j) &
189  + fstrdynamic%ray_m* hecmat%X(j) ) + fstrdynamic%ray_k*fstrdynamic%VEC3(j)
190  enddo
191 
192  !C ********************************************************************************
193  !C for couple analysis
194  call fstr_solve_dynamic_nlimplicit_couple_pre(hecmesh, hecmat, fstrsolid, &
195  & fstrparam, fstrdynamic, fstrcpl, restrt_step_num, pi, i)
196 
197  do j = 1 ,nn*hecmat%NP
198  hecmat%D(j) = c1* hecmat%D(j)
199  enddo
200  do j = 1 ,nn*hecmat%NPU
201  hecmat%AU(j) = c1* hecmat%AU(j)
202  enddo
203  do j = 1 ,nn*hecmat%NPL
204  hecmat%AL(j) = c1*hecmat%AL(j)
205  enddo
206  do j=1,nnod
207  do kk=1,ndof
208  idm = nn*(j-1)+1 + (ndof+1)*(kk-1)
209  imm = ndof*(j-1) + kk
210  hecmat%D(idm) = hecmat%D(idm) + c2*fstreig%mass(imm)
211  enddo
212  enddo
213 
214  !C-- geometrical boundary condition
215  call dynamic_mat_ass_bc (hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, heclagmat, iter)
216  call dynamic_mat_ass_bc_vl(hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, heclagmat, iter)
217  call dynamic_mat_ass_bc_ac(hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, heclagmat, iter)
218  call hecmw_mpc_mat_ass(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
219  call hecmw_mpc_trans_rhs(hecmesh, hecmat, hecmatmpc)
220 
221  !C-- RHS LOAD VECTOR CHECK
222  numnp=hecmatmpc%NP
223  call hecmw_innerproduct_r(hecmesh, ndof, hecmatmpc%B, hecmatmpc%B, bsize)
224 
225  if(iter == 1)then
226  resb = bsize
227  endif
228 
229  !res = dsqrt(bsize)/n_node_global
230  res = dsqrt(bsize/resb)
231  if( fstrparam%nlgeom .and. ndof /= 4 ) then
232  if(hecmesh%my_rank==0) write(*,'(a,i5,a,1pe12.4)')"iter: ",iter,", res: ",res
233  if(hecmesh%my_rank==0) write(ista,'(''iter='',I5,''- Residual'',E15.7)')iter,res
234  if( res<fstrsolid%step_ctrl(cstep)%converg ) exit
235  endif
236 
237  !C-- linear solver [A]{X} = {B}
238  hecmatmpc%X = 0.0d0
239  if( iexit .ne. 1 ) then
240  if( fstrparam%nlgeom ) then
241  if( iter == 1 ) then
242  hecmatmpc%Iarray(97) = 2 !Force numerical factorization
243  else
244  hecmatmpc%Iarray(97) = 1 !Need numerical factorization
245  endif
246  call fstr_set_current_config_to_mesh(hecmeshmpc,fstrsolid,coord)
247  endif
248  call solve_lineq(hecmeshmpc,hecmatmpc)
249  if( fstrparam%nlgeom ) then
250  call fstr_recover_initial_config_to_mesh(hecmeshmpc,fstrsolid,coord)
251  endif
252  endif
253  call hecmw_mpc_tback_sol(hecmesh, hecmat, hecmatmpc)
254 
255  do j=1,hecmesh%n_node*ndof
256  fstrsolid%dunode(j) = fstrsolid%dunode(j)+hecmat%X(j)
257  enddo
258  ! ----- update the strain, stress, and internal force
259  call fstr_updatenewton( hecmesh, hecmat, fstrsolid, fstrdynamic%t_curr, &
260  & fstrdynamic%t_delta, iter, fstrdynamic%strainEnergy )
261 
262  if(.not. fstrparam%nlgeom) exit
263  if(ndof == 4) exit
264  enddo
265 
266  ! ----- not convergence
267  if( iter>fstrsolid%step_ctrl(cstep)%max_iter ) then
268  if( hecmesh%my_rank==0) then
269  write(ilog,*) '### Fail to Converge : at step=', i
270  write(ista,*) '### Fail to Converge : at step=', i
271  write( *,*) ' ### Fail to Converge : at step=', i
272  endif
273  stop
274  endif
275 
276  !C *****************************************************
277  !C for couple analysis
278  call fstr_solve_dynamic_nlimplicit_couple_post(hecmesh, hecmat, fstrsolid, &
279  & fstrparam, fstrdynamic, fstrcpl, a1, a2, a3, b1, b2, b3, i, is_cycle)
280  if(is_cycle) cycle
281  exit
282  enddo
283  !C *****************************************************
284 
285  !C-- new displacement, velocity and acceleration
286  fstrdynamic%kineticEnergy = 0.0d0
287  do j = 1 ,ndof*nnod
288  fstrdynamic%ACC (j,2) = -a1*fstrdynamic%ACC(j,1) - a2*fstrdynamic%VEL(j,1) + &
289  a3*fstrsolid%dunode(j)
290  fstrdynamic%VEL (j,2) = -b1*fstrdynamic%ACC(j,1) - b2*fstrdynamic%VEL(j,1) + &
291  b3*fstrsolid%dunode(j)
292  fstrdynamic%ACC (j,1) = fstrdynamic%ACC (j,2)
293  fstrdynamic%VEL (j,1) = fstrdynamic%VEL (j,2)
294 
295  fstrsolid%unode(j) = fstrsolid%unode(j)+fstrsolid%dunode(j)
296  fstrdynamic%DISP(j,2) = fstrsolid%unode(j)
297 
298  fstrdynamic%kineticEnergy = fstrdynamic%kineticEnergy + &
299  0.5d0*fstreig%mass(j)*fstrdynamic%VEL(j,2)*fstrdynamic%VEL(j,2)
300  enddo
301 
302  !C-- output new displacement, velocity and acceleration
303  call fstr_dynamic_output(hecmesh, fstrsolid, fstrdynamic, fstrparam)
304 
305  !C-- output result of monitoring node
306  call dynamic_output_monit(hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
307  call fstr_updatestate( hecmesh, fstrsolid, fstrdynamic%t_delta )
308 
309  !--- Restart info
310  if( fstrdynamic%restart_nout > 0 ) then
311  if( mod(i,fstrdynamic%restart_nout).eq.0 .or. i.eq.fstrdynamic%n_step) then
312  call fstr_write_restart_dyna_nl(i,hecmesh,fstrsolid,fstrdynamic,fstrparam)
313  endif
314  endif
315 
316  enddo
317  !C-- end of time step loop
318  time_2 = hecmw_wtime()
319 
320  if( hecmesh%my_rank == 0 ) then
321  write(ista,'(a,f10.2,a)') ' solve (sec) :', time_2 - time_1, 's'
322  endif
323 
324  deallocate(coord)
325  call hecmw_mpc_mat_finalize(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
326  if (associated(hecmat0)) then
327  call hecmw_mat_finalize(hecmat0)
328  deallocate(hecmat0)
329  endif
330  end subroutine fstr_solve_dynamic_nlimplicit
331 
334  subroutine fstr_solve_dynamic_nlimplicit_contactslag(cstep, hecMESH,hecMAT,fstrSOLID,fstrEIG &
335  ,fstrDYNAMIC,fstrRESULT,fstrPARAM &
336  ,fstrCPL,hecLagMAT,restrt_step_num,infoCTChange &
337  ,conMAT )
338  implicit none
339  !C-- global variable
340  integer, intent(in) :: cstep
341  type(hecmwst_local_mesh) :: hecMESH
342  type(hecmwst_matrix) :: hecMAT
343  type(hecmwst_matrix), pointer :: hecMAT0
344  type(fstr_eigen) :: fstrEIG
345  type(fstr_solid) :: fstrSOLID
346  type(hecmwst_result_data) :: fstrRESULT
347  type(fstr_param) :: fstrPARAM
348  type(fstr_dynamic) :: fstrDYNAMIC
349  type(fstr_couple) :: fstrCPL !for COUPLE
350  type(hecmwst_matrix_lagrange) :: hecLagMAT
351  type(fstr_info_contactchange) :: infoCTChange
352  type(hecmwst_matrix) :: conMAT
353 
354  !C-- local variable
355  integer(kind=kint) :: nnod, ndof, numnp, nn
356  integer(kind=kint) :: i, j, ids, ide, ims, ime, kk, idm, imm
357  integer(kind=kint) :: iter
358  real(kind=kreal) :: a1, a2, a3, b1, b2, b3, c1, c2
359  real(kind=kreal) :: bsize, res, res1, rf
360  real(kind=kreal) :: res0, relres
361  real :: time_1, time_2
362  real(kind=kreal), parameter :: pi = 3.14159265358979323846d0
363 
364  integer(kind=kint) :: restrt_step_num
365  integer(kind=kint) :: ctAlgo
366  integer(kind=kint) :: max_iter_contact, count_step
367  integer(kind=kint) :: stepcnt
368  real(kind=kreal) :: maxdlag, converg_dlag
369 
370  integer(kind=kint) :: n_node_global
371  integer(kind=kint) :: contact_changed_global
372  integer(kind=kint) :: nndof,npdof
373  logical :: is_mat_symmetric
374  integer :: istat
375  logical :: is_cycle
376  real(kind=kreal),allocatable :: tmp_conb(:)
377  real(kind=kreal), allocatable :: coord(:)
378 
379  nullify(hecmat0)
380 
381  ! sum of n_node among all subdomains (to be used to calc res)
382  n_node_global = hecmesh%nn_internal
383  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
384 
385  ctalgo = fstrparam%contact_algo
386 
387  if( hecmat%Iarray(99)==4 .and. .not.fstr_is_matrixstruct_symmetric(fstrsolid,hecmesh) ) then
388  write(*,*) ' This type of direct solver is not yet available in such case ! '
389  write(*,*) ' Please use intel MKL direct solver !'
390  call hecmw_abort(hecmw_comm_get_comm())
391  endif
392 
393  hecmat%NDOF=hecmesh%n_dof
394 
395  nnod=hecmesh%n_node
396  ndof=hecmat%NDOF
397  nn=ndof*ndof
398 
399  allocate(coord(hecmesh%n_node*ndof))
400  if( associated( fstrsolid%contacts ) ) call initialize_contact_output_vectors(fstrsolid,hecmat)
401 
402  !!-- initial value
403  time_1 = hecmw_wtime()
404 
405  !C-- check parameters
406  if(dabs(fstrdynamic%beta) < 1.0e-20) then
407  if( hecmesh%my_rank == 0 ) then
408  write(imsg,*) 'stop due to Newmark-beta = 0'
409  endif
410  call hecmw_abort( hecmw_comm_get_comm())
411  endif
412 
413  !C-- matrix [M] lumped mass matrix
414  if(fstrdynamic%idx_mas == 1) then
415  call setmass(fstrsolid,hecmesh,hecmat,fstreig)
416 
417  !C-- consistent mass matrix
418  else if(fstrdynamic%idx_mas == 2) then
419  if( hecmesh%my_rank .eq. 0 ) then
420  write(imsg,*) 'stop: consistent mass matrix is not yet available !'
421  endif
422  call hecmw_abort( hecmw_comm_get_comm())
423  endif
424 
425  hecmat%Iarray(98) = 1 !Assembly complete
426  hecmat%Iarray(97) = 1 !Need numerical factorization
427 
428  !C-- initialize variables
429  if( restrt_step_num == 1 .and. fstrdynamic%VarInitialize .and. fstrdynamic%ray_m /= 0.0d0 ) &
430  call dynamic_init_varibles( hecmesh, hecmat, fstrsolid, fstreig, fstrdynamic, fstrparam )
431 
432  !C-- time step loop
433  a1 = .5d0/fstrdynamic%beta - 1.d0
434  a2 = 1.d0/(fstrdynamic%beta*fstrdynamic%t_delta)
435  a3 = 1.d0/(fstrdynamic%beta*fstrdynamic%t_delta*fstrdynamic%t_delta)
436  b1 = ( .5d0*fstrdynamic%gamma/fstrdynamic%beta - 1.d0 )*fstrdynamic%t_delta
437  b2 = fstrdynamic%gamma/fstrdynamic%beta - 1.d0
438  b3 = fstrdynamic%gamma/(fstrdynamic%beta*fstrdynamic%t_delta)
439  c1 = 1.d0 + fstrdynamic%ray_k*b3
440  c2 = a3 + fstrdynamic%ray_m*b3
441 
442  !C-- output of initial state
443  if( restrt_step_num == 1 ) then
444  call fstr_dynamic_output(hecmesh, fstrsolid, fstrdynamic, fstrparam)
445  call dynamic_output_monit(hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
446  endif
447 
448  fstrdynamic%VEC3(:) =0.d0
449  hecmat%X(:) =0.d0
450 
452  call fstr_scan_contact_state(cstep, restrt_step_num, 0, fstrdynamic%t_delta, ctalgo, hecmesh, fstrsolid, infoctchange)
453 
454  call hecmw_mat_copy_profile( hecmat, conmat )
455 
456  if ( fstr_is_contact_active() ) then
457  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
458  elseif( hecmat%Iarray(99)==4 ) then
459  write(*,*) ' This type of direct solver is not yet available in such case ! '
460  write(*,*) ' Please change solver type to intel MKL direct solver !'
461  call hecmw_abort(hecmw_comm_get_comm())
462  endif
463  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid,hecmesh)
464  call solve_lineq_contact_init(hecmesh,hecmat,heclagmat,is_mat_symmetric)
465 
466  max_iter_contact = fstrsolid%step_ctrl(cstep)%max_contiter
467  converg_dlag = fstrsolid%step_ctrl(cstep)%converg_lag
468 
469  !! step = 1,2,....,fstrDYNAMIC%n_step
470  do i = restrt_step_num, fstrdynamic%n_step
471  if (ndof == 4 .and. hecmesh%my_rank==0) write(*,'(a,i5)')"iter: ",i
472 
473  fstrdynamic%i_step = i
474  fstrdynamic%t_curr = fstrdynamic%t_delta * i
475 
476  if(hecmesh%my_rank==0) then
477  write(ista,'('' time step='',i10,'' time='',1pe13.4e3)') i,fstrdynamic%t_curr
478  write(*,'(A)')'-------------------------------------------------'
479  write(*,'('' time step='',i10,'' time='',1pe13.4e3)') i,fstrdynamic%t_curr
480  endif
481 
482  do j = 1 ,ndof*nnod
483  fstrdynamic%VEC1(j) = a1*fstrdynamic%ACC(j,1) + a2*fstrdynamic%VEL(j,1)
484  fstrdynamic%VEC2(j) = b1*fstrdynamic%ACC(j,1) + b2*fstrdynamic%VEL(j,1)
485  enddo
486 
487  count_step = 0
488  stepcnt = 0
489 
490  !C for couple analysis
491  do
492  fstrsolid%dunode(:) =0.d0
493  ! call fstr_UpdateEPState( hecMESH, fstrSOLID )
494  call fstr_solve_dynamic_nlimplicit_couple_init(fstrparam, fstrcpl)
495 
496  loopforcontactanalysis: do while( .true. )
497  count_step = count_step + 1
498 
499  ! ----- Inner Iteration
500  res0 = 0.d0
501  res1 = 0.d0
502  relres = 1.d0
503 
504  do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
505  stepcnt=stepcnt+1
506  if (fstrparam%nlgeom) then
507  call fstr_stiffmatrix( hecmesh, hecmat, fstrsolid, fstrdynamic%t_curr, fstrdynamic%t_delta )
508  else
509  if (.not. associated(hecmat0)) then
510  call fstr_stiffmatrix( hecmesh, hecmat, fstrsolid, fstrdynamic%t_curr, fstrdynamic%t_delta )
511  allocate(hecmat0)
512  call hecmw_mat_init(hecmat0)
513  call hecmw_mat_copy_profile(hecmat, hecmat0)
514  call hecmw_mat_copy_val(hecmat, hecmat0)
515  else
516  call hecmw_mat_copy_val(hecmat0, hecmat)
517  endif
518  endif
519 
520  if( fstrdynamic%ray_k/=0.d0 .or. fstrdynamic%ray_m/=0.d0 ) then
521  do j = 1 ,ndof*nnod
522  hecmat%X(j) = fstrdynamic%VEC2(j) - b3*fstrsolid%dunode(j)
523  enddo
524  endif
525  if( fstrdynamic%ray_k/=0.d0 ) then
526  if( hecmesh%n_dof == 3 ) then
527  call hecmw_matvec (hecmesh, hecmat, hecmat%X, fstrdynamic%VEC3)
528  else if( hecmesh%n_dof == 2 ) then
529  call hecmw_matvec (hecmesh, hecmat, hecmat%X, fstrdynamic%VEC3)
530  else if( hecmesh%n_dof == 6 ) then
531  call matvec(fstrdynamic%VEC3, hecmat%X, hecmat, ndof, hecmat%D, hecmat%AU, hecmat%AL)
532  endif
533  endif
534 
535  !C-- mechanical boundary condition
536  call dynamic_mat_ass_load (hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, iter)
537  do j=1, hecmesh%n_node* hecmesh%n_dof
538  hecmat%B(j)=hecmat%B(j)- fstrsolid%QFORCE(j) + fstreig%mass(j)*( fstrdynamic%VEC1(j)-a3*fstrsolid%dunode(j) &
539  + fstrdynamic%ray_m* hecmat%X(j) ) + fstrdynamic%ray_k*fstrdynamic%VEC3(j)
540  enddo
541 
542  !C for couple analysis
543  call fstr_solve_dynamic_nlimplicit_couple_pre(hecmesh, hecmat, fstrsolid, &
544  & fstrparam, fstrdynamic, fstrcpl, restrt_step_num, pi, i)
545 
546  do j = 1 ,nn*hecmat%NP
547  hecmat%D(j) = c1* hecmat%D(j)
548  enddo
549  do j = 1 ,nn*hecmat%NPU
550  hecmat%AU(j) = c1* hecmat%AU(j)
551  enddo
552  do j = 1 ,nn*hecmat%NPL
553  hecmat%AL(j) = c1*hecmat%AL(j)
554  enddo
555  do j=1,nnod
556  do kk=1,ndof
557  idm = nn*(j-1)+1 + (ndof+1)*(kk-1)
558  imm = ndof*(j-1) + kk
559  hecmat%D(idm) = hecmat%D(idm) + c2*fstreig%mass(imm)
560  enddo
561  enddo
562 
563  call hecmw_mat_clear( conmat )
564  call hecmw_mat_clear_b( conmat )
565  conmat%X = 0.0d0
566 
567  if( fstr_is_contact_active() ) then
568  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
569  call fstr_addcontactstiffness(cstep,ctalgo,iter,hecmesh,conmat,heclagmat,fstrsolid)
570  endif
571 
572  !C-- geometrical boundary condition
573  call dynamic_mat_ass_bc (hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, heclagmat, stepcnt, conmat=conmat)
574  call dynamic_mat_ass_bc_vl(hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, heclagmat, stepcnt, conmat=conmat)
575  call dynamic_mat_ass_bc_ac(hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, heclagmat, stepcnt, conmat=conmat)
576 
577  ! ----- check convergence
578  res = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
579 
580  if(iter == 1)then
581  res0 = res
582  endif
583 
584  ! ----- check convergence
585  if( .not.fstr_is_contact_active() ) then
586  maxdlag = 0.0d0
587  elseif( maxdlag == 0.0d0) then
588  maxdlag = 1.0d0
589  endif
590  call hecmw_allreduce_r1(hecmesh, maxdlag, hecmw_max)
591 
592  res = dsqrt(res/res0)
593  if( fstrparam%nlgeom .and. ndof /= 4 ) then
594  if( hecmesh%my_rank==0 ) then
595  write(*,'(a,i5,a,1pe12.4)')"iter: ",iter,", res: ",res
596  write(ista,'(''iter='',I5,''- Residual'',E15.7)')iter,res
597  write(*,'(a,1e15.7)') ' - MaxDLag =',maxdlag
598  write(ista,'(a,1e15.7)') ' - MaxDLag =',maxdlag
599  endif
600  if( res<fstrsolid%step_ctrl(cstep)%converg .and. maxdlag < converg_dlag ) exit
601  endif
602 
603  ! ---- For Parallel Contact with Multi-Partition Domains
604  hecmat%X = 0.0d0
605  call fstr_set_current_config_to_mesh(hecmesh,fstrsolid,coord)
606  call solve_lineq_contact(hecmesh,hecmat,heclagmat,conmat,istat,1.0d0,fstr_is_contact_active())
607  call fstr_recover_initial_config_to_mesh(hecmesh,fstrsolid,coord)
608 
609  ! ----- update external nodal displacement increments
610  call hecmw_update_r (hecmesh, hecmat%X, hecmat%NP, hecmat%NDOF)
611 
612  ! ----- update the strain, stress, and internal force
613  do j=1,hecmesh%n_node*ndof
614  fstrsolid%dunode(j) = fstrsolid%dunode(j)+hecmat%X(j)
615  enddo
616  call fstr_updatenewton( hecmesh, hecmat, fstrsolid, fstrdynamic%t_curr, &
617  & fstrdynamic%t_delta,iter, fstrdynamic%strainEnergy )
618 
619  if(.not. fstrparam%nlgeom) exit
620  if(ndof == 4) exit
621 
622  ! ----- update the Lagrange multipliers
623  if( fstr_is_contact_active() ) then
624  maxdlag = 0.0d0
625  do j=1,heclagmat%num_lagrange
626  heclagmat%lagrange(j) = heclagmat%lagrange(j) + hecmat%X(hecmesh%n_node*ndof+j)
627  if(dabs(hecmat%X(hecmesh%n_node*ndof+j))>maxdlag) maxdlag=dabs(hecmat%X(hecmesh%n_node*ndof+j))
628  ! write(*,*)'Lagrange:', j,hecLagMAT%lagrange(j),hecMAT%X(hecMESH%n_node*ndof+j)
629  enddo
630  endif
631  enddo
632 
633  ! ----- not convergence
634  if( iter>fstrsolid%step_ctrl(cstep)%max_iter ) then
635  if( hecmesh%my_rank==0) then
636  write(ilog,*) '### Fail to Converge : at step=', i
637  write(ista,*) '### Fail to Converge : at step=', i
638  write( *,*) ' ### Fail to Converge : at step=', i
639  endif
640  stop
641  endif
642 
643  ! ----- compute CONT_NFORCE/CONT_FRIC for output
644  if( fstr_is_contact_active() ) &
645  call fstr_calc_contact_output_force(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
646 
647  call fstr_scan_contact_state(cstep, i, count_step, fstrdynamic%t_delta, ctalgo, hecmesh, fstrsolid, infoctchange)
648 
649  if( hecmat%Iarray(99)==4 .and. .not. fstr_is_contact_active() ) then
650  write(*,*) ' This type of direct solver is not yet available in such case ! '
651  write(*,*) ' Please use intel MKL direct solver !'
652  call hecmw_abort(hecmw_comm_get_comm())
653  endif
654 
655  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid,hecmesh)
656  contact_changed_global=0
657  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) ) then
658  exit loopforcontactanalysis
659  elseif( fstr_is_matrixstructure_changed(infoctchange) ) then
660  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
661  contact_changed_global=1
662  endif
663  call hecmw_allreduce_i1(hecmesh,contact_changed_global,hecmw_max)
664  if (contact_changed_global > 0) then
665  call hecmw_mat_clear_b( hecmat )
666  call hecmw_mat_clear_b( conmat )
667  call solve_lineq_contact_init(hecmesh,hecmat,heclagmat,is_mat_symmetric)
668  endif
669 
670  if( count_step > max_iter_contact ) exit loopforcontactanalysis
671 
672  enddo loopforcontactanalysis
673 
674  !C for couple analysis
675  call fstr_solve_dynamic_nlimplicit_couple_post(hecmesh, hecmat, fstrsolid, &
676  & fstrparam, fstrdynamic, fstrcpl, a1, a2, a3, b1, b2, b3, i, is_cycle)
677  if(is_cycle) cycle
678  exit
679  enddo
680 
681  !C-- new displacement, velocity and acceleration
682  fstrdynamic%kineticEnergy = 0.0d0
683  do j = 1 ,ndof*nnod
684  fstrdynamic%ACC (j,2) = -a1*fstrdynamic%ACC(j,1) - a2*fstrdynamic%VEL(j,1) + &
685  a3*fstrsolid%dunode(j)
686  fstrdynamic%VEL (j,2) = -b1*fstrdynamic%ACC(j,1) - b2*fstrdynamic%VEL(j,1) + &
687  b3*fstrsolid%dunode(j)
688  fstrdynamic%ACC (j,1) = fstrdynamic%ACC (j,2)
689  fstrdynamic%VEL (j,1) = fstrdynamic%VEL (j,2)
690 
691  fstrsolid%unode(j) = fstrsolid%unode(j)+fstrsolid%dunode(j)
692  fstrdynamic%DISP(j,2) = fstrsolid%unode(j)
693 
694  fstrdynamic%kineticEnergy = fstrdynamic%kineticEnergy + &
695  0.5d0*fstreig%mass(j)*fstrdynamic%VEL(j,2)*fstrdynamic%VEL(j,2)
696  enddo
697 
698  !C-- output new displacement, velocity and acceleration
699  call fstr_dynamic_output(hecmesh, fstrsolid, fstrdynamic, fstrparam)
700 
701  !C-- output result of monitoring node
702  call dynamic_output_monit(hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
703 
704  call fstr_updatestate( hecmesh, fstrsolid, fstrdynamic%t_delta )
705 
706  !--- Restart info
707  if( fstrdynamic%restart_nout > 0 ) then
708  if( mod(i,fstrdynamic%restart_nout).eq.0 .or. i.eq.fstrdynamic%n_step ) then
709  call fstr_write_restart_dyna_nl(i,hecmesh,fstrsolid,fstrdynamic,fstrparam,&
710  infoctchange%contactNode_current)
711  endif
712  endif
713 
714  enddo
715  !C-- end of time step loop
716 
717  if (associated(hecmat0)) then
718  call hecmw_mat_finalize(hecmat0)
719  deallocate(hecmat0)
720  endif
721 
722  time_2 = hecmw_wtime()
723  if( hecmesh%my_rank == 0 ) then
724  write(ista,'(a,f10.2,a)') ' solve (sec) :', time_2 - time_1, 's'
725  endif
726 
727  deallocate(coord)
729 
730  subroutine fstr_solve_dynamic_nlimplicit_couple_init(fstrPARAM, fstrCPL)
731  implicit none
732  type(fstr_param) :: fstrparam
733  type(fstr_couple) :: fstrCPL
734  if( fstrparam%fg_couple == 1) then
735  if( fstrparam%fg_couple_type==1 .or. &
736  fstrparam%fg_couple_type==3 .or. &
737  fstrparam%fg_couple_type==5 ) call fstr_rcap_get( fstrcpl )
738  endif
740 
741  subroutine fstr_solve_dynamic_nlimplicit_couple_pre(hecMESH, hecMAT, fstrSOLID, &
742  & fstrPARAM, fstrDYNAMIC, fstrCPL, restrt_step_num, PI, i)
743  implicit none
744  type(hecmwst_local_mesh) :: hecMESH
745  type(hecmwst_matrix) :: hecMAT
746  type(fstr_solid) :: fstrSOLID
747  type(fstr_param) :: fstrPARAM
748  type(fstr_dynamic) :: fstrDYNAMIC
749  type(fstr_couple) :: fstrCPL
750  integer(kint) :: kkk0, kkk1, j, kk, i, restrt_step_num
751  real(kreal) :: bsize, PI
752 
753  if( fstrparam%fg_couple == 1) then
754  if( fstrparam%fg_couple_first /= 0 ) then
755  bsize = dfloat( i ) / dfloat( fstrparam%fg_couple_first )
756  if( bsize > 1.0 ) bsize = 1.0
757  do kkk0 = 1, fstrcpl%coupled_node_n
758  kkk1 = 3 * kkk0
759  fstrcpl%trac(kkk1-2) = bsize * fstrcpl%trac(kkk1-2)
760  fstrcpl%trac(kkk1-1) = bsize * fstrcpl%trac(kkk1-1)
761  fstrcpl%trac(kkk1 ) = bsize * fstrcpl%trac(kkk1 )
762  enddo
763  endif
764  if( fstrparam%fg_couple_window > 0 ) then
765  j = i - restrt_step_num + 1
766  kk = fstrdynamic%n_step - restrt_step_num + 1
767  bsize = 0.5*(1.0-cos(2.0*pi*dfloat(j)/dfloat(kk)))
768  do kkk0 = 1, fstrcpl%coupled_node_n
769  kkk1 = 3 * kkk0
770  fstrcpl%trac(kkk1-2) = bsize * fstrcpl%trac(kkk1-2)
771  fstrcpl%trac(kkk1-1) = bsize * fstrcpl%trac(kkk1-1)
772  fstrcpl%trac(kkk1 ) = bsize * fstrcpl%trac(kkk1 )
773  enddo
774  endif
775  call dynamic_mat_ass_couple( hecmesh, hecmat, fstrsolid, fstrcpl )
776  endif
778 
779  subroutine fstr_solve_dynamic_nlimplicit_couple_post(hecMESH, hecMAT, fstrSOLID, &
780  & fstrPARAM, fstrDYNAMIC, fstrCPL, a1, a2, a3, b1, b2, b3, i, is_cycle)
781  implicit none
782  type(hecmwst_local_mesh) :: hecMESH
783  type(hecmwst_matrix) :: hecMAT
784  type(fstr_solid) :: fstrSOLID
785  type(fstr_param) :: fstrPARAM
786  type(fstr_dynamic) :: fstrDYNAMIC
787  type(fstr_couple) :: fstrCPL
788  integer(kint) :: kkk0, kkk1, j, i, revocap_flag
789  real(kreal) :: bsize, a1, a2, a3, b1, b2, b3
790  logical :: is_cycle
791 
792  is_cycle = .false.
793 
794  if( fstrparam%fg_couple == 1 ) then
795  if( fstrparam%fg_couple_type>1 ) then
796  do j=1, fstrcpl%coupled_node_n
797  if( fstrcpl%dof == 3 ) then
798  kkk0 = j*3
799  kkk1 = fstrcpl%coupled_node(j)*3
800 
801  fstrcpl%disp (kkk0-2) = fstrsolid%unode(kkk1-2) + fstrsolid%dunode(kkk1-2)
802  fstrcpl%disp (kkk0-1) = fstrsolid%unode(kkk1-1) + fstrsolid%dunode(kkk1-1)
803  fstrcpl%disp (kkk0 ) = fstrsolid%unode(kkk1 ) + fstrsolid%dunode(kkk1 )
804 
805  fstrcpl%velo (kkk0-2) = -b1*fstrdynamic%ACC(kkk1-2,1) - b2*fstrdynamic%VEL(kkk1-2,1) + &
806  b3*fstrsolid%dunode(kkk1-2)
807  fstrcpl%velo (kkk0-1) = -b1*fstrdynamic%ACC(kkk1-1,1) - b2*fstrdynamic%VEL(kkk1-1,1) + &
808  b3*fstrsolid%dunode(kkk1-1)
809  fstrcpl%velo (kkk0 ) = -b1*fstrdynamic%ACC(kkk1,1) - b2*fstrdynamic%VEL(kkk1,1) + &
810  b3*fstrsolid%dunode(kkk1)
811  fstrcpl%accel(kkk0-2) = -a1*fstrdynamic%ACC(kkk1-2,1) - a2*fstrdynamic%VEL(kkk1-2,1) + &
812  a3*fstrsolid%dunode(kkk1-2)
813  fstrcpl%accel(kkk0-1) = -a1*fstrdynamic%ACC(kkk1-1,1) - a2*fstrdynamic%VEL(kkk1-1,1) + &
814  a3*fstrsolid%dunode(kkk1-1)
815  fstrcpl%accel(kkk0 ) = -a1*fstrdynamic%ACC(kkk1,1) - a2*fstrdynamic%VEL(kkk1,1) + &
816  a3*fstrsolid%dunode(kkk1)
817  else
818  kkk0 = j*2
819  kkk1 = fstrcpl%coupled_node(j)*2
820 
821  fstrcpl%disp (kkk0-1) = fstrsolid%unode(kkk1-1) + fstrsolid%dunode(kkk1-1)
822  fstrcpl%disp (kkk0 ) = fstrsolid%unode(kkk1 ) + fstrsolid%dunode(kkk1 )
823 
824  fstrcpl%velo (kkk0-1) = -b1*fstrdynamic%ACC(kkk1-1,1) - b2*fstrdynamic%VEL(kkk1-1,1) + &
825  b3*fstrsolid%dunode(kkk1-1)
826  fstrcpl%velo (kkk0 ) = -b1*fstrdynamic%ACC(kkk1,1) - b2*fstrdynamic%VEL(kkk1,1) + &
827  b3*fstrsolid%dunode(kkk1)
828  fstrcpl%accel(kkk0-1) = -a1*fstrdynamic%ACC(kkk1-1,1) - a2*fstrdynamic%VEL(kkk1-1,1) + &
829  a3*fstrsolid%dunode(kkk1-1)
830  fstrcpl%accel(kkk0 ) = -a1*fstrdynamic%ACC(kkk1,1) - a2*fstrdynamic%VEL(kkk1,1) + &
831  a3*fstrsolid%dunode(kkk1)
832  endif
833  enddo
834  call fstr_rcap_send( fstrcpl )
835  endif
836 
837  select case ( fstrparam%fg_couple_type )
838  case (4)
839  call fstr_rcap_get( fstrcpl )
840  case (5)
841  call fstr_get_convergence( revocap_flag )
842  if( revocap_flag==0 ) is_cycle = .true.
843  case (6)
844  call fstr_get_convergence( revocap_flag )
845  if( revocap_flag==0 ) then
846  call fstr_rcap_get( fstrcpl )
847  is_cycle = .true.
848  else
849  if( i /= fstrdynamic%n_step ) call fstr_rcap_get( fstrcpl )
850  endif
851  end select
852  endif
854 
855 end module fstr_dynamic_nlimplicit
This module contains subroutines for nonlinear implicit dynamic analysis.
subroutine fstr_solve_dynamic_nlimplicit_couple_init(fstrPARAM, fstrCPL)
subroutine fstr_solve_dynamic_nlimplicit_couple_post(hecMESH, hecMAT, fstrSOLID, fstrPARAM, fstrDYNAMIC, fstrCPL, a1, a2, a3, b1, b2, b3, i, is_cycle)
subroutine fstr_solve_dynamic_nlimplicit(cstep, hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrRESULT, fstrPARAM, fstrCPL, restrt_step_num)
This subroutine provides function of nonlinear implicit dynamic analysis using the Newmark method.
subroutine fstr_solve_dynamic_nlimplicit_couple_pre(hecMESH, hecMAT, fstrSOLID, fstrPARAM, fstrDYNAMIC, fstrCPL, restrt_step_num, PI, i)
subroutine fstr_solve_dynamic_nlimplicit_contactslag(cstep, hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrRESULT, fstrPARAM, fstrCPL, hecLagMAT, restrt_step_num, infoCTChange, conMAT)
This subroutine provides function of nonlinear implicit dynamic analysis using the Newmark method....
This module provides functions of reconstructing.
subroutine, public fstr_mat_con_contact(cstep, contact_algo, hecMAT, fstrSOLID, hecLagMAT, infoCTChange, conMAT, is_contact_active_flag)
this subroutine reconstructs node-based (stiffness) matrix structure \corresponding to contact state
subroutine, public fstr_save_originalmatrixstructure(hecMAT)
This subroutine saves original matrix structure constructed originally by hecMW_matrix.
logical function, public fstr_is_matrixstruct_symmetric(fstrSOLID, hecMESH)
this function judges whether sitiffness matrix is symmetric or not
This module provides functions to initialize variables when initial velocity or acceleration boundary...
subroutine dynamic_init_varibles(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrPARAM)
This module contains functions to set acceleration boundary condition in dynamic analysis.
subroutine dynamic_mat_ass_bc_ac(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, iter, conMAT)
This subrouitne set acceleration boundary condition in dynamic analysis.
This module contains functions to set velocity boundary condition in dynamic analysis.
subroutine dynamic_mat_ass_bc_vl(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, iter, conMAT)
This subrouitne set velocity boundary condition in dynamic analysis.
This module contains functions to set displacement boundary condition in dynamic analysis.
subroutine dynamic_mat_ass_bc(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, iter, conMAT)
This subroutine setup disp bundary condition.
This module contains functions relates to coupling analysis.
subroutine dynamic_mat_ass_couple(hecMESH, hecMAT, fstrSOLID, fstrCPL)
This module contains function to set boundary condition of external load in dynamic analysis.
subroutine dynamic_mat_ass_load(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, iter)
This function sets boundary condition of external load.
This module provides functions to output result.
subroutine dynamic_output_monit(hecMESH, fstrPARAM, fstrDYNAMIC, fstrEIG, fstrSOLID)
subroutine matvec(y, x, hecMAT, ndof, D, AU, AL)
subroutine fstr_dynamic_output(hecMESH, fstrSOLID, fstrDYNAMIC, fstrPARAM)
Output result.
Set up lumped mass matrix.
subroutine setmass(fstrSOLID, hecMESH, hecMAT, fstrEIG)
subroutine, public fstr_rcap_send(fstrCPL)
subroutine, public fstr_rcap_get(fstrCPL)
subroutine fstr_get_convergence(revocap_flag)
This module provides function to calculate residual of nodal force.
real(kind=kreal) function, public fstr_get_norm_para_contact(hecMAT, hecLagMAT, conMAT, hecMESH)
This module provides functions to read in and write out restart files.
Definition: fstr_Restart.f90:8
subroutine fstr_write_restart_dyna_nl(cstep, hecMESH, fstrSOLID, fstrDYNAMIC, fstrPARAM, contactNode)
write out restart file for nonlinear dynamic analysis
This module provides function to calculate tangent stiffness matrix.
subroutine, public fstr_stiffmatrix(hecMESH, hecMAT, fstrSOLID, time, tincr)
This subroutine creates tangential stiffness matrix.
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
subroutine fstr_recover_initial_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.F90:1141
integer(kind=kint), parameter imsg
Definition: m_fstr.F90:112
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.F90:109
subroutine fstr_set_current_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.F90:1128
integer(kind=kint), parameter ista
Definition: m_fstr.F90:110
This module provides functions to solve sparse system of \linear equitions in the case of contact ana...
subroutine, public solve_lineq_contact_init(hecMESH, hecMAT, hecLagMAT, is_sym)
This subroutine.
subroutine, public solve_lineq_contact(hecMESH, hecMAT, hecLagMAT, conMAT, istat, rf, is_contact_active)
This subroutine.
Top-level contact analysis module (System level)
subroutine fstr_update_ndforce_contact(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for residual vector (conMATB).
logical function fstr_is_matrixstructure_changed(infoCTChange)
logical function fstr_is_contact_conv(ctAlgo, infoCTChange, hecMESH)
subroutine fstr_addcontactstiffness(cstep, ctAlgo, iter, hecMESH, conMAT, hecLagMAT, fstrSOLID)
subroutine fstr_calc_contact_output_force(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for output (CONT_NFORCE/CONT_FRIC).
subroutine fstr_scan_contact_state(cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange)
Scanning contact state.
logical function fstr_is_contact_active()
Data for coupling analysis.
Definition: m_fstr.F90:622
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Definition: m_fstr.F90:515
Package of data used by Lanczos eigenvalue solver.
Definition: m_fstr.F90:604
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.F90:156