Neko  0.8.1
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 case
36  use gather_scatter
38  use file
39  use math
40  use logger
41  use device
42  use device_math
43  use jobctrl
44  use field, only : field_t
45  use profiler
46  use math, only : col2
50  implicit none
51  private
52 
53  public :: neko_solve
54 
55 contains
56 
58  subroutine neko_solve(C)
59  type(case_t), target, intent(inout) :: c
60  real(kind=rp) :: t, cfl
61  real(kind=dp) :: start_time_org, start_time, end_time
62  character(len=LOG_SIZE) :: log_buf
63  integer :: tstep, i
64  character(len=:), allocatable :: restart_file
65  logical :: output_at_end, found
66  ! for variable_tsteping
67  real(kind=rp) :: cfl_avrg = 0.0_rp
68  type(time_step_controller_t) :: dt_controller
69 
70  t = 0d0
71  tstep = 0
72  call neko_log%section('Starting simulation')
73  write(log_buf,'(A, E15.7,A,E15.7,A)') 'T : [', 0d0,',',c%end_time,')'
74  call neko_log%message(log_buf)
75  call dt_controller%init(c%params)
76  if (.not. dt_controller%if_variable_dt) then
77  write(log_buf,'(A, E15.7)') 'dt : ', c%dt
78  call neko_log%message(log_buf)
79  else
80  write(log_buf,'(A, E15.7)') 'CFL : ', dt_controller%set_cfl
81  call neko_log%message(log_buf)
82  end if
83 
84  call c%params%get('case.restart_file', restart_file, found)
85  if (found .and. len_trim(restart_file) .gt. 0) then
86  ! Restart the case
87  call simulation_restart(c, t)
88 
89  ! Restart the simulation components
90  call neko_simcomps%restart(t)
91  end if
92 
94  call neko_log%section('Postprocessing')
95  call c%q%eval(t, c%dt, tstep)
96  call c%s%sample(t, tstep)
97 
98  call c%usr%user_init_modules(t, c%fluid%u, c%fluid%v, c%fluid%w,&
99  c%fluid%p, c%fluid%c_Xh, c%params)
100  call neko_log%end_section()
101  call neko_log%newline()
102 
103  call profiler_start
104  cfl = c%fluid%compute_cfl(c%dt)
105  start_time_org = mpi_wtime()
106 
107  do while (t .lt. c%end_time .and. (.not. jobctrl_time_limit()))
108  call profiler_start_region('Time-Step')
109  tstep = tstep + 1
110  start_time = mpi_wtime()
111  if (dt_controller%dt_last_change .eq. 0) then
112  cfl_avrg = cfl
113  end if
114  call dt_controller%set_dt(c%dt, cfl, cfl_avrg, tstep)
115  !calculate the cfl after the possibly varied dt
116  cfl = c%fluid%compute_cfl(c%dt)
117 
118  call neko_log%status(t, c%end_time)
119  write(log_buf, '(A,I6)') 'Time-step: ', tstep
120  call neko_log%message(log_buf)
121  call neko_log%begin()
122 
123  write(log_buf, '(A,E15.7,1x,A,E15.7)') 'CFL:', cfl, 'dt:', c%dt
124  call neko_log%message(log_buf)
125  call simulation_settime(t, c%dt, c%ext_bdf, c%tlag, c%dtlag, tstep)
126 
127  call neko_log%section('Fluid')
128  call c%fluid%step(t, tstep, c%dt, c%ext_bdf, dt_controller)
129  end_time = mpi_wtime()
130  write(log_buf, '(A,E15.7,A,E15.7)') &
131  'Elapsed time (s):', end_time-start_time_org, ' Step time:', &
132  end_time-start_time
133  call neko_log%end_section(log_buf)
134 
135  ! Scalar step
136  if (allocated(c%scalar)) then
137  start_time = mpi_wtime()
138  call neko_log%section('Scalar')
139  call c%scalar%step(t, tstep, c%dt, c%ext_bdf, dt_controller)
140  end_time = mpi_wtime()
141  write(log_buf, '(A,E15.7,A,E15.7)') &
142  'Elapsed time (s):', end_time-start_time_org, ' Step time:', &
143  end_time-start_time
144  call neko_log%end_section(log_buf)
145  end if
146 
147  call neko_log%section('Postprocessing')
148  ! Execute all simulation components
149  call neko_simcomps%compute(t, tstep)
150 
151  call c%q%eval(t, c%dt, tstep)
152  call c%s%sample(t, tstep)
153 
154  ! Update material properties
155  call c%usr%material_properties(t, tstep, c%material_properties%rho,&
156  c%material_properties%mu, &
157  c%material_properties%cp, &
158  c%material_properties%lambda, &
159  c%params)
160 
161  call c%usr%user_check(t, tstep,&
162  c%fluid%u, c%fluid%v, c%fluid%w, c%fluid%p, c%fluid%c_Xh, c%params)
163 
164  call neko_log%end_section()
165 
166  call neko_log%end()
168  end do
169  call profiler_stop
170 
171  call json_get_or_default(c%params, 'case.output_at_end',&
172  output_at_end, .true.)
173  call c%s%sample(t, tstep, output_at_end)
174 
175  if (.not. (output_at_end) .and. t .lt. c%end_time) then
176  call simulation_joblimit_chkp(c, t)
177  end if
178 
179  call c%usr%user_finalize_modules(t, c%params)
180 
181  call neko_log%end_section('Normal end.')
182 
183  end subroutine neko_solve
184 
185  subroutine simulation_settime(t, dt, ext_bdf, tlag, dtlag, step)
186  real(kind=rp), intent(inout) :: t
187  real(kind=rp), intent(in) :: dt
188  type(time_scheme_controller_t), intent(inout) :: ext_bdf
189  real(kind=rp), dimension(10) :: tlag
190  real(kind=rp), dimension(10) :: dtlag
191  integer, intent(in) :: step
192  integer :: i
193 
194 
195  do i = 10, 2, -1
196  tlag(i) = tlag(i-1)
197  dtlag(i) = dtlag(i-1)
198  end do
199 
200  dtlag(1) = dt
201  tlag(1) = t
202  if (ext_bdf%ndiff .eq. 0) then
203  dtlag(2) = dt
204  tlag(2) = t
205  end if
206 
207  t = t + dt
208 
209  call ext_bdf%set_coeffs(dtlag)
210 
211  end subroutine simulation_settime
212 
214  subroutine simulation_restart(C, t)
215  implicit none
216  type(case_t), intent(inout) :: C
217  real(kind=rp), intent(inout) :: t
218  integer :: i
219  type(file_t) :: chkpf, previous_meshf
220  character(len=LOG_SIZE) :: log_buf
221  character(len=:), allocatable :: restart_file
222  character(len=:), allocatable :: restart_mesh_file
223  real(kind=rp) :: tol
224  logical :: found
225 
226  call c%params%get('case.restart_file', restart_file, found)
227  call c%params%get('case.restart_mesh_file', restart_mesh_file,&
228  found)
229 
230  if (found) then
231  previous_meshf = file_t(trim(restart_mesh_file))
232  call previous_meshf%read(c%fluid%chkp%previous_mesh)
233  end if
234 
235  call c%params%get('case.mesh2mesh_tolerance', tol,&
236  found)
237 
238  if (found) c%fluid%chkp%mesh2mesh_tol = tol
239 
240  chkpf = file_t(trim(restart_file))
241  call chkpf%read(c%fluid%chkp)
242  c%dtlag = c%fluid%chkp%dtlag
243  c%tlag = c%fluid%chkp%tlag
244 
245  !Free the previous mesh, dont need it anymore
246  call c%fluid%chkp%previous_mesh%free()
247  do i = 1, size(c%dtlag)
248  call c%ext_bdf%set_coeffs(c%dtlag)
249  end do
250 
251  call c%fluid%restart(c%dtlag, c%tlag)
252  if (allocated(c%scalar)) call c%scalar%restart( c%dtlag, c%tlag)
253 
254  t = c%fluid%chkp%restart_time()
255  call neko_log%section('Restarting from checkpoint')
256  write(log_buf,'(A,A)') 'File : ', &
257  trim(restart_file)
258  call neko_log%message(log_buf)
259  write(log_buf,'(A,E15.7)') 'Time : ', t
260  call neko_log%message(log_buf)
261  call neko_log%end_section()
262 
263  call c%s%set_counter(t)
264  end subroutine simulation_restart
265 
267  subroutine simulation_joblimit_chkp(C, t)
268  type(case_t), intent(inout) :: C
269  real(kind=rp), intent(inout) :: t
270  type(file_t) :: chkpf
271  character(len=LOG_SIZE) :: log_buf
272 
273  call c%fluid%chkp%sync_host()
274  chkpf = file_t('joblimit.chkp')
275  call chkpf%write(c%fluid%chkp, t)
276  write(log_buf, '(A)') '! saving checkpoint >>>'
277  call neko_log%message(log_buf)
278 
279  end subroutine simulation_joblimit_chkp
280 
281 end module simulation
282 
283 
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:53
Defines a simulation case.
Definition: case.f90:34
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Defines a field.
Definition: field.f90:34
Module for file I/O operations.
Definition: file.f90:34
Gather-scatter.
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:106
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:61
Definition: math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:645
Profiling interface.
Definition: profiler.F90:34
subroutine, public profiler_start
Start profiling.
Definition: profiler.F90:50
subroutine, public profiler_end_region
End the most recently started profiler region.
Definition: profiler.F90:95
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition: profiler.F90:76
subroutine, public profiler_stop
Stop profiling.
Definition: profiler.F90:63
Contains the simcomp_executor_t type.
type(simcomp_executor_t), 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:186
subroutine, public neko_solve(C)
Main driver to solve a case C.
Definition: simulation.f90:59
subroutine simulation_joblimit_chkp(C, t)
Write a checkpoint at joblimit.
Definition: simulation.f90:268
subroutine simulation_restart(C, t)
Restart a case C from a given checkpoint.
Definition: simulation.f90:215
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:53
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
Provides a tool to set time step dt.