FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_Spring.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  implicit none
9 contains
10 
11  subroutine fstr_addspring(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
12  use m_fstr
13  use m_static_lib
14  integer, intent(in) :: cstep
15  type (hecmwST_matrix),intent(inout) :: hecMAT
16  type (hecmwST_local_mesh),intent(in) :: hecMESH
17  type (fstr_solid),intent(inout) :: fstrSOLID
18  type (fstr_param),intent(inout) :: fstrPARAM
19 
20  integer(kind=kint) :: grpid, ndof, ig0, ig, ityp, iS0, iE0, ik, in, idx, num
21  real(kind=kreal) :: fval, factor
22 
23  factor = fstrsolid%factor(2)
24 
25  ndof = hecmat%NDOF
26  do ig0= 1, fstrsolid%SPRING_ngrp_tot
27  grpid = fstrsolid%SPRING_ngrp_GRPID(ig0)
28  if( .not. fstr_isloadactive( fstrsolid, grpid, cstep ) ) cycle
29  ig= fstrsolid%SPRING_ngrp_ID(ig0)
30  ityp= fstrsolid%SPRING_ngrp_DOF(ig0)
31  fval= fstrsolid%SPRING_ngrp_val(ig0)
32  if( fval < 0.d0 ) fval = -fval*(1.d0-factor)
33 
34  is0= hecmesh%node_group%grp_index(ig-1) + 1
35  ie0= hecmesh%node_group%grp_index(ig )
36  do ik= is0, ie0
37  in = hecmesh%node_group%grp_item(ik)
38  idx = ndof**2 * (in - 1) + ndof * (ityp - 1) + ityp
39  hecmat%D(idx) = hecmat%D(idx) + fval
40  enddo
41  enddo
42 
43  end subroutine fstr_addspring
44 
45  subroutine fstr_update_ndforce_spring( cstep, hecMESH, fstrSOLID, B )
46  use m_fstr
47  integer(kind=kint), intent(in) :: cstep
48  type (hecmwST_local_mesh),intent(in) :: hecMESH
49  type (fstr_solid), intent(in) :: fstrSOLID
50  real(kind=kreal), intent(inout) :: b(:)
51  ! Local variables
52  integer(kind=kint) ndof,ig0,ig,ityp,iS0,iE0,ik,in,idx,num
53  integer(kind=kint) :: grpid
54  real(kind=kreal) :: fval, factor
55 
56  factor = fstrsolid%factor(2)
57 
58  ndof = hecmesh%n_dof
59  do ig0= 1, fstrsolid%SPRING_ngrp_tot
60  grpid = fstrsolid%SPRING_ngrp_GRPID(ig0)
61  if( .not. fstr_isloadactive( fstrsolid, grpid, cstep ) ) cycle
62  ig= fstrsolid%SPRING_ngrp_ID(ig0)
63  ityp= fstrsolid%SPRING_ngrp_DOF(ig0)
64  fval= fstrsolid%SPRING_ngrp_val(ig0)
65  if( fval < 0.d0 ) fval = -fval*(1.d0-factor)
66 
67  is0= hecmesh%node_group%grp_index(ig-1) + 1
68  ie0= hecmesh%node_group%grp_index(ig )
69  do ik= is0, ie0
70  in = hecmesh%node_group%grp_item(ik)
71  idx = ndof * (in - 1) + ityp
72  b(idx) = b(idx) - fval * ( fstrsolid%dunode( idx ) + fstrsolid%unode( idx ) )
73  enddo
74  enddo
75  end subroutine fstr_update_ndforce_spring
76 
77 end module m_fstr_spring
m_fstr::fstr_isloadactive
logical function fstr_isloadactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.f90:1038
m_fstr_spring::fstr_addspring
subroutine fstr_addspring(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
Definition: fstr_Spring.f90:12
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_fstr_spring::fstr_update_ndforce_spring
subroutine fstr_update_ndforce_spring(cstep, hecMESH, fstrSOLID, B)
Definition: fstr_Spring.f90:46
m_static_lib
This modules just summarizes all modules used in static analysis.
Definition: static_LIB.f90:6
m_fstr_spring
This module provides functions to deal with spring force.
Definition: fstr_Spring.f90:7