FrontISTR  5.8.0
Large-scale structural analysis program with finit element method
heat_mat_ass_bc_DFLUX.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  use m_fstr
9 
10  implicit none
11 contains
12  !C
13  !C***
14  !C*** MAT_ASS_DFLUX
15  !C***
16  !C
17  subroutine heat_mat_ass_bc_dflux( hecMESH, hecMAT, fstrSOLID, fstrHEAT, CTIME, DTIME, beta )
20  use m_static_lib_3d
21  integer(kind=kint) :: k, icel, ic_type, isect, isuf, iamp, nn, is, j
22  real(kind=kreal) :: ctime, dtime, qq, val, asect, thick, vol, beta
23  type(fstr_solid) :: fstrsolid
24  type(fstr_heat) :: fstrHEAT
25  type(hecmwst_matrix) :: hecMAT
26  type(hecmwst_local_mesh) :: hecMESH
27  real(kind=kreal) :: xx(20), yy(20), zz(20)
28  real(kind=kreal) :: vect(20), ss(2000)
29  integer(kind=kint) :: ig0, ig, iS0, iE0, nodLocal(20)
30  real(kind=kreal), allocatable :: bbak(:)
31 
32  !C
33  !$omp parallel default(none), &
34  !$omp& private(k,icel,ic_type,isect,isuf,iamp,QQ,val,nn,is,j,nodLOCAL,xx,yy,zz,asect,vect,thick), &
35  !$omp& shared(fstrSOLID,fstrHEAT,CTIME,hecMAT,hecMESH)
36  !$omp do
37  do k = 1, fstrheat%Q_SUF_tot
38 
39  icel = fstrheat%Q_SUF_elem(k)
40  ic_type = hecmesh%elem_type(icel)
41  isect = hecmesh%section_ID(icel)
42  isuf = fstrheat%Q_SUF_surf(k)
43  iamp = fstrheat%Q_SUF_ampl(k)
44 
45  if( fstrsolid%elements(icel)%elemact_flag == kelact_inactive ) cycle
46 
47  call heat_get_amplitude( fstrheat,iamp,ctime,qq )
48  val = fstrheat%Q_SUF_val (k) * qq
49  if( dabs(val) < 1.d-16 ) cycle
50  !C**
51  nn = hecmw_get_max_node(ic_type)
52  !C**
53  is = hecmesh%elem_node_index(icel-1)
54  do j = 1, nn
55  nodlocal(j) = hecmesh%elem_node_item(is+j)
56  xx(j) = hecmesh%node( 3*nodlocal(j)-2 )
57  yy(j) = hecmesh%node( 3*nodlocal(j)-1 )
58  zz(j) = hecmesh%node( 3*nodlocal(j) )
59  enddo
60  !C**
61  if ( ic_type.eq.111 ) then
62  is = hecmesh%section%sect_R_index(isect)
63  asect = hecmesh%section%sect_R_item(is)
64  call heat_dflux_111(nn,xx,yy,zz,asect,isuf,val,vect)
65 
66  elseif( ic_type.eq.231 ) then
67  is = hecmesh%section%sect_R_index(isect)
68  thick = hecmesh%section%sect_R_item(is)
69  call heat_dflux_231(nn,xx,yy,zz,thick,isuf,val,vect)
70 
71  elseif( ic_type.eq.232 ) then
72  is = hecmesh%section%sect_R_index(isect)
73  thick = hecmesh%section%sect_R_item(is)
74  call heat_dflux_232(nn,xx,yy,zz,thick,isuf,val,vect)
75 
76  elseif( ic_type.eq.241 ) then
77  is = hecmesh%section%sect_R_index(isect)
78  thick = hecmesh%section%sect_R_item(is)
79  call heat_dflux_241(nn,xx,yy,zz,thick,isuf,val,vect)
80 
81  elseif( ic_type.eq.242 ) then
82  is = hecmesh%section%sect_R_index(isect)
83  thick = hecmesh%section%sect_R_item(is)
84  call heat_dflux_242(nn,xx,yy,zz,thick,isuf,val,vect)
85 
86  elseif( ic_type.eq.341 ) then
87  call heat_dflux_341(nn,xx,yy,zz,isuf,val,vect)
88 
89  elseif( ic_type.eq.342 ) then
90  call heat_dflux_342(nn,xx,yy,zz,isuf,val,vect)
91 
92  elseif( ic_type.eq.351 ) then
93  call heat_dflux_351(nn,xx,yy,zz,isuf,val,vect)
94 
95  elseif( ic_type.eq.352 ) then
96  call heat_dflux_352(nn,xx,yy,zz,isuf,val,vect)
97 
98  elseif( ic_type.eq.361 ) then
99  call heat_dflux_361(nn,xx,yy,zz,isuf,val,vect)
100 
101  elseif( ic_type.eq.362 ) then
102  call heat_dflux_362(nn,xx,yy,zz,isuf,val,vect)
103 
104  elseif( ic_type.eq.731 ) then
105  is = hecmesh%section%sect_R_index(isect)
106  thick = hecmesh%section%sect_R_item(is)
107  call heat_dflux_731(nn,xx,yy,zz,thick,isuf,val,vect)
108 
109  elseif( ic_type.eq.741 ) then
110  is = hecmesh%section%sect_R_index(isect)
111  thick = hecmesh%section%sect_R_item(is)
112  call heat_dflux_741(nn,xx,yy,zz,thick,isuf,val,vect)
113 
114  endif
115 
116  do j = 1, nn
117  !$omp atomic
118  hecmat%B( nodlocal(j) ) = hecmat%B( nodlocal(j) ) - vect(j)
119  enddo
120  enddo
121  !$omp end do
122  !$omp end parallel
123 
124  if( fstrheat%WL_tot>0 ) then
125  allocate( bbak(size(hecmat%B)) )
126 
127  do ig0 = 1, fstrheat%WL_tot
128  ig = fstrheat%weldline(ig0)%egrpid
129  is0= hecmesh%elem_group%grp_index(ig-1) + 1
130  ie0= hecmesh%elem_group%grp_index(ig )
131  vol =0.d0; bbak=0.d0
132  do k=is0, ie0
133  icel = hecmesh%elem_group%grp_item(k)
134  ic_type = hecmesh%elem_type(icel)
135  isect = hecmesh%section_ID(icel)
136  isuf = 0
137 
138  nn = hecmw_get_max_node(ic_type)
139  is = hecmesh%elem_node_index(icel-1)
140  do j = 1, nn
141  nodlocal(j) = hecmesh%elem_node_item(is+j)
142  xx(j) = hecmesh%node( 3*nodlocal(j)-2 )
143  yy(j) = hecmesh%node( 3*nodlocal(j)-1 )
144  zz(j) = hecmesh%node( 3*nodlocal(j) )
145  enddo
146 
147  if( .not. isweldactive(nn, xx,yy,zz, fstrheat%weldline(ig0), ctime-0.5d0*dtime) ) cycle
148 
149  vol = vol + volume_c3(ic_type,nn,xx,yy,zz)
150  val = fstrheat%weldline(ig0)%I * fstrheat%weldline(ig0)%U * &
151  fstrheat%weldline(ig0)%coe
152  write(idbg,*) "Element", hecmesh%global_elem_id(icel),"with dflux", val, vol
153 
154  if( ic_type.eq.341 ) then
155  call heat_dflux_341(nn,xx,yy,zz,isuf,val,vect)
156 
157  elseif( ic_type.eq.342 ) then
158  call heat_dflux_342(nn,xx,yy,zz,isuf,val,vect)
159 
160  elseif( ic_type.eq.351 ) then
161  call heat_dflux_351(nn,xx,yy,zz,isuf,val,vect)
162 
163  elseif( ic_type.eq.352 ) then
164  call heat_dflux_352(nn,xx,yy,zz,isuf,val,vect)
165 
166  elseif( ic_type.eq.361 ) then
167  call heat_dflux_361(nn,xx,yy,zz,isuf,val,vect)
168 
169  elseif( ic_type.eq.362 ) then
170  call heat_dflux_362(nn,xx,yy,zz,isuf,val,vect)
171 
172 
173  else
174  write(*,*) '###ERROR### : Element type not supported for heat analysis'
175  write(*,*) ' ic_type = ', ic_type
176  call hecmw_abort(hecmw_comm_get_comm())
177 
178  endif
179 
180  do j = 1, nn
181  bbak( nodlocal(j) ) = bbak( nodlocal(j) ) - vect(j)
182  enddo
183  enddo
184 
185  if( vol>0 ) then
186  bbak = bbak/vol
187  hecmat%B = hecmat%B + bbak
188  endif
189  enddo
190 
191  deallocate(bbak)
192  endif
193 
194  end subroutine heat_mat_ass_bc_dflux
195 
196  logical function isweldactive(nn, xx,yy,zz, weldline, ctime)
197  integer, intent(in) :: nn
198  real(kind=kreal), intent(in) :: xx(:), yy(:), zz(:)
199  type(tweldline), intent(in) :: weldline
200  real(kind=kreal), intent(in) :: ctime
201 
202  integer :: i
203  real(kind=kreal) :: cpos, wpos, tend
204 
205  isweldactive = .false.
206  tend = weldline%tstart+ (weldline%n2-weldline%n1)/weldline%v
207  if( ctime<weldline%tstart .or. ctime>tend ) return
208  cpos = 0
209  do i=1,nn
210  if( weldline%xyz==1 ) then
211  cpos = cpos+xx(i)
212  elseif( weldline%xyz==2 ) then
213  cpos = cpos+yy(i)
214  else
215  cpos = cpos+zz(i)
216  endif
217  enddo
218  cpos = cpos/nn
219  wpos = weldline%n1 + weldline%v*(ctime-weldline%tstart)
220  if( dabs(cpos-wpos)<weldline%distol ) isweldactive = .true.
221  end function
222 end module m_heat_mat_ass_bc_dflux
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
integer(kind=kint), parameter idbg
Definition: m_fstr.f90:112
This module provide a function to get amplitude definition.
subroutine heat_get_amplitude(fstrHEAT, id, TT, QQ, OutOfRange)
This module provides subroutines for calculating distributed heat flux for various elements.
subroutine heat_dflux_242(NN, XX, YY, ZZ, THICK, LTYPE, val, VECT)
subroutine heat_dflux_351(NN, XX, YY, ZZ, LTYPE, val, VECT)
subroutine heat_dflux_341(NN, XX, YY, ZZ, LTYPE, val, VECT)
subroutine heat_dflux_352(NN, XX, YY, ZZ, LTYPE, val, VECT)
subroutine heat_dflux_231(NN, XX, YY, ZZ, THICK, LTYPE, val, VECT)
subroutine heat_dflux_741(NN, XX, YY, ZZ, THICK, LTYPE, val, VECT)
subroutine heat_dflux_111(NN, XX, YY, ZZ, ASECT, LTYPE, val, VECT)
subroutine heat_dflux_342(NN, XX, YY, ZZ, LTYPE, val, VECT)
subroutine heat_dflux_232(NN, XX, YY, ZZ, THICK, LTYPE, val, VECT)
subroutine heat_dflux_731(NN, XX, YY, ZZ, THICK, LTYPE, val, VECT)
subroutine heat_dflux_362(NN, XX, YY, ZZ, LTYPE, val, VECT)
subroutine heat_dflux_361(NN, XX, YY, ZZ, LTYPE, val, VECT)
subroutine heat_dflux_241(NN, XX, YY, ZZ, THICK, LTYPE, val, VECT)
This module provides a subroutine for setting distributed heat flux boundary conditions.
subroutine heat_mat_ass_bc_dflux(hecMESH, hecMAT, fstrSOLID, fstrHEAT, CTIME, DTIME, beta)
logical function isweldactive(nn, xx, yy, zz, weldline, ctime)
This module provides common functions of Solid elements.
real(kind=kreal) function volume_c3(etype, nn, XX, YY, ZZ)
Volume of element.
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:431
-1:not relation, >1:index of coupled_node
Definition: m_fstr.f90:637