FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_StiffMatrix.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  use m_fstr
9  implicit none
10 
11  private
12  public :: fstr_stiffmatrix
13 
14 contains
15 
16  !---------------------------------------------------------------------*
18  subroutine fstr_stiffmatrix( hecMESH, hecMAT, fstrSOLID, time, tincr)
19  !---------------------------------------------------------------------*
20  use m_static_lib
21  use mmechgauss
22 
23  type (hecmwst_local_mesh) :: hecmesh
24  type (hecmwst_matrix) :: hecmat
25  type (fstr_solid) :: fstrsolid
26  real(kind=kreal),intent(in) :: time
27  real(kind=kreal),intent(in) :: tincr
28 
29  type( tmaterial ), pointer :: material
30 
31  real(kind=kreal) :: stiffness(fstrsolid%max_ncon_stf*6, fstrsolid%max_ncon_stf*6)
32  integer(kind=kint) :: nodlocal(fstrsolid%max_ncon)
33  real(kind=kreal) :: tt(fstrsolid%max_ncon), ecoord(3,fstrsolid%max_ncon)
34  real(kind=kreal) :: thick
35  integer(kind=kint) :: ndof, itype, is, ie, ic_type, nn, icel, iis, i, j
36  real(kind=kreal) :: u(6,fstrsolid%max_ncon), du(6,fstrsolid%max_ncon), coords(3,3), u_prev(6,fstrsolid%max_ncon)
37  integer :: isect, ihead, cdsys_id
38 
39  ! ----- initialize
40  call hecmw_mat_clear( hecmat )
41 
42  ndof = hecmat%NDOF
43  tt(:) = 0.d0
44 
45  do itype= 1, hecmesh%n_elem_type
46  is= hecmesh%elem_type_index(itype-1) + 1
47  ie= hecmesh%elem_type_index(itype )
48  ic_type= hecmesh%elem_type_item(itype)
49  ! ----- Ignore link and patch elements
50  if (hecmw_is_etype_link(ic_type)) cycle
51  if (hecmw_is_etype_patch(ic_type)) cycle
52  ! ----- Set number of nodes
53 
54  ! ----- element loop
55  !$omp parallel default(none), &
56  !$omp& private(icel,iiS,nn,j,nodLOCAL,i,ecoord,du,u,u_prev,tt,cdsys_ID,coords, &
57  !$omp& material,thick,stiffness,isect,ihead), &
58  !$omp& shared(iS,iE,hecMESH,ndof,fstrSOLID,ic_type,hecMAT,time,tincr)
59  !$omp do
60  do icel= is, ie
61 
62  ! ----- nodal coordinate & displacement
63  iis= hecmesh%elem_node_index(icel-1)
64  nn = hecmesh%elem_node_index(icel)-iis
65  ! if( nn>150 ) stop "elemental nodes > 150!"
66 
67  do j=1,nn
68  nodlocal(j)= hecmesh%elem_node_item (iis+j)
69  do i=1, 3
70  ecoord(i,j) = hecmesh%node(3*nodlocal(j)+i-3)
71  enddo
72  do i=1,ndof
73  du(i,j) = fstrsolid%dunode(ndof*nodlocal(j)+i-ndof)
74  u(i,j) = fstrsolid%unode(ndof*nodlocal(j)+i-ndof) + du(i,j)
75  u_prev(i,j) = fstrsolid%unode(ndof*nodlocal(j)+i-ndof)
76  enddo
77  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres >0 ) &
78  tt(j)=fstrsolid%temperature( nodlocal(j) )
79  enddo
80 
81  isect = hecmesh%section_ID(icel)
82  ihead = hecmesh%section%sect_R_index(isect-1)
83  cdsys_id = hecmesh%section%sect_orien_ID(isect)
84  if( cdsys_id > 0 ) call get_coordsys(cdsys_id, hecmesh, fstrsolid, coords)
85  thick = hecmesh%section%sect_R_item(ihead+1)
86  if( getspacedimension( ic_type )==2 ) thick =1.d0
87  material => fstrsolid%elements(icel)%gausses(1)%pMaterial
88 
89  if( ic_type==241 .or. ic_type==242 .or. ic_type==231 .or. ic_type==232 .or. ic_type==2322) then
90  if( material%nlgeom_flag /= infinitesimal ) call stiffmat_abort( ic_type, 2 )
91  call stf_c2( ic_type,nn,ecoord(1:2,1:nn),fstrsolid%elements(icel)%gausses(:),thick, &
92  stiffness(1:nn*ndof,1:nn*ndof), fstrsolid%elements(icel)%iset, &
93  u(1:2,1:nn) )
94 
95  elseif ( ic_type==301 ) then
96  call stf_c1( ic_type,nn,ecoord(:,1:nn),thick,fstrsolid%elements(icel)%gausses(:), &
97  stiffness(1:nn*ndof,1:nn*ndof), u(1:3,1:nn) )
98 
99  elseif ( ic_type==361 ) then
100  if( fstrsolid%sections(isect)%elemopt361 == kel361fi ) then ! full integration element
101  call stf_c3 &
102  ( ic_type, nn, ecoord(:, 1:nn), fstrsolid%elements(icel)%gausses(:), &
103  stiffness(1:nn*ndof, 1:nn*ndof), cdsys_id, coords, time, tincr, u(1:3,1:nn), tt(1:nn) )
104  else if( fstrsolid%sections(isect)%elemopt361 == kel361bbar ) then ! B-bar element
105  call stf_c3d8bbar &
106  ( ic_type, nn, ecoord(:, 1:nn), fstrsolid%elements(icel)%gausses(:), &
107  stiffness(1:nn*ndof,1:nn*ndof), cdsys_id, coords, time, tincr, u(1:3, 1:nn), tt(1:nn) )
108  else if( fstrsolid%sections(isect)%elemopt361 == kel361ic ) then ! incompatible element
109  call stf_c3d8ic &
110  ( ic_type, nn, ecoord(:,1:nn), fstrsolid%elements(icel)%gausses(:), &
111  stiffness(1:nn*ndof, 1:nn*ndof), cdsys_id, coords, time, tincr, u(1:3,1:nn), &
112  fstrsolid%elements(icel)%aux, tt(1:nn) )
113  else if( fstrsolid%sections(isect)%elemopt361 == kel361fbar ) then ! F-bar element
114  call stf_c3d8fbar &
115  ( ic_type, nn, ecoord(:, 1:nn), fstrsolid%elements(icel)%gausses(:), &
116  stiffness(1:nn*ndof,1:nn*ndof), cdsys_id, coords, time, tincr, u(1:3, 1:nn), tt(1:nn) )
117  endif
118 
119  elseif (ic_type==341 .or. ic_type==351 .or. ic_type==342 .or. ic_type==352 .or. ic_type==362 ) then
120  if( ic_type==341 .and. fstrsolid%sections(isect)%elemopt341 == kel341sesns ) cycle ! skip smoothed fem
121  call stf_c3 &
122  ( ic_type, nn, ecoord(:, 1:nn), fstrsolid%elements(icel)%gausses(:), &
123  stiffness(1:nn*ndof, 1:nn*ndof), cdsys_id, coords, time, tincr, u(1:3,1:nn), tt(1:nn) )
124 
125  else if( ic_type == 511) then
126  call stf_connector( ic_type,nn,ecoord(:,1:nn),fstrsolid%elements(icel)%gausses(:), &
127  stiffness(1:nn*ndof,1:nn*ndof), u(1:3,1:nn), tincr, tt(1:nn) )
128 
129  else if( ic_type == 611) then
130  if( material%nlgeom_flag /= infinitesimal ) call stiffmat_abort( ic_type, 2 )
131  call stf_beam(ic_type, nn, ecoord, hecmesh%section%sect_R_item(ihead+1:), &
132  & material%variables(m_youngs), material%variables(m_poisson), stiffness(1:nn*ndof,1:nn*ndof))
133 
134  else if( ic_type == 641 ) then
135  if( material%nlgeom_flag /= infinitesimal ) call stiffmat_abort( ic_type, 2 )
136  call stf_beam_641(ic_type, nn, ecoord, fstrsolid%elements(icel)%gausses(:), &
137  & hecmesh%section%sect_R_item(ihead+1:), stiffness(1:nn*ndof,1:nn*ndof))
138 
139  else if( ( ic_type == 741 ) .or. ( ic_type == 743 ) .or. ( ic_type == 731 ) ) then
140  if( material%nlgeom_flag /= infinitesimal ) call stiffmat_abort( ic_type, 2 )
141  call stf_shell_mitc(ic_type, nn, ndof, ecoord(1:3, 1:nn), fstrsolid%elements(icel)%gausses(:), &
142  & stiffness(1:nn*ndof, 1:nn*ndof), thick, 0)
143 
144  else if( ic_type == 761 ) then !for shell-solid mixed analysis
145  if( material%nlgeom_flag /= infinitesimal ) call stiffmat_abort( ic_type, 2 )
146  call stf_shell_mitc(731, 3, 6, ecoord(1:3, 1:3), fstrsolid%elements(icel)%gausses(:), &
147  & stiffness(1:nn*ndof, 1:nn*ndof), thick, 2)
148 
149  else if( ic_type == 781 ) then !for shell-solid mixed analysis
150  if( material%nlgeom_flag /= infinitesimal ) call stiffmat_abort( ic_type, 2 )
151  call stf_shell_mitc(741, 4, 6, ecoord(1:3, 1:4), fstrsolid%elements(icel)%gausses(:), &
152  & stiffness(1:nn*ndof, 1:nn*ndof), thick, 1)
153 
154  elseif ( ic_type==3414 ) then
155  if( material%mtype /= incomp_newtonian) call stiffmat_abort( ic_type, 3, material%mtype )
156  call stf_c3_vp &
157  ( ic_type, nn, ecoord(:, 1:nn),fstrsolid%elements(icel)%gausses(:), &
158  stiffness(1:nn*ndof, 1:nn*ndof), tincr, u_prev(1:4, 1:nn) )
159  else if ( ic_type == 881 .or. ic_type == 891 ) then !for selective es/ns smoothed fem
160  call stf_c3d4_sesns &
161  ( ic_type,nn,nodlocal,ecoord(:, 1:nn), fstrsolid%elements(icel)%gausses(:), &
162  stiffness, cdsys_id, coords, time, tincr, u(1:3,1:nn), tt(1:nn) )
163  else
164  call stiffmat_abort( ic_type, 1 )
165  endif
166  !
167  ! ----- CONSTRUCT the GLOBAL MATRIX STARTED
168  call hecmw_mat_ass_elem(hecmat, nn, nodlocal, stiffness)
169 
170  enddo ! icel
171  !$omp end do
172  !$omp end parallel
173  enddo ! itype
174 
175  end subroutine fstr_stiffmatrix
176 
177  subroutine stiffmat_abort( ic_type, flag, mtype )
178  integer(kind=kint), intent(in) :: ic_type
179  integer(kind=kint), intent(in) :: flag
180  integer(kind=kint), intent(in), optional :: mtype
181 
182  if( flag == 1 ) then
183  write(*,*) '###ERROR### : Element type not supported for static analysis'
184  else if( flag == 2 ) then
185  write(*,*) '###ERROR### : Element type not supported for nonlinear static analysis'
186  else if( flag == 3 ) then
187  write(*,*) '###ERROR### : This element is not supported for this material'
188  endif
189  write(*,*) ' ic_type = ', ic_type
190  if( present(mtype) ) write(*,*) ' mtype = ', mtype
191  call hecmw_abort(hecmw_comm_get_comm())
192  end subroutine
193 
194 end module m_fstr_stiffmatrix
m_fstr::kel361fi
integer(kind=kint), parameter kel361fi
Definition: m_fstr.f90:77
m_fstr::kel361bbar
integer(kind=kint), parameter kel361bbar
Definition: m_fstr.f90:78
m_fstr_stiffmatrix::stiffmat_abort
subroutine stiffmat_abort(ic_type, flag, mtype)
Definition: fstr_StiffMatrix.f90:178
m_fstr_stiffmatrix::fstr_stiffmatrix
subroutine, public fstr_stiffmatrix(hecMESH, hecMAT, fstrSOLID, time, tincr)
This subroutine creates tangential stiffness matrix.
Definition: fstr_StiffMatrix.f90:19
m_fstr::fstr_solid
Definition: m_fstr.f90:238
mmechgauss
This modules defines a structure to record history dependent parameter in static analysis.
Definition: mechgauss.f90:6
m_fstr::kel341sesns
integer(kind=kint), parameter kel341sesns
Definition: m_fstr.f90:75
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_fstr::get_coordsys
subroutine get_coordsys(cdsys_ID, hecMESH, fstrSOLID, coords)
This subroutine fetch coords defined by local coordinate system.
Definition: m_fstr.f90:1073
m_fstr::kel361fbar
integer(kind=kint), parameter kel361fbar
Definition: m_fstr.f90:80
m_fstr::kel361ic
integer(kind=kint), parameter kel361ic
Definition: m_fstr.f90:79
m_fstr_stiffmatrix
This module provides function to calculate tangent stiffness matrix.
Definition: fstr_StiffMatrix.f90:7
m_static_lib
This modules just summarizes all modules used in static analysis.
Definition: static_LIB.f90:6