FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
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 !-------------------------------------------------------------------------------
8 contains
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
30 end module m_heat_mat_ass_bc_fixt
m_fstr::fstr_heat
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:425
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_heat_mat_ass_bc_fixt::heat_mat_ass_bc_fixt
subroutine heat_mat_ass_bc_fixt(hecMAT, fstrHEAT, CTIME, DTIME, beta)
Definition: heat_mat_ass_bc_FIXT.f90:11
m_heat_mat_ass_bc_fixt
This module provides a subroutine for setting fixed temperature boundary conditions.
Definition: heat_mat_ass_bc_FIXT.f90:7
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