Neko  0.9.99
A portable framework for high-order spectral element flow simulations
simulation.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2021, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module simulation
35  use mpi_f08
36  use case, only : case_t
37  use num_types, only : rp, dp
39  use file, only : file_t
40  use logger, only : log_size, neko_log
41  use jobctrl, only : jobctrl_time_limit
42  use field, only : field_t
43  use profiler, only : profiler_start, profiler_stop, &
48  implicit none
49  private
50 
51  public :: neko_solve
52 
53 contains
54 
56  subroutine neko_solve(C)
57  type(case_t), target, intent(inout) :: c
58  real(kind=rp) :: t, cfl
59  real(kind=dp) :: start_time_org, start_time, end_time, tstep_start_time
60  character(len=LOG_SIZE) :: log_buf
61  integer :: tstep
62  character(len=:), allocatable :: restart_file
63  logical :: output_at_end, found
64  ! for variable_tsteping
65  real(kind=rp) :: cfl_avrg = 0.0_rp
66  type(time_step_controller_t) :: dt_controller
67  real(kind=rp) :: rho, mu, cp, lambda
68 
69  t = 0d0
70  tstep = 0
71  call neko_log%section('Starting simulation')
72  write(log_buf, '(A, E15.7,A,E15.7,A)') 'T : [', 0d0, ',', c%end_time, ')'
73  call neko_log%message(log_buf)
74  call dt_controller%init(c%params)
75  if (.not. dt_controller%if_variable_dt) then
76  write(log_buf, '(A, E15.7)') 'dt : ', c%dt
77  call neko_log%message(log_buf)
78  else
79  write(log_buf, '(A, E15.7)') 'CFL : ', dt_controller%set_cfl
80  call neko_log%message(log_buf)
81  end if
82 
83  call c%params%get('case.restart_file', restart_file, found)
84  if (found .and. len_trim(restart_file) .gt. 0) then
85  ! Restart the case
86  call simulation_restart(c, t)
87 
88  ! Restart the simulation components
89  call neko_simcomps%restart(t)
90  end if
91 
93  call neko_log%section('Postprocessing')
94  call c%output_controller%execute(t, tstep)
95 
96  call c%usr%user_init_modules(t, c%fluid%u, c%fluid%v, c%fluid%w,&
97  c%fluid%p, c%fluid%c_Xh, c%params)
98  call neko_log%end_section()
99  call neko_log%newline()
100 
101  call profiler_start
102  cfl = c%fluid%compute_cfl(c%dt)
103  start_time_org = mpi_wtime()
104 
105  do while (t .lt. c%end_time .and. (.not. jobctrl_time_limit()))
106  call profiler_start_region('Time-Step')
107  tstep = tstep + 1
108  start_time = mpi_wtime()
109  tstep_start_time = start_time
110  if (dt_controller%dt_last_change .eq. 0) then
111  cfl_avrg = cfl
112  end if
113  call dt_controller%set_dt(c%dt, cfl, cfl_avrg, tstep)
114  !calculate the cfl after the possibly varied dt
115  cfl = c%fluid%compute_cfl(c%dt)
116 
117  call neko_log%status(t, c%end_time)
118  write(log_buf, '(A,I6)') 'Time-step: ', tstep
119  call neko_log%message(log_buf)
120  call neko_log%begin()
121 
122  write(log_buf, '(A,E15.7,1x,A,E15.7)') 'CFL:', cfl, 'dt:', c%dt
123  call neko_log%message(log_buf)
124  call simulation_settime(t, c%dt, c%ext_bdf, c%tlag, c%dtlag, tstep)
125 
126  ! Run the preprocessing
127  call neko_log%section('Preprocessing')
128  call neko_simcomps%preprocess(t, tstep)
129  call neko_log%end_section()
130 
131  call neko_log%section('Fluid')
132  call c%fluid%step(t, tstep, c%dt, c%ext_bdf, dt_controller)
133  end_time = mpi_wtime()
134  write(log_buf, '(A,E15.7)') &
135  'Fluid step time (s): ', end_time-start_time
136  call neko_log%message(log_buf)
137  write(log_buf, '(A,E15.7)') &
138  'Total elapsed time (s):', end_time-start_time_org
139  call neko_log%end_section(log_buf)
140 
141  ! Scalar step
142  if (allocated(c%scalar)) then
143  start_time = mpi_wtime()
144  call neko_log%section('Scalar')
145  call c%scalar%step(t, tstep, c%dt, c%ext_bdf, dt_controller)
146  end_time = mpi_wtime()
147  write(log_buf, '(A,E15.7)') &
148  'Scalar step time: ', end_time-start_time
149  call neko_log%message(log_buf)
150  write(log_buf, '(A,E15.7)') &
151  'Total elapsed time (s):', end_time-start_time_org
152  call neko_log%end_section(log_buf)
153 
155  cp = c%scalar%cp
156  lambda = c%scalar%lambda
157  end if
158 
159  call neko_log%section('Postprocessing')
160  ! Execute all simulation components
161  call neko_simcomps%compute(t, tstep)
162 
163  call c%output_controller%execute(t, tstep)
164 
166  rho = c%fluid%rho
167  mu = c%fluid%mu
168 
169  ! Update material properties
170  call c%usr%material_properties(t, tstep, rho, mu, cp, lambda, c%params)
171 
173  c%fluid%rho = rho
174  c%fluid%mu = mu
175  call c%fluid%update_material_properties()
176 
177  if (allocated(c%scalar)) then
178  c%scalar%cp = cp
179  c%scalar%lambda = lambda
180  call c%scalar%update_material_properties()
181  end if
182 
183  call c%usr%user_check(t, tstep, c%fluid%u, c%fluid%v, c%fluid%w, &
184  c%fluid%p, c%fluid%c_Xh, c%params)
185 
186  call neko_log%end_section()
187  end_time = mpi_wtime()
188  call neko_log%section('Step summary')
189  write(log_buf, '(A,I8,A,E15.7)') &
190  'Total time for step ', tstep, ' (s): ', end_time-tstep_start_time
191  call neko_log%message(log_buf)
192  write(log_buf, '(A,E15.7)') &
193  'Total elapsed time (s): ', end_time-start_time_org
194  call neko_log%message(log_buf)
195  call neko_log%end_section()
196  call neko_log%end()
198  end do
199  call profiler_stop
200 
201  call json_get_or_default(c%params, 'case.output_at_end',&
202  output_at_end, .true.)
203  call c%output_controller%execute(t, tstep, output_at_end)
204 
205  if (.not. (output_at_end) .and. t .lt. c%end_time) then
206  call simulation_joblimit_chkp(c, t)
207  end if
208 
209  call c%usr%user_finalize_modules(t, c%params)
210 
211  call neko_log%end_section('Normal end.')
212 
213  end subroutine neko_solve
214 
215  subroutine simulation_settime(t, dt, ext_bdf, tlag, dtlag, step)
216  real(kind=rp), intent(inout) :: t
217  real(kind=rp), intent(in) :: dt
218  type(time_scheme_controller_t), intent(inout) :: ext_bdf
219  real(kind=rp), dimension(10) :: tlag
220  real(kind=rp), dimension(10) :: dtlag
221  integer, intent(in) :: step
222  integer :: i
223 
224 
225  do i = 10, 2, -1
226  tlag(i) = tlag(i-1)
227  dtlag(i) = dtlag(i-1)
228  end do
229 
230  dtlag(1) = dt
231  tlag(1) = t
232  if (ext_bdf%ndiff .eq. 0) then
233  dtlag(2) = dt
234  tlag(2) = t
235  end if
236 
237  t = t + dt
238 
239  call ext_bdf%set_coeffs(dtlag)
240 
241  end subroutine simulation_settime
242 
244  subroutine simulation_restart(C, t)
245  implicit none
246  type(case_t), intent(inout) :: C
247  real(kind=rp), intent(inout) :: t
248  integer :: i
249  type(file_t) :: chkpf, previous_meshf
250  character(len=LOG_SIZE) :: log_buf
251  character(len=:), allocatable :: restart_file
252  character(len=:), allocatable :: restart_mesh_file
253  real(kind=rp) :: tol
254  logical :: found, check_cont
255 
256  call c%params%get('case.restart_file', restart_file, found)
257  call c%params%get('case.restart_mesh_file', restart_mesh_file,&
258  found)
259 
260  if (found) then
261  previous_meshf = file_t(trim(restart_mesh_file))
262  call previous_meshf%read(c%fluid%chkp%previous_mesh)
263  end if
264 
265  call c%params%get('case.mesh2mesh_tolerance', tol,&
266  found)
267 
268  if (found) c%fluid%chkp%mesh2mesh_tol = tol
269 
270  chkpf = file_t(trim(restart_file))
271  call chkpf%read(c%fluid%chkp)
272  c%dtlag = c%fluid%chkp%dtlag
273  c%tlag = c%fluid%chkp%tlag
274 
275  !Free the previous mesh, dont need it anymore
276  do i = 1, size(c%dtlag)
277  call c%ext_bdf%set_coeffs(c%dtlag)
278  end do
279 
280  call c%fluid%restart(c%dtlag, c%tlag)
281  call c%fluid%chkp%previous_mesh%free()
282  if (allocated(c%scalar)) &
283  call c%scalar%restart( c%dtlag, c%tlag)
284 
285  t = c%fluid%chkp%restart_time()
286  call neko_log%section('Restarting from checkpoint')
287  write(log_buf, '(A,A)') 'File : ', &
288  trim(restart_file)
289  call neko_log%message(log_buf)
290  write(log_buf, '(A,E15.7)') 'Time : ', t
291  call neko_log%message(log_buf)
292  call neko_log%end_section()
293 
294  call c%output_controller%set_counter(t)
295  end subroutine simulation_restart
296 
298  subroutine simulation_joblimit_chkp(C, t)
299  type(case_t), intent(inout) :: C
300  real(kind=rp), intent(inout) :: t
301  type(file_t) :: chkpf
302  character(len=:), allocatable :: chkp_format
303  character(len=LOG_SIZE) :: log_buf
304  character(len=10) :: format_str
305  logical :: found
306 
307  call c%params%get('case.checkpoint_format', chkp_format, found)
308  call c%fluid%chkp%sync_host()
309  format_str = '.chkp'
310  if (found) then
311  if (chkp_format .eq. 'hdf5') then
312  format_str = '.h5'
313  end if
314  end if
315  chkpf = file_t(c%output_directory // 'joblimit'//trim(format_str))
316  call chkpf%write(c%fluid%chkp, t)
317  write(log_buf, '(A)') '! saving checkpoint >>>'
318  call neko_log%message(log_buf)
319 
320  end subroutine simulation_joblimit_chkp
321 
322 end module simulation
323 
324 
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:54
Defines a simulation case.
Definition: case.f90:34
Defines a field.
Definition: field.f90:34
Module for file I/O operations.
Definition: file.f90:34
Job control.
Definition: jobctrl.f90:34
logical function, public jobctrl_time_limit()
Check if the job's time limit has been reached.
Definition: jobctrl.f90:107
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Logging routines.
Definition: log.f90:34
type(log_t), public neko_log
Global log stream.
Definition: log.f90:65
integer, parameter, public log_size
Definition: log.f90:42
integer, parameter, public dp
Definition: num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Profiling interface.
Definition: profiler.F90:34
subroutine, public profiler_start
Start profiling.
Definition: profiler.F90:52
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition: profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition: profiler.F90:109
subroutine, public profiler_stop
Stop profiling.
Definition: profiler.F90:65
Contains the simcomp_executor_t type.
type(simcomp_executor_t), target, public neko_simcomps
Global variable for the simulation component driver.
Simulation driver.
Definition: simulation.f90:34
subroutine simulation_settime(t, dt, ext_bdf, tlag, dtlag, step)
Definition: simulation.f90:216
subroutine, public neko_solve(C)
Main driver to solve a case C.
Definition: simulation.f90:57
subroutine simulation_joblimit_chkp(C, t)
Write a checkpoint at joblimit.
Definition: simulation.f90:299
subroutine simulation_restart(C, t)
Restart a case C from a given checkpoint.
Definition: simulation.f90:245
Compound scheme for the advection and diffusion operators in a transport equation.
Implements type time_step_controller.
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition: file.f90:54
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
Provides a tool to set time step dt.