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
62 character(len=:),
allocatable :: restart_file
63 logical :: output_at_end, found
65 real(kind=
rp) :: cfl_avrg = 0.0_rp
67 real(kind=
rp) :: rho, mu, cp, lambda
71 call neko_log%section(
'Starting simulation')
72 write(log_buf,
'(A, E15.7,A,E15.7,A)')
'T : [', 0d0,
',', c%end_time,
')'
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
79 write(log_buf,
'(A, E15.7)')
'CFL : ', dt_controller%set_cfl
83 call c%params%get(
'case.restart_file', restart_file, found)
84 if (found .and. len_trim(restart_file) .gt. 0)
then
93 call neko_log%section(
'Postprocessing')
94 call c%output_controller%execute(t, tstep)
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)
102 cfl = c%fluid%compute_cfl(c%dt)
103 start_time_org = mpi_wtime()
108 start_time = mpi_wtime()
109 tstep_start_time = start_time
110 if (dt_controller%dt_last_change .eq. 0)
then
113 call dt_controller%set_dt(c%dt, cfl, cfl_avrg, tstep)
115 cfl = c%fluid%compute_cfl(c%dt)
118 write(log_buf,
'(A,I6)')
'Time-step: ', tstep
122 write(log_buf,
'(A,E15.7,1x,A,E15.7)')
'CFL:', cfl,
'dt:', c%dt
127 call neko_log%section(
'Preprocessing')
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
137 write(log_buf,
'(A,E15.7)') &
138 'Total elapsed time (s):', end_time-start_time_org
142 if (
allocated(c%scalar))
then
143 start_time = mpi_wtime()
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
150 write(log_buf,
'(A,E15.7)') &
151 'Total elapsed time (s):', end_time-start_time_org
156 lambda = c%scalar%lambda
159 call neko_log%section(
'Postprocessing')
163 call c%output_controller%execute(t, tstep)
170 call c%usr%material_properties(t, tstep, rho, mu, cp, lambda, c%params)
175 call c%fluid%update_material_properties()
177 if (
allocated(c%scalar))
then
179 c%scalar%lambda = lambda
180 call c%scalar%update_material_properties()
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)
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
192 write(log_buf,
'(A,E15.7)') &
193 'Total elapsed time (s): ', end_time-start_time_org
202 output_at_end, .true.)
203 call c%output_controller%execute(t, tstep, output_at_end)
205 if (.not. (output_at_end) .and. t .lt. c%end_time)
then
209 call c%usr%user_finalize_modules(t, c%params)
211 call neko_log%end_section(
'Normal end.')
216 real(kind=
rp),
intent(inout) :: t
217 real(kind=
rp),
intent(in) :: dt
219 real(kind=
rp),
dimension(10) :: tlag
220 real(kind=
rp),
dimension(10) :: dtlag
221 integer,
intent(in) :: step
227 dtlag(i) = dtlag(i-1)
232 if (ext_bdf%ndiff .eq. 0)
then
239 call ext_bdf%set_coeffs(dtlag)
246 type(
case_t),
intent(inout) :: C
247 real(kind=
rp),
intent(inout) :: t
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
254 logical :: found, check_cont
256 call c%params%get(
'case.restart_file', restart_file, found)
257 call c%params%get(
'case.restart_mesh_file', restart_mesh_file,&
261 previous_meshf =
file_t(trim(restart_mesh_file))
262 call previous_meshf%read(c%fluid%chkp%previous_mesh)
265 call c%params%get(
'case.mesh2mesh_tolerance', tol,&
268 if (found) c%fluid%chkp%mesh2mesh_tol = tol
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
276 do i = 1,
size(c%dtlag)
277 call c%ext_bdf%set_coeffs(c%dtlag)
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)
285 t = c%fluid%chkp%restart_time()
286 call neko_log%section(
'Restarting from checkpoint')
287 write(log_buf,
'(A,A)')
'File : ', &
290 write(log_buf,
'(A,E15.7)')
'Time : ', t
294 call c%output_controller%set_counter(t)
299 type(
case_t),
intent(inout) :: C
300 real(kind=
rp),
intent(inout) :: t
302 character(len=:),
allocatable :: chkp_format
303 character(len=LOG_SIZE) :: log_buf
304 character(len=10) :: format_str
307 call c%params%get(
'case.checkpoint_format', chkp_format, found)
308 call c%fluid%chkp%sync_host()
311 if (chkp_format .eq.
'hdf5')
then
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 >>>'
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Defines a simulation case.
Module for file I/O operations.
logical function, public jobctrl_time_limit()
Check if the job's time limit has been reached.
Utilities for retrieving parameters from the case files.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter, public dp
integer, parameter, public rp
Global precision used in computations.
subroutine, public profiler_start
Start profiling.
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
subroutine, public profiler_stop
Stop profiling.
Contains the simcomp_executor_t type.
type(simcomp_executor_t), target, public neko_simcomps
Global variable for the simulation component driver.
subroutine simulation_settime(t, dt, ext_bdf, tlag, dtlag, step)
subroutine, public neko_solve(C)
Main driver to solve a case C.
subroutine simulation_joblimit_chkp(C, t)
Write a checkpoint at joblimit.
subroutine simulation_restart(C, t)
Restart a case C from a given checkpoint.
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...
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
Provides a tool to set time step dt.