FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
heat_solve_TRAN.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 heat_solve_tran ( hecMESH,hecMAT,fstrRESULT,fstrPARAM,fstrHEAT,ISTEP,total_step,current_time )
10  use m_fstr
14  use m_heat_init
16  use m_solve_lineq
17  use m_heat_io
18  implicit none
19  integer(kind=kint) :: ISTEP, iterALL, i, inod, mnod
20  integer(kind=kint) :: total_step
21  real(kind=kreal) :: start_time, delta_time_base, delta_time, current_time, next_time, total_time, end_time, delmax, delmin
22  real(kind=kreal) :: tmpmax, dltmp, tmpmax_myrank, remain_time
23  type(hecmwst_local_mesh) :: hecmesh
24  type(hecmwst_matrix) :: hecMAT
25  type(hecmwst_result_data) :: fstrRESULT
26  type(fstr_param) :: fstrPARAM
27  type(fstr_heat) :: fstrHEAT
28  type(hecmwst_local_mesh), pointer :: hecMESHmpc
29  type(hecmwst_matrix), pointer :: hecMATmpc
30  integer(kind=kint), parameter :: miniter = 4
31  logical :: is_end, outflag
32 
33  call hecmw_mpc_mat_init(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
34 
35  start_time = 0.0d0
36  do i = 1, istep - 1
37  start_time = start_time + fstrheat%STEP_EETIME(i)
38  enddo
39 
40  total_time = start_time + current_time
41 
42  delta_time_base = fstrheat%STEP_DLTIME(istep)
43  end_time = fstrheat%STEP_EETIME(istep)
44  delmin = fstrheat%STEP_DELMIN(istep)
45  delmax = fstrheat%STEP_DELMAX(istep)
46 
47  is_end = .false.
48  hecmat%NDOF = 1
49  hecmat%Iarray(98) = 1 !Assembly complete
50  hecmat%X = 0.0d0
51 
52  if(fstrheat%beta == -1.0d0)then
53  if(fstrheat%is_steady == 1)then
54  fstrheat%beta = 1.0d0
55  else
56  fstrheat%beta = 0.5d0
57  endif
58  endif
59 
60  if(fstrheat%is_steady /= 1 .and. total_step == 1) then
61  call heat_output_result(hecmesh, fstrheat, 0, total_time, .true.)
62  call heat_output_visual(hecmesh, fstrresult, fstrheat, 0, total_time, .true.)
63  endif
64 
65  !C-------------------- START TRANSIENT LOOP ------------------------
66  tr_loop: do
67 
68  if(end_time <= current_time + delta_time_base + delta_time_base*1.0d-6) then
69  delta_time_base = end_time - current_time
70  endif
71  if( 0.0d0 < delmin .and. fstrheat%timepoint_id > 0 ) then
72  remain_time = get_remain_to_next_timepoints(total_time, 0.0d0, fstrparam%timepoints(fstrheat%timepoint_id))
73  delta_time = dmin1(delta_time_base, remain_time)
74  else
75  delta_time = delta_time_base
76  endif
77  next_time = current_time + delta_time
78  total_time = start_time + next_time
79 
80  if( fstrheat%is_steady == 1 ) then
81  is_end = .true.
82  else
83  if( (end_time - next_time) / end_time < 1.d-12 ) is_end = .true.
84  endif
85 
86  if( 0.0d0 < delmin .and. fstrheat%timepoint_id > 0 ) then
87  outflag = is_end .or. is_at_timepoints(total_time, 0.0d0, fstrparam%timepoints(fstrheat%timepoint_id))
88  else
89  outflag = is_end
90  endif
91 
92  if( hecmesh%my_rank.eq.0 ) then
93  write(imsg,"(a,i8,a,1pe12.5,a,1pe12.5)") " ** Increment No. :", total_step, ", total time: ", &
94  & total_time, ", delta t: ", delta_time
95  write(*, "(a,i8,a,1pe12.5,a,1pe12.5)") " ** Increment No. :", total_step, ", total time: ", &
96  & total_time, ", delta t: ", delta_time
97  endif
98 
99  if(delta_time_base < delmin .and. (.not. is_end))then
100  if(hecmesh%my_rank == 0) write(imsg,*) ' !!! DELTA TIME EXCEEDED TOLERANCE OF TIME INCREMENT'
101  call hecmw_abort(hecmw_comm_get_comm())
102  endif
103 
104  call heat_solve_main(hecmesh, hecmat, hecmeshmpc, hecmatmpc, fstrparam, fstrheat, istep, iterall, total_time, delta_time)
105 
106  if(0.0d0 < delmin)then
107  tmpmax = 0.0d0
108  do i = 1, hecmesh%nn_internal
109  inod = fstrparam%global_local_id(1,i)
110  dltmp = fstrheat%TEMP0(i) - fstrheat%TEMP(i)
111  if(tmpmax < dabs(dltmp)) then
112  mnod = inod
113  tmpmax = dabs(dltmp)
114  endif
115  enddo
116  tmpmax_myrank = tmpmax
117  call hecmw_allreduce_r1(hecmesh, tmpmax, hecmw_max)
118 
119  if(tmpmax_myrank < tmpmax) mnod = -1
120  call hecmw_allreduce_i1(hecmesh, mnod, hecmw_max)
121 
122  if(delmax < tmpmax .or. fstrheat%is_iter_max_limit)then
123  if(hecmesh%my_rank == 0)then
124  write(*,*) ' *** EXCEEDED TOLERANCE OF VARIATION IN TEMPERATUTE.'
125  write(*,*) ' : NODE NUMBER = ', mnod, ' : DELTA TEMP = ', tmpmax
126  endif
127  delta_time_base = 0.5d0*delta_time_base
128  cycle tr_loop
129  endif
130 
131  if(iterall <= miniter) delta_time_base = delta_time_base*1.5d0
132  else
133  if(fstrheat%is_iter_max_limit) call hecmw_abort( hecmw_comm_get_comm() )
134  endif
135 
136  do i = 1, hecmesh%n_node
137  fstrheat%TEMP0(i) = fstrheat%TEMP(i)
138  enddo
139 
140  call heat_output_log(hecmesh, fstrparam, fstrheat, total_step, total_time)
141  call heat_output_result(hecmesh, fstrheat, total_step, total_time, outflag)
142  call heat_output_visual(hecmesh, fstrresult, fstrheat, total_step, total_time, outflag)
143  call heat_output_restart(hecmesh, fstrheat, istep, total_step, next_time, is_end)
144 
145  total_step = total_step + 1
146  current_time = next_time
147  if( is_end ) exit
148  enddo tr_loop
149  !C-------------------- END TRANSIENT LOOP ------------------------
150 
151  call hecmw_mpc_mat_finalize(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
152 
153  end subroutine heat_solve_tran
154 end module m_heat_solve_tran
m_heat_io
This module provides a function to control heat analysis.
Definition: heat_io.f90:6
m_heat_init
This module provides functions to initialize heat analysis.
Definition: heat_init.f90:6
m_heat_solve_main::heat_solve_main
subroutine heat_solve_main(hecMESH, hecMAT, hecMESHmpc, hecMATmpc, fstrPARAM, fstrHEAT, ISTEP, iterALL, next_time, delta_time)
Definition: heat_solve_main.f90:10
m_heat_solve_tran
This module provide a function to control nonsteady heat analysis.
Definition: heat_solve_TRAN.f90:6
m_heat_io::heat_output_restart
subroutine heat_output_restart(hecMESH, fstrHEAT, istep, tstep, current_time, outflag)
Definition: heat_io.f90:157
m_heat_solve_main
This module provides a function for stationary heat analysis.
Definition: heat_solve_main.f90:6
m_fstr::fstr_heat
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:425
m_fstr::fstr_param
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:154
m_heat_mat_ass_boundary
This module provides a subroutine for all boundary conditions needed in heat analysis.
Definition: heat_mat_ass_boundary.f90:7
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_heat_io::heat_output_visual
subroutine heat_output_visual(hecMESH, fstrRESULT, fstrHEAT, tstep, ctime, outflag)
Definition: heat_io.f90:120
m_heat_io::heat_output_log
subroutine heat_output_log(hecMESH, fstrPARAM, fstrHEAT, tstep, ctime)
Definition: heat_io.f90:42
m_heat_solve_tran::heat_solve_tran
subroutine heat_solve_tran(hecMESH, hecMAT, fstrRESULT, fstrPARAM, fstrHEAT, ISTEP, total_step, current_time)
Definition: heat_solve_TRAN.f90:10
m_heat_mat_ass_capacity
This module provides a subroutine to assemble heat capacity matrix.
Definition: heat_mat_ass_capacity.f90:6
m_heat_mat_ass_conductivity
Definition: heat_mat_ass_conductivity.f90:5
m_solve_lineq
Definition: solve_LINEQ.f90:6
m_heat_io::heat_output_result
subroutine heat_output_result(hecMESH, fstrHEAT, tstep, ctime, outflag)
Definition: heat_io.f90:90
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110