FrontISTR 5.2.0
Large-scale structural analysis program with finit element method
Loading...
Searching...
No Matches
heat_mat_ass_bc_FIXT.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!-------------------------------------------------------------------------------
8contains
9
10 subroutine heat_mat_ass_bc_fixt(hecMAT, fstrHEAT, CTIME, DTIME, beta)
11 use m_fstr
13 implicit none
14 type(fstr_heat) :: fstrHEAT
15 type(hecmwst_matrix) :: hecMAT
16 integer(kind=kint) :: ib, ii, id
17 real(kind=kreal) :: ctime, dtime, qq, beta
18 logical :: OutOfRange
19
20 do ib = 1, fstrheat%T_FIX_tot
21 ii = fstrheat%T_FIX_node(ib)
22 id = fstrheat%T_FIX_ampl(ib)
23 call heat_get_amplitude(fstrheat, id, ctime, qq, outofrange)
24
25 if(outofrange) cycle
26
27 call hecmw_mat_ass_bc(hecmat, ii, 1, fstrheat%T_FIX_VAL(ib)*qq)
28 enddo
29 end subroutine heat_mat_ass_bc_fixt
This module defined coomon data and basic structures for analysis.
Definition: m_fstr.f90:15
This moudle provide a function to get amplitude definition.
subroutine heat_get_amplitude(fstrheat, id, tt, qq, outofrange)
This module provides a subroutine for setting fixed temperature boundary conditions.
subroutine heat_mat_ass_bc_fixt(hecmat, fstrheat, ctime, dtime, beta)
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:394