FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
heat_get_amplitude.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 !-------------------------------------------------------------------------------
7 contains
8  !C***
9  !C*** GET_AMPLITUDE
10  !C***
11  subroutine heat_get_amplitude ( fstrHEAT,id,TT,QQ,OutOfRange )
12 
13  use m_fstr
14 
15  implicit none
16  integer(kind=kint) :: id, nn, ii, ikk
17  real(kind=kreal) :: tt, qq
18  type(fstr_heat) :: fstrheat
19  logical, optional :: outofrange
20 
21  qq = 1.0
22  if (present(outofrange)) outofrange = .false.
23  if( id>0 ) then
24  nn = fstrheat%AMPLtab(id)
25  if ( tt < fstrheat%AMPLtime(id,1) ) then
26  ii = 1
27  if (present(outofrange)) outofrange = .true.
28  elseif( tt >= fstrheat%AMPLtime(id,nn) ) then
29  ii = nn + 1
30  if (present(outofrange)) outofrange = .true.
31  else
32  ii = 2
33  if (present(outofrange)) outofrange = .false.
34  do ikk= 1, nn - 1
35  if( tt .ge. fstrheat%AMPLtime(id,ikk) &
36  .and. tt .lt. fstrheat%AMPLtime(id,ikk+1) ) then
37  ii = ikk + 1
38  exit
39  endif
40  enddo
41  endif
42  qq = fstrheat%AMPLfuncA(id,ii) * tt + fstrheat%AMPLfuncB(id,ii)
43  endif
44 
45  end subroutine heat_get_amplitude
46 end module m_heat_get_amplitude
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_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