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