FrontISTR  5.8.0
Large-scale structural analysis program with finit element method
fstr_solve_heat.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 
9  subroutine fstr_solve_heat(hecMESH, hecMAT, fstrSOLID, fstrRESULT, fstrPARAM, fstrHEAT)
10  use m_fstr
11  use m_heat_init
13  use m_heat_io
14  implicit none
15  integer(kind=kint) :: ISTEP
16  type(hecmwst_local_mesh) :: hecMESH
17  type(hecmwst_matrix) :: hecMAT
18  type(fstr_solid) :: fstrSOLID
19  type(hecmwst_result_data) :: fstrRESULT
20  type(fstr_param) :: fstrPARAM
21  type(fstr_heat) :: fstrHEAT
22 
23  integer(kind=kint) :: total_step, restart_step_num
24  real(kind=kreal) :: start_time, end_time
25 
26  call heat_init(hecmesh, fstrheat)
27  call heat_init_log(hecmesh)
28 
29  if(fstrheat%restart_nout < 0) then
30  call heat_input_restart(fstrheat, hecmesh, restart_step_num, total_step, start_time)
31  end_time = fstrheat%STEP_EETIME(restart_step_num)
32  if( (end_time - start_time) / end_time < 1.d-12 ) then
33  restart_step_num = restart_step_num + 1
34  start_time = 0.0d0
35  endif
36  else
37  restart_step_num = 1
38  total_step = 1
39  start_time = 0.0d0
40  endif
41 
42  do istep = restart_step_num, fstrheat%STEPtot
43  fstrheat%is_steady = 0
44  if(fstrheat%STEP_DLTIME(istep) <= 0.0d0) fstrheat%is_steady = 1
45 
46  if(hecmesh%my_rank == 0)then
47  write(imsg,"(a,i8,a,i8)")"* Current step / Total step: ", istep, "/", fstrheat%STEPtot
48  write(imsg,"(a,i8)")"* max iteration at each step: ", fstrparam%ITMAX(istep)
49  if(fstrheat%is_steady == 1)then
50  write(imsg,"(a)")"* Steady state analysis"
51  else
52  write(imsg,"(a)")"* Transient anslysis"
53  endif
54  endif
55 
56  call heat_solve_tran(hecmesh, hecmat, fstrsolid, fstrresult, fstrparam, fstrheat, istep, total_step, start_time)
57 
58  start_time = 0.0d0
59  enddo
60 
61  call heat_finalize(fstrheat)
62  end subroutine fstr_solve_heat
63 end module m_fstr_solve_heat
This module provides a function to control heat analysis.
subroutine fstr_solve_heat(hecMESH, hecMAT, fstrSOLID, fstrRESULT, fstrPARAM, fstrHEAT)
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:111
This module provides functions to initialize heat analysis.
Definition: heat_init.f90:6
subroutine heat_finalize(fstrHEAT)
Definition: heat_init.f90:63
subroutine heat_init(hecMESH, fstrHEAT)
Definition: heat_init.f90:10
subroutine heat_init_log(hecMESH)
Definition: heat_init.f90:48
This module provides a function to control heat analysis.
Definition: heat_io.f90:6
subroutine heat_input_restart(fstrHEAT, hecMESH, istep, tstep, tt)
Definition: heat_io.f90:11
This module provide a function to control nonsteady heat analysis.
subroutine heat_solve_tran(hecMESH, hecMAT, fstrSOLID, fstrRESULT, fstrPARAM, fstrHEAT, ISTEP, total_step, current_time)
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:431
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:155