42submodule(simulation_checkpoint) checkpoint_linear
 
   43  use simulation, 
only: simulation_step, simulation_restart
 
   44  use file, 
only: file_t, file_free
 
   45  use time_step_controller, 
only: time_step_controller_t
 
   52  module subroutine checkpoint_save_linear(this, neko_case)
 
   53    class(simulation_checkpoint_t), 
intent(inout) :: this
 
   54    class(case_t), 
intent(inout) :: neko_case
 
   58    save_disc = modulo(neko_case%time%tstep, this%n_saves_memory) .eq. 0
 
   61    if (save_disc .or. neko_case%time%tstep .le. this%first_valid_timestep) 
then 
   63       call this%chkp_output%set_counter(neko_case%time%tstep)
 
   64       call this%chkp_output%sample(neko_case%time%t)
 
   65       this%n_saves_disc = this%n_saves_disc + 1
 
   67  end subroutine checkpoint_save_linear
 
   73  module subroutine checkpoint_restore_linear(this, neko_case, tstep)
 
   74    class(simulation_checkpoint_t), 
intent(inout) :: this
 
   75    class(case_t), 
target, 
intent(inout) :: neko_case
 
   76    integer, 
intent(in) :: tstep
 
   77    type(time_step_controller_t) :: dt_controller
 
   78    real(kind=dp) :: loop_start
 
   79    integer :: j, k, previous_save, next_save
 
   81    type(field_t), 
pointer :: u, v, w, p, s
 
   83    loop_start = mpi_wtime()
 
   85    u => neko_case%fluid%u
 
   86    v => neko_case%fluid%v
 
   87    w => neko_case%fluid%w
 
   88    p => neko_case%fluid%p
 
   92    previous_save = tstep - modulo(tstep, this%n_saves_memory)
 
   93    next_save = previous_save + this%n_saves_memory
 
   97    if (tstep .lt. this%first_valid_timestep) 
then 
   99       next_save = previous_save + 1
 
  100    else if (previous_save .lt. this%first_valid_timestep) 
then 
  101       previous_save = this%first_valid_timestep
 
  105    if (this%loaded_checkpoint .ne. previous_save) 
then 
  108       call this%chkp_output%set_counter(previous_save)
 
  109       call this%chkp_output%file_%read(neko_case%chkp)
 
  110       call simulation_restart(neko_case, neko_case%chkp)
 
  113       call dt_controller%init(neko_case%params)
 
  114       neko_case%time%tstep = previous_save
 
  115       this%loaded_checkpoint = neko_case%time%tstep
 
  118       do k = previous_save, min(next_save - 1, this%n_timesteps)
 
  121          if (k .ne. previous_save) 
then 
  122             if (neko_case%time%t .ge. neko_case%time%end_time) 
exit 
  123             call simulation_step(neko_case, dt_controller, loop_start)
 
  126          local_idx = modulo(k, this%n_saves_memory) + 1
 
  127          call field_copy(this%p_list(local_idx), p)
 
  128          call field_copy(this%u_list(local_idx), u)
 
  129          call field_copy(this%v_list(local_idx), v)
 
  130          call field_copy(this%w_list(local_idx), w)
 
  131          do i_scalars = 1, this%n_scalars
 
  132             j = (local_idx - 1) * this%n_scalars + i_scalars
 
  133             s => neko_case%scalars%scalar_fields(i_scalars)%s
 
  134             call field_copy(this%s_list(j), s)
 
  140    local_idx = modulo(tstep, this%n_saves_memory) + 1
 
  141    call field_copy(p, this%p_list(local_idx))
 
  142    call field_copy(u, this%u_list(local_idx))
 
  143    call field_copy(v, this%v_list(local_idx))
 
  144    call field_copy(w, this%w_list(local_idx))
 
  145    do i_scalars = 1, this%n_scalars
 
  146       j = (local_idx - 1) * this%n_scalars + i_scalars
 
  147       s => neko_case%scalars%scalar_fields(i_scalars)%s
 
  148       call field_copy(s, this%s_list(j))
 
  152  end subroutine checkpoint_restore_linear
 
  154end submodule checkpoint_linear