14 subroutine dynamic_mat_ass_bc_ac(cstep, hecMESH, hecMAT, fstrSOLID ,fstrDYNAMIC, fstrPARAM, hecLagMAT, t_curr, iter, conMAT)
20 integer(kind=kint) :: cstep
21 type(hecmwst_matrix) :: hecMAT
22 type(hecmwst_local_mesh) :: hecMESH
26 type(hecmwst_matrix_lagrange) :: hecLagMAT
27 real(kind=kreal) :: t_curr
28 integer,
optional :: iter
29 type(hecmwst_matrix),
optional :: conMAT
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
36 if( fstrsolid%ACCELERATION_type ==
kbcinitial )
return
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
50 if( fstrdynamic%idx_eqa == 1 )
then
52 do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
53 ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
54 grpid = fstrsolid%ACCELERATION_ngrp_GRPID(ig0)
56 rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
58 call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
62 ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
65 idofe = ityp - idofs*10
67 is0 = hecmesh%node_group%grp_index(ig-1) + 1
68 ie0 = hecmesh%node_group%grp_index(ig )
71 in = hecmesh%node_group%grp_item(ik)
72 do idof = idofs, idofe
74 if(
present(iter) )
then
79 + b2*fstrdynamic%VEL (ndof*in-(ndof-idof),1) &
80 + b3*fstrdynamic%ACC (ndof*in-(ndof-idof),1) &
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) &
89 if(
present(conmat))
then
90 call hecmw_mat_ass_bc(hecmat, in, idof, rhs, conmat)
92 call hecmw_mat_ass_bc(hecmat, in, idof, rhs)
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)
99 call hecmw_mat_ass_bc_contactlag(hecmat,heclagmat,in,idof,rhs)
113 else if( fstrdynamic%idx_eqa == 11 )
then
115 do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
116 ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
117 rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
119 call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
123 ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
125 is0 = hecmesh%node_group%grp_index(ig-1) + 1
126 ie0 = hecmesh%node_group%grp_index(ig )
128 idofe = ityp - idofs*10
131 in = hecmesh%node_group%grp_item(ik)
133 do idof = idofs, idofe
134 rhs = 2.0*fstrdynamic%DISP(ndof*in-(ndof-idof),1) &
135 - fstrdynamic%DISP(ndof*in-(ndof-idof),3) &
137 hecmat%B (ndof*in-(ndof-idof)) = rhs
138 fstrdynamic%VEC1(ndof*in-(ndof-idof)) = 1.0d0
159 type(hecmwst_matrix) :: hecmat
160 type(hecmwst_local_mesh) :: hecMESH
163 real(kind=kreal) :: t_curr
165 integer(kind=kint) :: NDOF, ig0, ig, ityp, iS0, iE0, ik, in, idofS, idofE, idof
167 integer(kind=kint) :: flag_u, grpid
168 real(kind=kreal) :: rhs, f_t
170 if( fstrsolid%ACCELERATION_type ==
kbctransit )
return
177 do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
178 ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
179 grpid = fstrsolid%ACCELERATION_ngrp_GRPID(ig0)
181 rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
184 call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
187 ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
189 is0 = hecmesh%node_group%grp_index(ig-1) + 1
190 ie0 = hecmesh%node_group%grp_index(ig )
192 idofe = ityp - idofs*10
195 in = hecmesh%node_group%grp_item(ik)
197 do idof = idofs, idofe
198 fstrdynamic%ACC (ndof*in-(ndof-idof),1) = rhs
211 type(hecmwst_matrix) :: hecmat
212 type(hecmwst_local_mesh) :: hecMESH
215 real(kind=kreal) :: t_curr
216 integer,
optional :: iter
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
223 if( fstrsolid%ACCELERATION_type ==
kbcinitial )
return
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
234 do ig0 = 1, fstrsolid%ACCELERATION_ngrp_tot
235 ig = fstrsolid%ACCELERATION_ngrp_ID(ig0)
236 rhs = fstrsolid%ACCELERATION_ngrp_val(ig0)
238 call table_dyn(hecmesh, fstrsolid, fstrdynamic, ig0, t_curr, f_t, flag_u)
242 ityp = fstrsolid%ACCELERATION_ngrp_type(ig0)
244 is0 = hecmesh%node_group%grp_index(ig-1) + 1
245 ie0 = hecmesh%node_group%grp_index(ig )
247 idofe = ityp - idofs*10
250 in = hecmesh%node_group%grp_item(ik)
252 do idof = idofs, idofe
253 rhs = 2.0*fstrdynamic%DISP(ndof*in-(ndof-idof),1) &
254 - fstrdynamic%DISP(ndof*in-(ndof-idof),3) &
256 hecmat%B(ndof*in-(ndof-idof)) = rhs* fstrdynamic%VEC1(ndof*in-(ndof-idof))
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.
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
integer(kind=kint), parameter kbcinitial
integer(kind=kint), parameter kbctransit
logical function fstr_isboundaryactive(fstrSOLID, nbc, cstep)
Table of lading step in dynamic analysis.
subroutine table_dyn(hecMESH, fstrSOLID, fstrDYNAMIC, ig0, t_curr, f_t, flag_u)
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
FSTR INNER CONTROL PARAMETERS (fstrPARAM)