FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
heat_mat_ass_bc_FILM.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 !-------------------------------------------------------------------------------
8 contains
9  !C
10  !C***
11  !C*** MAT_ASS_FILM
12  !C***
13  !C
14  subroutine heat_mat_ass_bc_film( hecMESH, hecMAT, fstrHEAT, CTIME, DTIME, beta )
15 
16  use m_fstr
18  use m_heat_lib_film
19 
20  implicit none
21  integer(kind=kint) :: k, icel, isuf, iam1, iam2, ic_type, isect, nn, is, j, mm, m, ic, ip
22  integer(kind=kint) :: inod, jp, jnod, isu, ieu, ik, isl, iel
23  real(kind=kreal) :: ctime, dtime, qq, hh, sink, thick, beta
24  type(fstr_heat) :: fstrheat
25  type(hecmwst_matrix) :: hecMAT
26  type(hecmwst_local_mesh) :: hecMESH
27 
28  real(kind=kreal) :: xx(20), yy(20), zz(20)
29  real(kind=kreal) :: term1(64), term2(20), stiff(8,8)
30  integer(kind=kint) :: nodLocal(20), nsuf(8), nodSurf(8)
31 
32 !C
33  !$omp parallel default(none), &
34  !$omp& private(k,icel,isuf,iam1,iam2,QQ,HH,SINK,ic_type,isect,nn,is,j,nodLocal, &
35  !$omp& xx,yy,zz,thick,mm,term1,term2,stiff,nsuf,nodSurf,ip,inod,jnod,ic,isU,ieU,ik,jp,isL,ieL), &
36  !$omp& shared(fstrHEAT,CTIME,hecMAT,hecMESH)
37  !$omp do
38  do k = 1, fstrheat%H_SUF_tot
39  icel = fstrheat%H_SUF_elem(k)
40  isuf = fstrheat%H_SUF_surf(k)
41  iam1 = fstrheat%H_SUF_ampl(k,1)
42  iam2 = fstrheat%H_SUF_ampl(k,2)
43  call heat_get_amplitude ( fstrheat, iam1, ctime, qq )
44  hh = fstrheat%H_SUF_val (k,1) * qq
45  call heat_get_amplitude ( fstrheat, iam2, ctime, qq )
46  sink = fstrheat%H_SUF_val (k,2) * qq
47  ic_type = hecmesh%elem_type(icel)
48  isect = hecmesh%section_ID(icel)
49  !C**
50  nn = hecmw_get_max_node(ic_type)
51  !C**
52  is = hecmesh%elem_node_index(icel-1)
53  do j = 1, nn
54  nodlocal(j) = hecmesh%elem_node_item(is+j)
55  xx(j) = hecmesh%node( 3*nodlocal(j)-2 )
56  yy(j) = hecmesh%node( 3*nodlocal(j)-1 )
57  zz(j) = hecmesh%node( 3*nodlocal(j) )
58  enddo
59  !C**
60  if ( ic_type.eq.231 ) then
61  is = hecmesh%section%sect_R_index(isect)
62  thick = hecmesh%section%sect_R_item(is)
63  mm = 2
64  call heat_film_231(nn,xx,yy,zz,thick,isuf,hh,sink,mm,term1,term2,nsuf)
65  elseif( ic_type.eq.232 ) then
66  is = hecmesh%section%sect_R_index(isect)
67  thick = hecmesh%section%sect_R_item(is)
68  mm = 3
69  call heat_film_232(nn,xx,yy,zz,thick,isuf,hh,sink,mm,term1,term2,nsuf)
70  elseif( ic_type.eq.241 ) then
71  is = hecmesh%section%sect_R_index(isect)
72  thick = hecmesh%section%sect_R_item(is)
73  mm = 2
74  call heat_film_241(nn,xx,yy,zz,thick,isuf,hh,sink,mm,term1,term2,nsuf)
75  elseif( ic_type.eq.242 ) then
76  is = hecmesh%section%sect_R_index(isect)
77  thick = hecmesh%section%sect_R_item(is)
78  mm = 3
79  call heat_film_242(nn,xx,yy,zz,thick,isuf,hh,sink,mm,term1,term2,nsuf)
80  elseif( ic_type.eq.341 ) then
81  mm = 3
82  call heat_film_341(nn,xx,yy,zz,isuf,hh,sink,mm,term1,term2,nsuf)
83  elseif( ic_type.eq.342 ) then
84  mm = 6
85  call heat_film_342(nn,xx,yy,zz,isuf,hh,sink,mm,term1,term2,nsuf)
86  elseif( ic_type.eq.351 ) then
87  mm = 4
88  if( isuf.eq.1 .or. isuf.eq.2 ) mm = 3
89  call heat_film_351(nn,xx,yy,zz,isuf,hh,sink,mm,term1,term2,nsuf)
90  elseif( ic_type.eq.352 ) then
91  mm = 8
92  if( isuf.eq.1 .or. isuf.eq.2 ) mm = 6
93  call heat_film_352(nn,xx,yy,zz,isuf,hh,sink,mm,term1,term2,nsuf)
94  elseif( ic_type.eq.361 ) then
95  mm = 4
96  call heat_film_361(nn,xx,yy,zz,isuf,hh,sink,mm,term1,term2,nsuf)
97  elseif( ic_type.eq.362 ) then
98  mm = 8
99  call heat_film_362(nn,xx,yy,zz,isuf,hh,sink,mm,term1,term2,nsuf)
100  elseif( ic_type.eq.731 ) then
101  mm = 3
102  call heat_film_731(nn,xx,yy,zz,hh,sink,term1,term2)
103  do m = 1, mm
104  nsuf(m) = m
105  enddo
106 
107  elseif( ic_type.eq.741 ) then
108  mm = 4
109  call heat_film_741(nn,xx,yy,zz,hh,sink,term1,term2)
110  do m = 1, mm
111  nsuf(m) = m
112  enddo
113 
114  endif
115  !C
116 
117  do ip = 1, mm
118  nodsurf(ip) = nodlocal(nsuf(ip))
119  end do
120 
121  ic = 0
122  stiff = 0.d0
123  do ip = 1, mm
124  do jp = 1, mm
125  ic = ic + 1
126  stiff(ip,jp) = -term1(ic)
127  enddo
128  enddo
129 
130  call hecmw_mat_ass_elem(hecmat, mm, nodsurf, stiff)
131 
132  do ip = 1, mm
133  !$omp atomic
134  hecmat%B(nodsurf(ip)) = hecmat%B(nodsurf(ip)) - term2(ip)
135  end do
136 
137  !C
138  !C
139  enddo
140  !$omp end do
141  !$omp end parallel
142 
143  end subroutine heat_mat_ass_bc_film
144 end module m_heat_mat_ass_bc_film
m_heat_lib_film::heat_film_741
subroutine heat_film_741(NN, XX, YY, ZZ, HH, SINK, TERM1, TERM2)
Definition: heat_LIB_FILM.f90:1245
m_heat_mat_ass_bc_film
This module provides a subroutine for setting heat transfer boundary conditions.
Definition: heat_mat_ass_bc_FILM.f90:7
m_heat_lib_film::heat_film_361
subroutine heat_film_361(NN, XX, YY, ZZ, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:925
m_heat_lib_film::heat_film_731
subroutine heat_film_731(NN, XX, YY, ZZ, HH, SINK, TERM1, TERM2)
Definition: heat_LIB_FILM.f90:1209
m_fstr::fstr_heat
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:425
m_heat_lib_film::heat_film_351
subroutine heat_film_351(NN, XX, YY, ZZ, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:537
m_heat_lib_film::heat_film_341
subroutine heat_film_341(NN, XX, YY, ZZ, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:317
m_heat_lib_film::heat_film_231
subroutine heat_film_231(NN, XX, YY, ZZ, THICK, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:24
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_heat_lib_film
This module provides subroutines to generate heat transfer boundary.
Definition: heat_LIB_FILM.f90:7
m_heat_lib_film::heat_film_352
subroutine heat_film_352(NN, XX, YY, ZZ, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:677
m_heat_lib_film::heat_film_241
subroutine heat_film_241(NN, XX, YY, ZZ, THICK, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:166
m_heat_lib_film::heat_film_362
subroutine heat_film_362(NN, XX, YY, ZZ, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:1042
m_heat_mat_ass_bc_film::heat_mat_ass_bc_film
subroutine heat_mat_ass_bc_film(hecMESH, hecMAT, fstrHEAT, CTIME, DTIME, beta)
Definition: heat_mat_ass_bc_FILM.f90:15
m_heat_lib_film::heat_film_242
subroutine heat_film_242(NN, XX, YY, ZZ, THICK, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:235
m_heat_lib_film::heat_film_232
subroutine heat_film_232(NN, XX, YY, ZZ, THICK, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:89
m_heat_lib_film::heat_film_342
subroutine heat_film_342(NN, XX, YY, ZZ, LTYPE, HH, SINK, MM, TERM1, TERM2, NOD)
Definition: heat_LIB_FILM.f90:378
m_heat_get_amplitude
This module provide a function to get amplitude definition.
Definition: heat_get_amplitude.f90:6
m_heat_get_amplitude::heat_get_amplitude
subroutine heat_get_amplitude(fstrHEAT, id, TT, QQ, OutOfRange)
Definition: heat_get_amplitude.f90:12