FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
readtemp.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 !-------------------------------------------------------------------------------
5 module mreadtemp
6 
7  private
9 
10 contains
11 
13  subroutine read_temperature_result(hecMESH, nstep, sstep, rtype, interval, factor, ctime, temp, temp_bak)
14  use m_fstr
15  type(hecmwst_local_mesh), intent(in) :: hecmesh
16  integer(kind=kint), intent(in) :: nstep
17  integer(kind=kint), intent(in) :: sstep
18  integer(kind=kint), intent(in) :: rtype
19  integer(kind=kint), intent(in) :: interval
20  real(kind=kreal), intent(in) :: factor
21  real(kind=kreal), intent(in) :: ctime
22  real(kind=kreal), pointer :: temp(:)
23  real(kind=kreal), pointer :: temp_bak(:)
24  if( rtype == 1 ) then
25  call read_temperature_result_by_step(hecmesh, nstep, sstep, interval, factor, temp, temp_bak)
26  elseif( rtype == 2 ) then
27  call read_temperature_result_by_time(hecmesh, nstep, sstep, ctime, temp)
28  endif
29  write(idbg,*) ' Read temperature from result file : OK'
30  end subroutine read_temperature_result
31 
32  subroutine read_temperature_result_by_step(hecMESH, nstep, sstep, interval, factor, temp, temp_bak)
34  use m_fstr
35  type(hecmwst_local_mesh), intent(in) :: hecmesh
36  integer(kind=kint), intent(in) :: nstep
37  integer(kind=kint), intent(in) :: sstep
38  integer(kind=kint), intent(in) :: interval
39  real(kind=kreal), intent(in) :: factor
40  real(kind=kreal), pointer :: temp(:)
41  real(kind=kreal), pointer :: temp_bak(:)
42  type(hecmwst_result_data) :: result
43  character(len=HECMW_NAME_LEN) :: name_id
44  integer(kind=kint) :: i, k0, kt, nstep_active, tstep
45  real(kind=kreal) :: w
46 
47  real(kind=kreal), pointer :: temp0(:,:)
48 
49  allocate(temp0(2,size(temp)))
50 
51  nstep_active = (nstep-sstep)/interval + 1
52  kt = floor(factor*dble(nstep_active)-1.d-10)
53  w = factor*dble(nstep_active)-dble(kt)
54 
55  do k0=1,2
56  tstep = sstep+(kt+k0-2)*interval
57  if( tstep == 0 ) then
58  do i=1, hecmesh%n_node
59  temp0(k0,i) = temp_bak(i)
60  enddo
61  cycle
62  end if
63 
64  call hecmw_nullify_result_data( result )
65  name_id = 'fstrTEMP'
66  call hecmw_result_read_by_name(hecmesh, name_id, tstep, result)
67 
68  if (result%nn_component /= 1 .or. result%nn_dof(1) /= 1) then
69  write(*,*) ' Read temperature result failed; not heat analysis result'
70  endif
71 
72  do i=1, hecmesh%n_node
73  temp0(k0,i) = result%node_val_item(i)
74  enddo
75 
76  call hecmw_result_free(result)
77  enddo
78 
79  do i=1, hecmesh%n_node
80  temp(i) = (1.d0-w)*temp0(1,i)+w*temp0(2,i)
81  enddo
82 
83  deallocate(temp0)
84  end subroutine read_temperature_result_by_step
85 
87  subroutine read_temperature_result_by_time(hecMESH, nstep, sstep, ctime, temp)
89  use m_fstr
90  type(hecmwst_local_mesh), intent(in) :: hecmesh
91  integer(kind=kint), intent(in) :: nstep
92  integer(kind=kint), intent(in) :: sstep
93  real(kind=kreal), intent(in) :: ctime
94  real(kind=kreal), pointer :: temp(:)
95 
96  logical, save :: flg_first_call = .true.
97  logical, save :: flg_constant = .false.
98  logical, save :: flg_finished = .false.
99  integer(kind=kint), save :: tstep1, tstep2
100  real(kind=kreal), save :: ttime1, ttime2
101 
102  integer(kind=kint) :: ierr, i
103  real(kind=kreal) :: w
104  logical :: flg_read1, flg_read2
105  real(kind=kreal), allocatable :: temp0(:,:)
106 
107  allocate(temp0(2,size(temp)))
108  flg_read1 = .false.
109  flg_read2 = .false.
110  ! read first 2 results at first call
111  if (flg_first_call) then
112  call read_next_result(hecmesh, nstep, sstep, tstep1, ttime1, temp0(1,:), ierr)
113  ! if no more result; error: no result file
114  if (ierr /= 0) then
115  write(*,*) ' Read temperature result failed; no result file'
116  call hecmw_abort(hecmw_comm_get_comm())
117  endif
118  flg_read1 = .true.
119  call read_next_result(hecmesh, nstep, tstep1+1, tstep2, ttime2, temp0(2,:), ierr)
120  ! if no more result (only 1 result exists); assume temperature is constant
121  if (ierr /= 0) flg_constant = .true.
122  flg_read2 = .true.
123  flg_first_call = .false.
124  endif
125  ! when temperature is assumed constant
126  if (flg_constant) then
127  if (.not. flg_read1) call read_result_at(hecmesh, tstep1, ttime1, temp0(1,:))
128  do i=1, hecmesh%n_node
129  temp(i) = temp0(1,i)
130  enddo
131  return
132  endif
133  ! when ctime is out-of-range (left side)
134  if (ctime <= ttime1) then
135  if (.not. flg_read1) call read_result_at(hecmesh, tstep1, ttime1, temp0(1,:))
136  do i=1, hecmesh%n_node
137  temp(i) = temp0(1,i)
138  enddo
139  return
140  endif
141  ! proceed until ctime is within [ttime1, ttime2]
142  do while (ttime2 < ctime)
143  tstep1 = tstep2
144  ttime1 = ttime2
145  if (flg_read2) then
146  do i=1, hecmesh%n_node
147  temp0(1,i) = temp0(2,i)
148  enddo
149  flg_read1 = .true.
150  endif
151  flg_read2 = .false.
152  call read_next_result(hecmesh, nstep, tstep1+1, tstep2, ttime2, temp0(2,:), ierr)
153  ! if no more result; finished reading (ctime is out-of-range)
154  if (ierr /= 0) then
155  flg_finished = .true.
156  exit
157  endif
158  flg_read2 = .true.
159  enddo
160  ! when ctime is out-of-range (right side)
161  if (flg_finished) then
162  if (.not. flg_read1) call read_result_at(hecmesh, tstep1, ttime1, temp0(1,:))
163  do i=1, hecmesh%n_node
164  temp(i) = temp0(1,i)
165  enddo
166  return
167  endif
168  ! now, ttime1 <= ctime <= ttime2
169  if (.not. flg_read1) call read_result_at(hecmesh, tstep1, ttime1, temp0(1,:))
170  if (.not. flg_read2) call read_result_at(hecmesh, tstep2, ttime2, temp0(2,:))
171  w = (ctime - ttime1) / (ttime2 - ttime1)
172  do i=1, hecmesh%n_node
173  temp(i) = (1-w)*temp0(1,i) + w*temp0(2,i)
174  enddo
175  !write(IDBG,'(a,f10.2,a,i0,a,f10.2,a,i0,a,f10.2,a)') &
176  ! ' Interpolated temp at',ctime,' from step ',tstep1,' (time',ttime1,') and step ',tstep2,' (time',ttime2,')'
177  return
178  end subroutine read_temperature_result_by_time
179 
180  subroutine read_result_at(hecMESH, tstep, ttime, temp)
181  use m_fstr
182  use hecmw_result
183  implicit none
184  type(hecmwst_local_mesh), intent(in) :: hecmesh
185  integer(kind=kint), intent(in) :: tstep
186  real(kind=kreal), intent(out) :: ttime
187  real(kind=kreal), intent(out) :: temp(:)
188  type(hecmwst_result_data) :: result
189  character(len=HECMW_NAME_LEN) :: name_id
190  integer(kind=kint) :: i
191  call hecmw_nullify_result_data( result )
192  name_id = 'fstrTEMP'
193  call hecmw_result_read_by_name(hecmesh, name_id, tstep, result)
194  if (result%nn_component /= 1 .or. result%nn_dof(1) /= 1) then
195  write(*,*) ' Read temperature result failed; not heat analysis result'
196  endif
197  ttime = result%global_val_item(1)
198  do i=1, hecmesh%n_node
199  temp(i) = result%node_val_item(i)
200  enddo
201  call hecmw_result_free(result)
202  end subroutine read_result_at
203 
204  subroutine read_next_result(hecMESH, nstep, sstep, tstep, ttime, temp, ierr)
205  use m_fstr
206  use hecmw_result
207  implicit none
208  type(hecmwst_local_mesh), intent(in) :: hecMESH
209  integer(kind=kint), intent(in) :: nstep
210  integer(kind=kint), intent(in) :: sstep
211  integer(kind=kint), intent(out) :: tstep
212  real(kind=kreal), intent(out) :: ttime
213  real(kind=kreal), intent(out) :: temp(:)
214  integer(kind=kint), intent(out) :: ierr
215  character(len=HECMW_NAME_LEN) :: name_ID
216  integer(kind=kint) :: istep
217  name_id = 'fstrTEMP'
218  tstep = -1
219  do istep = sstep, nstep
220  call hecmw_result_checkfile_by_name(name_id, istep, ierr)
221  if (ierr == 0) then
222  tstep = istep
223  exit
224  endif
225  enddo
226  if (tstep < 0) then
227  ierr = 1
228  return
229  endif
230  call read_result_at(hecmesh, tstep, ttime, temp)
231  end subroutine read_next_result
232 
233 end module mreadtemp
hecmw_result::hecmw_result_checkfile_by_name
subroutine, public hecmw_result_checkfile_by_name(name_ID, i_step, ierr)
Definition: hecmw_result_f.f90:243
hecmw_result::hecmw_result_read_by_name
subroutine, public hecmw_result_read_by_name(hecMESH, name_ID, i_step, result)
Definition: hecmw_result_f.f90:252
mreadtemp::read_temperature_result_by_step
subroutine read_temperature_result_by_step(hecMESH, nstep, sstep, interval, factor, temp, temp_bak)
Definition: readtemp.f90:33
m_fstr::idbg
integer(kind=kint), parameter idbg
Definition: m_fstr.f90:111
mreadtemp::read_result_at
subroutine read_result_at(hecMESH, tstep, ttime, temp)
Definition: readtemp.f90:181
mreadtemp::read_next_result
subroutine read_next_result(hecMESH, nstep, sstep, tstep, ttime, temp, ierr)
Definition: readtemp.f90:205
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
mreadtemp::read_temperature_result_by_time
subroutine read_temperature_result_by_time(hecMESH, nstep, sstep, ctime, temp)
Read in temperature distribution from external file.
Definition: readtemp.f90:88
mreadtemp::read_temperature_result
subroutine, public read_temperature_result(hecMESH, nstep, sstep, rtype, interval, factor, ctime, temp, temp_bak)
Read in temperature distribution from external file.
Definition: readtemp.f90:14
mreadtemp
Definition: readtemp.f90:5
hecmw_result::hecmw_nullify_result_data
subroutine, public hecmw_nullify_result_data(P)
Definition: hecmw_result_f.f90:64
hecmw_result
I/O and Utility.
Definition: hecmw_result_f.f90:7
hecmw_result::hecmw_result_free
subroutine, public hecmw_result_free(result_data)
Definition: hecmw_result_f.f90:423