FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
heat_mat_ass_capacity.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 contains
8 
9  subroutine heat_mat_ass_capacity(hecMESH, hecMAT, fstrHEAT, delta_time)
10  use m_fstr
12  implicit none
13  type(fstr_heat) :: fstrheat
14  type(hecmwst_matrix) :: hecMAT
15  type(hecmwst_local_mesh) :: hecMESH
16  integer(kind=kint) :: i, in, j, nodLOCAL(20), ip, inod
17  integer(kind=kint) :: itype, iS, iE, ic_type, icel, isect, IMAT, in0, nn
18  real(kind=kreal) :: temp(20), lumped(120), mass(20*6, 20*6), ecoord(3,20)
19  real(kind=kreal) :: delta_time, surf, thick, alfa, beta
20 
21  beta = fstrheat%beta
22 
23  do itype = 1, hecmesh%n_elem_type
24  is = hecmesh%elem_type_index(itype-1) + 1
25  ie = hecmesh%elem_type_index(itype )
26  ic_type = hecmesh%elem_type_item(itype)
27 
28  if (hecmw_is_etype_link(ic_type)) cycle
29  if (hecmw_is_etype_patch(ic_type)) cycle
30  if(ic_type == 3414) cycle
31 
32  !$omp parallel default(none), &
33  !$omp& private(icel,isect,IMAT,nn,temp,in0,i,j,nodLOCAL,ecoord,in,thick,surf,inod,lumped,mass), &
34  !$omp& shared(iS,iE,hecMESH,ic_type,hecMAT,fstrHEAT,delta_time)
35  !$omp do
36  do icel = is, ie
37  isect = hecmesh%section_ID(icel)
38  imat = hecmesh%section%sect_mat_ID_item(isect)
39  in0 = hecmesh%elem_node_index(icel-1)
40  nn = hecmw_get_max_node(ic_type)
41  do i = 1, nn
42  nodlocal(i) = hecmesh%elem_node_item(in0+i)
43  do j = 1, 3
44  ecoord(j,i) = hecmesh%node(3*(nodlocal(i)-1)+j)
45  enddo
46  temp(i) = fstrheat%TEMP0(nodlocal(i))
47  enddo
48 
49  lumped = 0.0d0
50  if(ic_type == 231 .or. ic_type == 232 .or. ic_type == 241 .or. ic_type == 242)then
51  in = hecmesh%section%sect_R_index(isect)
52  thick = hecmesh%section%sect_R_item(in)
53  call heat_capacity_c2(ic_type, nn, ecoord(1:2,1:nn), temp, imat, thick, lumped, mass, &
54  fstrheat%CPtab(imat), fstrheat%CPtemp(imat,:), fstrheat%CPfuncA(imat,:) ,fstrheat%CPfuncB(imat,:), &
55  fstrheat%RHOtab(imat), fstrheat%RHOtemp(imat,:), fstrheat%RHOfuncA(imat,:), fstrheat%RHOfuncB(imat,:))
56 
57  elseif(ic_type == 341 .or. ic_type == 342 .or. ic_type == 351 .or. ic_type == 352 .or. &
58  & ic_type == 361 .or. ic_type == 362)then
59  call heat_capacity_c3(ic_type, nn, ecoord(1:3,1:nn), temp, imat, lumped, mass, &
60  fstrheat%CPtab(imat), fstrheat%CPtemp(imat,:), fstrheat%CPfuncA(imat,:) ,fstrheat%CPfuncB(imat,:), &
61  fstrheat%RHOtab(imat), fstrheat%RHOtemp(imat,:), fstrheat%RHOfuncA(imat,:), fstrheat%RHOfuncB(imat,:))
62 
63  elseif (ic_type == 541) then
64  ! skip warning
65 
66  elseif(ic_type == 731)then
67  in = hecmesh%section%sect_R_index(isect)
68  thick = hecmesh%section%sect_R_item(in)
69  call heat_capacity_shell_731(ic_type, nn, ecoord(1:3,1:nn), temp, imat, thick, lumped, mass, &
70  fstrheat%CPtab(imat), fstrheat%CPtemp(imat,:), fstrheat%CPfuncA(imat,:) ,fstrheat%CPfuncB(imat,:), &
71  fstrheat%RHOtab(imat), fstrheat%RHOtemp(imat,:), fstrheat%RHOfuncA(imat,:), fstrheat%RHOfuncB(imat,:))
72 
73  elseif(ic_type == 741)then
74  in = hecmesh%section%sect_R_index(isect)
75  thick = hecmesh%section%sect_R_item(in)
76  call heat_capacity_shell_741(ic_type, nn, ecoord(1:3,1:nn), temp, imat, thick, lumped, mass, &
77  fstrheat%CPtab(imat), fstrheat%CPtemp(imat,:), fstrheat%CPfuncA(imat,:) ,fstrheat%CPfuncB(imat,:), &
78  fstrheat%RHOtab(imat), fstrheat%RHOtemp(imat,:), fstrheat%RHOfuncA(imat,:), fstrheat%RHOfuncB(imat,:))
79 
80  elseif(ic_type == 111 .or. ic_type == 301 .or. ic_type == 611 .or. ic_type == 641)then
81  in = hecmesh%section%sect_R_index(isect)
82  surf = hecmesh%section%sect_R_item(in)
83  call heat_capacity_c1(ic_type, nn, ecoord(1:3,1:nn), temp, imat, surf, lumped, mass, &
84  fstrheat%CPtab(imat), fstrheat%CPtemp(imat,:), fstrheat%CPfuncA(imat,:) ,fstrheat%CPfuncB(imat,:), &
85  fstrheat%RHOtab(imat), fstrheat%RHOtemp(imat,:), fstrheat%RHOfuncA(imat,:), fstrheat%RHOfuncB(imat,:))
86 
87  else
88  write(*,*)"** error setMASS"
89  endif
90 
91  do ip = 1, nn
92  inod = nodlocal(ip)
93  !$omp atomic
94  hecmat%D(inod) = hecmat%D(inod) + lumped(ip) / delta_time
95  !$omp atomic
96  hecmat%B(inod) = hecmat%B(inod) + lumped(ip)*temp(ip) / delta_time
97  enddo
98  enddo
99  !$omp end do
100  !$omp end parallel
101  enddo
102  end subroutine heat_mat_ass_capacity
103 end module m_heat_mat_ass_capacity
m_heat_mat_ass_capacity::heat_mat_ass_capacity
subroutine heat_mat_ass_capacity(hecMESH, hecMAT, fstrHEAT, delta_time)
Definition: heat_mat_ass_capacity.f90:10
m_heat_lib_capacity::heat_capacity_c2
subroutine heat_capacity_c2(etype, nn, ecoord, temperature, IMAT, thick, lumped, mass, ntab1, temp1, funcA1, funcB1, ntab2, temp2, funcA2, funcB2)
Definition: heat_LIB_CAPACITY.f90:74
m_fstr::fstr_heat
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:425
m_heat_lib_capacity::heat_capacity_shell_741
subroutine heat_capacity_shell_741(etype, nn, ecoord, temperature, IMAT, thick, lumped, mass, ntab1, temp1, funcA1, funcB1, ntab2, temp2, funcA2, funcB2)
Definition: heat_LIB_CAPACITY.f90:232
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_heat_lib_capacity::heat_capacity_c3
subroutine heat_capacity_c3(etype, nn, ecoord, temperature, IMAT, lumped, mass, ntab1, temp1, funcA1, funcB1, ntab2, temp2, funcA2, funcB2)
Definition: heat_LIB_CAPACITY.f90:132
m_heat_lib_capacity
Definition: heat_LIB_CAPACITY.f90:5
m_heat_lib_capacity::heat_capacity_shell_731
subroutine heat_capacity_shell_731(etype, nn, ecoord, temperature, IMAT, thick, lumped, mass, ntab1, temp1, funcA1, funcB1, ntab2, temp2, funcA2, funcB2)
Definition: heat_LIB_CAPACITY.f90:191
m_heat_mat_ass_capacity
This module provides a subroutine to assemble heat capacity matrix.
Definition: heat_mat_ass_capacity.f90:6
m_heat_lib_capacity::heat_capacity_c1
subroutine heat_capacity_c1(etype, nn, ecoord, temperature, IMAT, surf, lumped, mass, ntab1, temp1, funcA1, funcB1, ntab2, temp2, funcA2, funcB2)
Definition: heat_LIB_CAPACITY.f90:38