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