FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
dynamic_mat_ass_bc_ac.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 contains
9 
10  !C***
12  !C***
13 
14  subroutine dynamic_mat_ass_bc_ac(cstep, hecMESH, hecMAT, fstrSOLID ,fstrDYNAMIC, fstrPARAM, hecLagMAT, t_curr, iter, conMAT)
15  use m_fstr
16  use m_table_dyn
17  use mcontact
18 
19  implicit none
20  integer(kind=kint) :: cstep
21  type(hecmwst_matrix) :: hecMAT
22  type(hecmwst_local_mesh) :: hecMESH
23  type(fstr_solid) :: fstrSOLID
24  type(fstr_dynamic) :: fstrDYNAMIC
25  type(fstr_param) :: fstrPARAM
26  type(hecmwst_matrix_lagrange) :: hecLagMAT
27  real(kind=kreal) :: t_curr
28  integer, optional :: iter
29  type(hecmwst_matrix), optional :: conMAT
30 
31  integer(kind=kint) :: ig0, ig, ityp, NDOF, iS0, iE0, ik, in, idofS, idofE, idof
32  integer(kind=kint) :: flag_u, grpid
33  real(kind=kreal) :: b2, b3, b4, c1
34  real(kind=kreal) :: rhs, rhs0, f_t
35 
36  if( fstrsolid%ACCELERATION_type == kbcinitial )return
37 
38  flag_u = 3
39 
40  b2 = fstrdynamic%t_delta
41  b3 = fstrdynamic%t_delta**2*(0.5d0-fstrdynamic%beta)
42  b4 = fstrdynamic%t_delta**2*fstrdynamic%beta
43  c1 = fstrdynamic%t_delta**2
44 
45  ndof = hecmat%NDOF
46 
47  !C=============================C
48  !C-- implicit dynamic analysis
49  !C=============================C
50  if( fstrdynamic%idx_eqa == 1 ) then
51 
52  do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
53  ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
54  grpid = fstrsolid%ACCELERATION_ngrp_GRPID(ig0)
55  if( .not. fstr_isboundaryactive( fstrsolid, grpid, cstep ) ) cycle
56  rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
57 
58  call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
59  rhs = rhs * f_t
60  rhs0 = rhs
61 
62  ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
63 
64  idofs = ityp/10
65  idofe = ityp - idofs*10
66 
67  is0 = hecmesh%node_group%grp_index(ig-1) + 1
68  ie0 = hecmesh%node_group%grp_index(ig )
69 
70  do ik = is0, ie0
71  in = hecmesh%node_group%grp_item(ik)
72  do idof = idofs, idofe
73 
74  if( present(iter) ) then ! increment
75  if( iter>1 ) then
76  rhs = 0.d0
77  else
78  rhs = &
79  + b2*fstrdynamic%VEL (ndof*in-(ndof-idof),1) &
80  + b3*fstrdynamic%ACC (ndof*in-(ndof-idof),1) &
81  + b4*rhs0
82  endif
83  else
84  rhs = fstrdynamic%DISP(ndof*in-(ndof-idof),1) &
85  + b2*fstrdynamic%VEL (ndof*in-(ndof-idof),1) &
86  + b3*fstrdynamic%ACC (ndof*in-(ndof-idof),1) &
87  + b4*rhs0
88  endif
89  if(present(conmat)) then
90  call hecmw_mat_ass_bc(hecmat, in, idof, rhs, conmat)
91  else
92  call hecmw_mat_ass_bc(hecmat, in, idof, rhs)
93  endif
94  if( fstr_is_contact_active() .and. fstrparam%contact_algo == kcaslagrange &
95  .and. fstrparam%nlgeom .and. fstrdynamic%idx_resp == 1 ) then
96  if(present(conmat)) then
97  call hecmw_mat_ass_bc_contactlag(conmat,heclagmat,in,idof,rhs)
98  else
99  call hecmw_mat_ass_bc_contactlag(hecmat,heclagmat,in,idof,rhs)
100  endif
101  endif
102  enddo
103  enddo
104 
105  enddo
106  !C
107  !C-- end of implicit dynamic analysis
108  !C
109 
110  !C=============================C
111  !C-- explicit dynamic analysis
112  !C=============================C
113  else if( fstrdynamic%idx_eqa == 11 ) then
114  !C
115  do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
116  ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
117  rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
118 
119  call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
120  rhs = rhs * f_t
121  rhs0 = rhs
122 
123  ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
124 
125  is0 = hecmesh%node_group%grp_index(ig-1) + 1
126  ie0 = hecmesh%node_group%grp_index(ig )
127  idofs = ityp/10
128  idofe = ityp - idofs*10
129 
130  do ik = is0, ie0
131  in = hecmesh%node_group%grp_item(ik)
132 
133  do idof = idofs, idofe
134  rhs = 2.0*fstrdynamic%DISP(ndof*in-(ndof-idof),1) &
135  - fstrdynamic%DISP(ndof*in-(ndof-idof),3) &
136  + c1*rhs0
137  hecmat%B (ndof*in-(ndof-idof)) = rhs
138  fstrdynamic%VEC1(ndof*in-(ndof-idof)) = 1.0d0
139  end do
140  enddo
141  enddo
142  !C
143  !C-- end of explicit dynamic analysis
144  !C
145  end if
146  !
147  return
148  end subroutine dynamic_mat_ass_bc_ac
149 
150 
151  !C***
153  !C***
154  subroutine dynamic_bc_init_ac(hecMESH, hecMAT, fstrSOLID ,fstrDYNAMIC, t_curr)
155  use m_fstr
156  use m_table_dyn
157 
158  implicit none
159  type(hecmwst_matrix) :: hecmat
160  type(hecmwst_local_mesh) :: hecMESH
161  type(fstr_solid) :: fstrSOLID
162  type(fstr_dynamic) :: fstrDYNAMIC
163  real(kind=kreal) :: t_curr
164 
165  integer(kind=kint) :: NDOF, ig0, ig, ityp, iS0, iE0, ik, in, idofS, idofE, idof
166  !!!
167  integer(kind=kint) :: flag_u, grpid
168  real(kind=kreal) :: rhs, f_t
169 
170  if( fstrsolid%ACCELERATION_type == kbctransit )return
171 
172  !!!
173  flag_u = 3
174  !!!
175 
176  ndof = hecmat%NDOF
177  do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
178  ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
179  grpid = fstrsolid%ACCELERATION_ngrp_GRPID(ig0)
180  if( .not. fstr_isboundaryactive( fstrsolid, grpid, 1 ) ) cycle
181  rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
182 
183  !!!!!! time history
184  call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
185  rhs = rhs * f_t
186  !!!!!!
187  ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
188 
189  is0 = hecmesh%node_group%grp_index(ig-1) + 1
190  ie0 = hecmesh%node_group%grp_index(ig )
191  idofs = ityp/10
192  idofe = ityp - idofs*10
193 
194  do ik = is0, ie0
195  in = hecmesh%node_group%grp_item(ik)
196  !C
197  do idof = idofs, idofe
198  fstrdynamic%ACC (ndof*in-(ndof-idof),1) = rhs
199  end do
200  enddo
201  !*End ig0 (main) loop
202  enddo
203 
204  return
205  end subroutine dynamic_bc_init_ac
206 
207  subroutine dynamic_explicit_ass_ac(hecMESH, hecMAT, fstrSOLID ,fstrDYNAMIC, t_curr, iter)
208  use m_fstr
209  use m_table_dyn
210  use mcontact
211  type(hecmwst_matrix) :: hecmat
212  type(hecmwst_local_mesh) :: hecMESH
213  type(fstr_solid) :: fstrSOLID
214  type(fstr_dynamic) :: fstrDYNAMIC
215  real(kind=kreal) :: t_curr
216  integer, optional :: iter
217 
218  integer(kind=kint) :: ig0, ig, ityp, NDOF, iS0, iE0, ik, in, idofS, idofE, idof
219  integer(kind=kint) :: flag_u
220  real(kind=kreal) :: b2, b3, b4, c1
221  real(kind=kreal) :: rhs, rhs0, f_t
222 
223  if( fstrsolid%ACCELERATION_type == kbcinitial )return
224 
225  flag_u = 3
226 
227  b2 = fstrdynamic%t_delta
228  b3 = fstrdynamic%t_delta**2*(0.5d0-fstrdynamic%beta)
229  b4 = fstrdynamic%t_delta**2*fstrdynamic%beta
230  c1 = fstrdynamic%t_delta**2
231 
232  ndof = hecmat%NDOF
233 
234  do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
235  ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
236  rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
237 
238  call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
239  rhs = rhs * f_t
240  rhs0 = rhs
241 
242  ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
243 
244  is0 = hecmesh%node_group%grp_index(ig-1) + 1
245  ie0 = hecmesh%node_group%grp_index(ig )
246  idofs = ityp/10
247  idofe = ityp - idofs*10
248 
249  do ik = is0, ie0
250  in = hecmesh%node_group%grp_item(ik)
251 
252  do idof = idofs, idofe
253  rhs = 2.0*fstrdynamic%DISP(ndof*in-(ndof-idof),1) &
254  - fstrdynamic%DISP(ndof*in-(ndof-idof),3) &
255  + c1*rhs0
256  hecmat%B(ndof*in-(ndof-idof)) = rhs* fstrdynamic%VEC1(ndof*in-(ndof-idof))
257  ! fstrDYNAMIC%VEC1(NDOF*in-(NDOF-idof)) = 1.0d0
258  end do
259  enddo
260  enddo
261  end subroutine dynamic_explicit_ass_ac
262 
263 end module m_dynamic_mat_ass_bc_ac
This module contains functions to set acceleration boundary condition in dynamic analysis.
subroutine dynamic_bc_init_ac(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, t_curr)
This function sets initial condition of acceleration.
subroutine dynamic_explicit_ass_ac(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, t_curr, iter)
subroutine dynamic_mat_ass_bc_ac(cstep, hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, t_curr, iter, conMAT)
This subrouitne set acceleration boundary condition in dynamic analysis.
This module defines common data and basic structures for analysis.
Definition: m_fstr.F90:15
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.F90:61
integer(kind=kint), parameter kbcinitial
Definition: m_fstr.F90:68
integer(kind=kint), parameter kbctransit
Definition: m_fstr.F90:69
logical function fstr_isboundaryactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.F90:1053
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
Top-level contact analysis module (System level)
logical function fstr_is_contact_active()
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Definition: m_fstr.F90:530
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.F90:156