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
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(:)
37 do k = 1, fstrheat%Q_SUF_tot
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)
45 if( fstrsolid%elements(icel)%elemact_flag == kelact_inactive ) cycle
48 val = fstrheat%Q_SUF_val (k) * qq
49 if( dabs(val) < 1.d-16 ) cycle
51 nn = hecmw_get_max_node(ic_type)
53 is = hecmesh%elem_node_index(icel-1)
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) )
61 if ( ic_type.eq.111 )
then
62 is = hecmesh%section%sect_R_index(isect)
63 asect = hecmesh%section%sect_R_item(is)
66 elseif( ic_type.eq.231 )
then
67 is = hecmesh%section%sect_R_index(isect)
68 thick = hecmesh%section%sect_R_item(is)
71 elseif( ic_type.eq.232 )
then
72 is = hecmesh%section%sect_R_index(isect)
73 thick = hecmesh%section%sect_R_item(is)
76 elseif( ic_type.eq.241 )
then
77 is = hecmesh%section%sect_R_index(isect)
78 thick = hecmesh%section%sect_R_item(is)
81 elseif( ic_type.eq.242 )
then
82 is = hecmesh%section%sect_R_index(isect)
83 thick = hecmesh%section%sect_R_item(is)
86 elseif( ic_type.eq.341 )
then
89 elseif( ic_type.eq.342 )
then
92 elseif( ic_type.eq.351 )
then
95 elseif( ic_type.eq.352 )
then
98 elseif( ic_type.eq.361 )
then
101 elseif( ic_type.eq.362 )
then
104 elseif( ic_type.eq.731 )
then
105 is = hecmesh%section%sect_R_index(isect)
106 thick = hecmesh%section%sect_R_item(is)
109 elseif( ic_type.eq.741 )
then
110 is = hecmesh%section%sect_R_index(isect)
111 thick = hecmesh%section%sect_R_item(is)
118 hecmat%B( nodlocal(j) ) = hecmat%B( nodlocal(j) ) - vect(j)
124 if( fstrheat%WL_tot>0 )
then
125 allocate( bbak(
size(hecmat%B)) )
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 )
133 icel = hecmesh%elem_group%grp_item(k)
134 ic_type = hecmesh%elem_type(icel)
135 isect = hecmesh%section_ID(icel)
138 nn = hecmw_get_max_node(ic_type)
139 is = hecmesh%elem_node_index(icel-1)
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) )
147 if( .not.
isweldactive(nn, xx,yy,zz, fstrheat%weldline(ig0), ctime-0.5d0*dtime) ) cycle
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
154 if( ic_type.eq.341 )
then
157 elseif( ic_type.eq.342 )
then
160 elseif( ic_type.eq.351 )
then
163 elseif( ic_type.eq.352 )
then
166 elseif( ic_type.eq.361 )
then
169 elseif( ic_type.eq.362 )
then
174 write(*,*)
'###ERROR### : Element type not supported for heat analysis'
175 write(*,*)
' ic_type = ', ic_type
176 call hecmw_abort(hecmw_comm_get_comm())
181 bbak( nodlocal(j) ) = bbak( nodlocal(j) ) - vect(j)
187 hecmat%B = hecmat%B + bbak
197 integer,
intent(in) :: nn
198 real(kind=kreal),
intent(in) :: xx(:), yy(:), zz(:)
200 real(kind=kreal),
intent(in) :: ctime
203 real(kind=kreal) :: cpos, wpos, tend
206 tend = weldline%tstart+ (weldline%n2-weldline%n1)/weldline%v
207 if( ctime<weldline%tstart .or. ctime>tend )
return
210 if( weldline%xyz==1 )
then
212 elseif( weldline%xyz==2 )
then
219 wpos = weldline%n1 + weldline%v*(ctime-weldline%tstart)
220 if( dabs(cpos-wpos)<weldline%distol )
isweldactive = .true.
This module defines common data and basic structures for analysis.
integer(kind=kint), parameter idbg
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)
-1:not relation, >1:index of coupled_node