FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_EIG_setMASS.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 !-------------------------------------------------------------------------------
7 
8  use m_eigen_lib
10  use m_dynamic_mass
11  use m_fstr
12 
13 contains
14 
15  subroutine setmass(fstrSOLID, hecMESH, hecMAT, fstrEIG)
17  use m_fstr
18  use mmaterial
20  implicit none
21  type(hecmwst_matrix) :: hecmat
22  type(hecmwst_local_mesh) :: hecMESH
23  type(fstr_solid) :: fstrSOLID
24  type(fstr_eigen) :: fstrEIG
25  integer(kind=kint) :: i, iS, iE, ii, nn, j, jn, jS, k
26  integer(kind=kint) :: N, NP, NDOF
27  integer(kind=kint) :: icel, ic_type, itype, isect, ihead, sec_opt, cid
28  integer(kind=kint) :: nodLOCAL(20)
29  real(kind=kreal) :: val, surf, chkmass
30  real(kind=kreal) :: rho, thick, length
31  real(kind=kreal) :: lumped(120), mass(20*6, 20*6), ecoord(3,20)
32 
33  if(myrank == 0)then
34  write(imsg,*)"* mass matrix generation"
35  endif
36 
37  n = hecmat%N
38  np = hecmat%NP
39  ndof = hecmesh%n_dof
40 
41  allocate(fstreig%mass(np*ndof))
42  fstreig%mass = 0.0d0
43  val = 0.0d0
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  nn = hecmw_get_max_node(ic_type)
50 
51  if(hecmw_is_etype_link(ic_type)) cycle
52  if(hecmw_is_etype_patch(ic_type)) cycle
53  if(ic_type == 3414) cycle
54 
55  do icel = is, ie
56  js = hecmesh%elem_node_index(icel-1)
57  do j = 1, nn
58  nodlocal(j) = hecmesh%elem_node_item(js+j)
59  do i = 1, 3
60  ecoord(i,j) = hecmesh%node(3*(nodlocal(j)-1)+i)
61  enddo
62  enddo
63 
64  isect = hecmesh%section_ID(icel)
65  ihead = hecmesh%section%sect_R_index(isect-1)
66  cid = hecmesh%section%sect_mat_ID_item(isect)
67  sec_opt = hecmesh%section%sect_opt(isect)
68  rho = fstrsolid%materials(cid)%variables(m_density)
69  thick = fstrsolid%materials(cid)%variables(m_thick)
70 
71  lumped = 0.0d0
72  if(ic_type == 231 .or. ic_type == 232 .or. ic_type == 241 .or. ic_type == 242)then
73  call mass_c2(ic_type, nn, ecoord(1:2,1:nn), fstrsolid%elements(icel)%gausses, sec_opt, thick, mass, lumped)
74 
75  elseif(ic_type == 341 .or. ic_type == 342 .or. ic_type == 351 .or. ic_type == 352 .or. &
76  & ic_type == 361 .or. ic_type == 362 )then
77  call mass_c3(ic_type, nn, ecoord(1:3,1:nn), fstrsolid%elements(icel)%gausses, mass, lumped)
78 
79  elseif(ic_type==731 .or. ic_type==741 .or. ic_type==743) then
80  rho = fstrsolid%materials(cid)%variables(m_density)
81  thick = fstrsolid%materials(cid)%variables(m_thick)
82  call mass_shell(ic_type, nn, ecoord(1:3,1:nn), rho, thick, fstrsolid%elements(icel)%gausses, mass, lumped)
83 
84  elseif(ic_type == 761)then
85  surf = get_face3(ecoord(1:3,1:nn))
86  rho = fstrsolid%materials(cid)%variables(m_density)
87  thick = fstrsolid%materials(cid)%variables(m_thick)
88  val = surf*thick*rho/3.0d0
89 
90  elseif(ic_type == 781)then
91  surf = get_face4(ecoord(1:3,1:nn))
92  rho = fstrsolid%materials(cid)%variables(m_density)
93  thick = fstrsolid%materials(cid)%variables(m_thick)
94  val = surf*thick*rho/4.0d0
95 
96  elseif(ic_type == 611 .or. ic_type == 641)then
97  surf = hecmesh%section%sect_R_item(ihead+4)
98  length = get_length(ecoord(1:3,1:nn))
99  rho = fstrsolid%materials(cid)%variables(m_density)
100  val = 0.5d0*surf*length*rho
101 
102  elseif(ic_type == 301)then
103  surf = hecmesh%section%sect_R_item(ihead+1)
104  length = get_length(ecoord(1:3,1:nn))
105  rho = fstrsolid%materials(cid)%variables(m_density)
106  val = 0.5d0*surf*length*rho
107 
108  elseif( ic_type/100 == 5 )then !skip interface element
109  else
110  write(*,*)"** error setMASS"
111  endif
112 
113  do j = 1,nn
114  jn = nodlocal(j)
115  js = ndof*(jn-1)
116 
117  if(ic_type == 611)then
118  fstreig%mass(js+1) = fstreig%mass(js+1) + val
119  fstreig%mass(js+2) = fstreig%mass(js+2) + val
120  fstreig%mass(js+3) = fstreig%mass(js+3) + val
121  fstreig%mass(js+4) = fstreig%mass(js+4) + 0.0d0
122  fstreig%mass(js+5) = fstreig%mass(js+5) + 0.0d0
123  fstreig%mass(js+6) = fstreig%mass(js+6) + 0.0d0
124 
125  elseif(ic_type == 641)then
126  if(j == 1 .or. j == 2)then
127  fstreig%mass(js+1) = fstreig%mass(js+1) + val
128  fstreig%mass(js+2) = fstreig%mass(js+2) + val
129  fstreig%mass(js+3) = fstreig%mass(js+3) + val
130  elseif(j == 3 .or. j == 4)then
131  fstreig%mass(js+1) = fstreig%mass(js+1) + 0.0d0
132  fstreig%mass(js+2) = fstreig%mass(js+2) + 0.0d0
133  fstreig%mass(js+3) = fstreig%mass(js+3) + 0.0d0
134  end if
135 
136  elseif(ic_type == 761)then
137  if(j == 1 .or. j == 2 .or. j == 3)then
138  fstreig%mass(js+1) = fstreig%mass(js+1) + val
139  fstreig%mass(js+2) = fstreig%mass(js+2) + val
140  fstreig%mass(js+3) = fstreig%mass(js+3) + val
141  elseif(j == 4 .or. j == 5 .or. j == 6)then
142  fstreig%mass(js+1) = fstreig%mass(js+1) + 0.0d0
143  fstreig%mass(js+2) = fstreig%mass(js+2) + 0.0d0
144  fstreig%mass(js+3) = fstreig%mass(js+3) + 0.0d0
145  endif
146 
147  else if( ic_type == 781 ) then
148  if(j == 1 .or. j == 2 .or. j == 3 .or. j == 4)then
149  fstreig%mass(js+1) = fstreig%mass(js+1) + val
150  fstreig%mass(js+2) = fstreig%mass(js+2) + val
151  fstreig%mass(js+3) = fstreig%mass(js+3) + val
152  elseif(j == 5 .or. j == 6 .or. j == 7 .or. j == 8)then
153  fstreig%mass(js+1) = fstreig%mass(js+1) + 0.0d0
154  fstreig%mass(js+2) = fstreig%mass(js+2) + 0.0d0
155  fstreig%mass(js+3) = fstreig%mass(js+3) + 0.0d0
156  endif
157 
158  elseif(ic_type == 301) then
159  fstreig%mass(js+1) = fstreig%mass(js+1) + val
160  fstreig%mass(js+2) = fstreig%mass(js+2) + val
161  fstreig%mass(js+3) = fstreig%mass(js+3) + val
162 
163  else
164  do k = 1, ndof
165  fstreig%mass(js+k) = fstreig%mass(js+k) + lumped(ndof*(j-1)+k)
166  enddo
167  endif
168  enddo
169  enddo
170  enddo
171 
172  call hecmw_update_r(hecmesh, fstreig%mass, np, ndof)
173 
174  chkmass = 0.0d0
175  do i = 1, n
176  ii = ndof*(i-1) + 1
177  chkmass = chkmass + fstreig%mass(ii)
178  end do
179  call hecmw_allreduce_r1(hecmesh, chkmass, hecmw_sum)
180  fstreig%totalmass = chkmass
181 
182  if(myrank == 0)then
183  write(imsg,"(a,1pe12.5)")"** Total mass: ", chkmass
184  endif
185  end subroutine setmass
186 
187 end module m_fstr_eig_setmass
188 
mmaterial::m_density
integer(kind=kint), parameter m_density
Definition: material.f90:94
m_static_lib_shell::mass_shell
subroutine mass_shell(etype, nn, elem, rho, thick, gausses, mass, lumped)
Definition: static_LIB_shell.f90:3201
m_fstr_eig_setmass::setmass
subroutine setmass(fstrSOLID, hecMESH, hecMAT, fstrEIG)
Definition: fstr_EIG_setMASS.f90:16
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
hecmw_util::hecmw_sum
integer(kind=kint), parameter hecmw_sum
Definition: hecmw_util_f.F90:23
m_fstr::fstr_eigen
Package of data used by Lanczos eigenvalue solver.
Definition: m_fstr.f90:593
m_fstr::fstr_solid
Definition: m_fstr.f90:238
m_dynamic_mass::get_face4
real(kind=kreal) function get_face4(ecoord)
Definition: dynamic_mass.f90:220
m_dynamic_mass::get_face3
real(kind=kreal) function get_face3(ecoord)
Definition: dynamic_mass.f90:200
m_fstr_eig_setmass
Set up lumped mass matrix.
Definition: fstr_EIG_setMASS.f90:6
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_util::hecmwst_local_mesh
Definition: hecmw_util_f.F90:234
m_dynamic_mass::mass_c3
subroutine mass_c3(etype, nn, ecoord, gausses, mass, lumped, temperature)
Definition: dynamic_mass.f90:85
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
mmaterial::m_thick
integer(kind=kint), parameter m_thick
Definition: material.f90:95
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
m_eigen_lib
This modules just summarizes all modules used in eigen analysis.
Definition: eigen_LIB.f90:6
m_dynamic_mass::mass_c2
subroutine mass_c2(etype, nn, ecoord, gausses, sec_opt, thick, mass, lumped, temperature)
Definition: dynamic_mass.f90:11
m_static_lib_shell
This module ...
Definition: static_LIB_shell.f90:6
mmaterial
This module summarizes all information of material properties.
Definition: material.f90:6
m_static_get_prop
This module provides a function to fetch material properties from hecmw.
Definition: fstr_get_prop.f90:6
m_dynamic_mass
This module contains subroutines used in 3d eigen analysis for.
Definition: dynamic_mass.f90:6
m_dynamic_mass::get_length
real(kind=kreal) function get_length(ecoord)
Definition: dynamic_mass.f90:189
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110