FrontISTR  5.7.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(hecMESH, hecMAT, fstrSOLID ,fstrDYNAMIC, fstrPARAM, hecLagMAT, iter, conMAT)
15  use m_fstr
16  use m_table_dyn
17  use mcontact
18 
19  implicit none
20  type(hecmwst_matrix) :: hecmat
21  type(hecmwst_local_mesh) :: hecMESH
22  type(fstr_solid) :: fstrSOLID
23  type(fstr_dynamic) :: fstrDYNAMIC
24  type(fstr_param) :: fstrPARAM
25  type(hecmwst_matrix_lagrange) :: hecLagMAT
26  integer, optional :: iter
27  type(hecmwst_matrix), optional :: conMAT
28 
29  integer(kind=kint) :: ig0, ig, ityp, NDOF, iS0, iE0, ik, in, idofS, idofE, idof
30  integer(kind=kint) :: dyn_step, flag_u
31  real(kind=kreal) :: b2, b3, b4, c1
32  real(kind=kreal) :: rhs, rhs0, f_t
33 
34  if( fstrsolid%ACCELERATION_type == kbcinitial )return
35 
36  dyn_step = fstrdynamic%i_step
37  flag_u = 3
38 
39  b2 = fstrdynamic%t_delta
40  b3 = fstrdynamic%t_delta**2*(0.5d0-fstrdynamic%beta)
41  b4 = fstrdynamic%t_delta**2*fstrdynamic%beta
42  c1 = fstrdynamic%t_delta**2
43 
44  ndof = hecmat%NDOF
45 
46  !C=============================C
47  !C-- implicit dynamic analysis
48  !C=============================C
49  if( fstrdynamic%idx_eqa == 1 ) then
50 
51  do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
52  ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
53  rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
54 
55  call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, f_t, flag_u)
56  rhs = rhs * f_t
57  rhs0 = rhs
58 
59  ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
60 
61  idofs = ityp/10
62  idofe = ityp - idofs*10
63 
64  is0 = hecmesh%node_group%grp_index(ig-1) + 1
65  ie0 = hecmesh%node_group%grp_index(ig )
66 
67  do ik = is0, ie0
68  in = hecmesh%node_group%grp_item(ik)
69  do idof = idofs, idofe
70 
71  if( present(iter) ) then ! increment
72  if( iter>1 ) then
73  rhs = 0.d0
74  else
75  rhs = &
76  + b2*fstrdynamic%VEL (ndof*in-(ndof-idof),1) &
77  + b3*fstrdynamic%ACC (ndof*in-(ndof-idof),1) &
78  + b4*rhs0
79  endif
80  else
81  rhs = fstrdynamic%DISP(ndof*in-(ndof-idof),1) &
82  + b2*fstrdynamic%VEL (ndof*in-(ndof-idof),1) &
83  + b3*fstrdynamic%ACC (ndof*in-(ndof-idof),1) &
84  + b4*rhs0
85  endif
86  if(present(conmat)) then
87  call hecmw_mat_ass_bc(hecmat, in, idof, rhs, conmat)
88  else
89  call hecmw_mat_ass_bc(hecmat, in, idof, rhs)
90  endif
91  if( fstr_is_contact_active() .and. fstrparam%contact_algo == kcaslagrange &
92  .and. fstrparam%nlgeom .and. fstrdynamic%idx_resp == 1 ) then
93  if(present(conmat)) then
94  call hecmw_mat_ass_bc_contactlag(conmat,heclagmat,in,idof,rhs)
95  else
96  call hecmw_mat_ass_bc_contactlag(hecmat,heclagmat,in,idof,rhs)
97  endif
98  endif
99 
100  !for output reaction force
101  fstrsolid%REACTION(ndof*(in-1)+idof) = fstrsolid%QFORCE(ndof*(in-1)+idof)
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, 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 
140  !for output reaction force
141  fstrsolid%REACTION(ndof*(in-1)+idof) = fstrsolid%QFORCE(ndof*(in-1)+idof)
142  end do
143  enddo
144  enddo
145  !C
146  !C-- end of explicit dynamic analysis
147  !C
148  end if
149  !
150  return
151  end subroutine dynamic_mat_ass_bc_ac
152 
153 
154  !C***
156  !C***
157  subroutine dynamic_bc_init_ac(hecMESH, hecMAT, fstrSOLID ,fstrDYNAMIC)
158  use m_fstr
159  use m_table_dyn
160 
161  implicit none
162  type(hecmwst_matrix) :: hecmat
163  type(hecmwst_local_mesh) :: hecMESH
164  type(fstr_solid) :: fstrSOLID
165  type(fstr_dynamic) :: fstrDYNAMIC
166 
167  integer(kind=kint) :: NDOF, ig0, ig, ityp, iS0, iE0, ik, in, idofS, idofE, idof
168  !!!
169  integer(kind=kint) :: flag_u
170  real(kind=kreal) :: rhs, f_t
171 
172  if( fstrsolid%ACCELERATION_type == kbctransit )return
173 
174  !!!
175  flag_u = 3
176  !!!
177 
178  ndof = hecmat%NDOF
179  do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
180  ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
181  rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
182 
183  !!!!!! time history
184  call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, 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, 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  integer, optional :: iter
216 
217  integer(kind=kint) :: ig0, ig, ityp, NDOF, iS0, iE0, ik, in, idofS, idofE, idof
218  integer(kind=kint) :: dyn_step, flag_u
219  real(kind=kreal) :: b2, b3, b4, c1
220  real(kind=kreal) :: rhs, rhs0, f_t
221 
222  if( fstrsolid%ACCELERATION_type == kbcinitial )return
223 
224  dyn_step = fstrdynamic%i_step
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, 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 
259  !for output reaction force
260  fstrsolid%REACTION(ndof*(in-1)+idof) = fstrsolid%QFORCE(ndof*(in-1)+idof)
261  end do
262  enddo
263  enddo
264  end subroutine dynamic_explicit_ass_ac
265 
266 end module m_dynamic_mat_ass_bc_ac
m_table_dyn
Table of lading step in dynamic analysis.
Definition: table_dyn.f90:6
m_fstr::kbcinitial
integer(kind=kint), parameter kbcinitial
Definition: m_fstr.f90:66
m_fstr::kbctransit
integer(kind=kint), parameter kbctransit
Definition: m_fstr.f90:67
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::fstr_param
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:154
m_fstr::kcaslagrange
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.f90:59
m_dynamic_mat_ass_bc_ac::dynamic_bc_init_ac
subroutine dynamic_bc_init_ac(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC)
This function sets initial condition of acceleration.
Definition: dynamic_mat_ass_bc_ac.f90:158
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_dynamic_mat_ass_bc_ac::dynamic_mat_ass_bc_ac
subroutine dynamic_mat_ass_bc_ac(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, iter, conMAT)
This subrouitne set acceleration boundary condition in dynamic analysis.
Definition: dynamic_mat_ass_bc_ac.f90:15
mcontact::fstr_is_contact_active
logical function fstr_is_contact_active()
Definition: fstr_contact.f90:52
m_dynamic_mat_ass_bc_ac
This module contains functions to set acceleration boundary condition in dynamic analysis.
Definition: dynamic_mat_ass_bc_ac.f90:7
m_dynamic_mat_ass_bc_ac::dynamic_explicit_ass_ac
subroutine dynamic_explicit_ass_ac(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, iter)
Definition: dynamic_mat_ass_bc_ac.f90:208
m_table_dyn::table_dyn
subroutine table_dyn(hecMESH, fstrSOLID, fstrDYNAMIC, ig0, f_t, flag_u)
Definition: table_dyn.f90:19
mcontact
This module provides functions to calculate contact stiff matrix.
Definition: fstr_contact.f90:6