FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
table_dyn.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 module m_table_dyn
7 
8  use m_fstr
9  use hecmw
10 
11  implicit none
12 
13 contains
14 
15  !C================================================================C
16  !C-- subroutine table_dyn
17  !C================================================================C
18  subroutine table_dyn(hecMESH, fstrSOLID, fstrDYNAMIC, ig0, f_t, flag_u)
19  type(hecmwst_local_mesh) :: hecMESH
20  type(fstr_solid) :: fstrSOLID
21  type(fstr_dynamic) :: fstrDYNAMIC
22 
23  integer(kind=kint) :: i, ig0
24  integer(kind=kint) :: jj_n_amp, jj1, jj2
25  integer(kind=kint) :: s1, s2, flag_u
26  real(kind=kreal) :: t_1, t_2, t_t, f_1, f_2, f_t
27 
28  jj_n_amp = 0
29  s1 = 0; s2 = 0
30  t_1 = 0.0d0; t_2 = 0.0d0; t_t = 0.0d0; f_1 = 0.0d0; f_2 = 0.0d0; f_t = 0.0d0
31 
32  if( flag_u .eq. 1 ) then
33  jj_n_amp = fstrsolid%BOUNDARY_ngrp_amp(ig0)
34  else if( flag_u .eq. 2 ) then
35  jj_n_amp = fstrsolid%VELOCITY_ngrp_amp(ig0)
36  else if( flag_u .eq. 3 ) then
37  jj_n_amp = fstrsolid%ACCELERATION_ngrp_amp(ig0)
38  else if( flag_u .eq. 0 ) then
39  jj_n_amp = fstrsolid%CLOAD_ngrp_amp(ig0)
40  else if( flag_u .eq. 10 ) then
41  jj_n_amp = fstrsolid%DLOAD_ngrp_amp(ig0)
42  end if
43 
44  if( jj_n_amp == 0 ) then
45  f_t = 1.d0
46  else
47 
48  jj1 = hecmesh%amp%amp_index(jj_n_amp - 1)
49  jj2 = hecmesh%amp%amp_index(jj_n_amp)
50 
51  jj1 = jj1 + 2
52  if( fstrdynamic%idx_eqa == 1 ) then
53  t_t = fstrdynamic%t_curr
54 
55  else if( fstrdynamic%idx_eqa == 11 ) then
56  select case (flag_u)
57  case (0)
58  t_t = fstrdynamic%t_curr - fstrdynamic%t_delta
59  case (10)
60  t_t = fstrdynamic%t_curr - fstrdynamic%t_delta
61  case (1)
62  t_t = fstrdynamic%t_curr - fstrdynamic%t_delta
63  case (2)
64  t_t = fstrdynamic%t_curr - fstrdynamic%t_delta
65  case (3)
66  t_t = fstrdynamic%t_curr - fstrdynamic%t_delta
67  end select
68  end if
69 
70  if( fstrdynamic%i_step == 0 ) then
71  t_t = fstrdynamic%t_delta*fstrdynamic%i_step
72  end if
73 
74  ! if(jj2 .eq. 0) then
75  ! f_t = 1.0
76 
77  if(t_t .gt. hecmesh%amp%amp_table(jj2)) then
78  f_t = hecmesh%amp%amp_val(jj2)
79  else if(t_t .le. hecmesh%amp%amp_table(jj2)) then
80  do i = jj1, jj2
81  if(t_t .le. hecmesh%amp%amp_table(i)) then
82  s2 = i
83  s1 = i - 1
84  exit
85  end if
86  end do
87 
88  t_2 = hecmesh%amp%amp_table(s2)
89  t_1 = hecmesh%amp%amp_table(s1)
90  f_2 = hecmesh%amp%amp_val(s2)
91  f_1 = hecmesh%amp%amp_val(s1)
92  if( t_2-t_1 .lt. 1.0e-20) then
93  if( hecmesh%my_rank.eq.0) then
94  write(imsg,*) 'stop due to t_2-t_1 <= 0'
95  end if
96  call hecmw_abort( hecmw_comm_get_comm())
97  end if
98  f_t = ((t_2*f_1 - t_1*f_2) + (f_2 - f_1)*t_t) / (t_2 - t_1)
99  end if
100 
101  end if
102 
103  end subroutine table_dyn
104 
105 end module m_table_dyn
m_table_dyn
Table of lading step in dynamic analysis.
Definition: table_dyn.f90:6
m_fstr::fstr_solid
Definition: m_fstr.f90:238
m_fstr::fstr_dynamic
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Definition: m_fstr.f90:504
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
hecmw
Definition: hecmw.f90:6
m_table_dyn::table_dyn
subroutine table_dyn(hecMESH, fstrSOLID, fstrDYNAMIC, ig0, f_t, flag_u)
Definition: table_dyn.f90:19
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110