Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
time_step_controller.f90
Go to the documentation of this file.
1! Copyright (c) 2022, 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!
35 use num_types, only : rp
36 use logger, only : neko_log, log_size
37 use json_module, only : json_file
39 use time_state, only : time_state_t
41 use mpi_f08, only : mpi_min, mpi_in_place
42 implicit none
43 private
44
46 type, public :: time_step_controller_t
47 logical :: is_variable_dt
48 real(kind=rp) :: cfl_trg = 0.0_rp
49 real(kind=rp) :: cfl_avg = 0.0_rp
50 real(kind=rp) :: init_dt = huge(0.0_rp)
51 real(kind=rp) :: max_dt = 0.0_rp
52 real(kind=rp) :: min_dt = 0.0_rp
53 integer :: max_update_frequency = 0
54 integer :: min_update_frequency = 0
55 integer :: dt_last_change = 0
56 real(kind=rp) :: alpha = 0.0_rp
57 real(kind=rp) :: max_dt_increase_factor = 0.0_rp
58 real(kind=rp) :: min_dt_decrease_factor = 0.0_rp
59 real(kind=rp) :: dev_tol = 0.0_rp
60 contains
62 procedure, pass(this) :: init => time_step_controller_init
64 procedure, pass(this) :: set_dt => time_step_controller_set_dt
65
67
68contains
69
72 subroutine time_step_controller_init(this, params)
73 class(time_step_controller_t), intent(inout) :: this
74 type(json_file), intent(inout) :: params
75
76 this%dt_last_change = -1
77 call json_get_or_default(params, 'variable_timestep', &
78 this%is_variable_dt, .false.)
79 if (this%is_variable_dt) then
80 call json_get_or_lookup_or_default(params, 'target_cfl', &
81 this%cfl_trg, 0.4_rp)
82 call json_get_or_lookup_or_default(params, 'timestep', &
83 this%init_dt, huge(0.0_rp))
84 call json_get_or_lookup_or_default(params, 'max_timestep', &
85 this%max_dt, huge(0.0_rp))
86 call json_get_or_lookup_or_default(params, 'min_timestep', &
87 this%min_dt, 0.0_rp)
88 call json_get_or_lookup_or_default(params, 'max_update_frequency',&
89 this%max_update_frequency, 0)
90 call json_get_or_lookup_or_default(params, 'min_update_frequency',&
91 this%min_update_frequency, huge(0))
92 call json_get_or_lookup_or_default(params, 'cfl_running_avg_coeff', &
93 this%alpha, 0.5_rp)
94 call json_get_or_lookup_or_default(params, 'max_dt_increase_factor', &
95 this%max_dt_increase_factor, 1.2_rp)
96 call json_get_or_lookup_or_default(params, 'min_dt_decrease_factor', &
97 this%min_dt_decrease_factor, 0.5_rp)
98 call json_get_or_lookup_or_default(params, 'cfl_deviation_tolerance', &
99 this%dev_tol, 0.2_rp)
100 end if
101
102 end subroutine time_step_controller_init
103
111 subroutine time_step_controller_set_dt(this, time, cfl)
112 class(time_step_controller_t), intent(inout) :: this
113 type(time_state_t), intent(inout) :: time
114 real(kind=rp), intent(in) :: cfl
115 real(kind=rp) :: dt_old, scaling_factor, global_min_dt
116 character(len=LOG_SIZE) :: log_buf
117 integer :: ierr
118
119 ! Check if variable dt is requested
120 if (.not. this%is_variable_dt) return
121
122 ! Reset the average cfl if it is the first time step since the last change
123 if (this%dt_last_change .eq. 0) then
124 this%cfl_avg = cfl
125 end if
126
127 if (this%dt_last_change .eq. -1) then
128
129 ! Set the first dt for desired cfl, or use the provided initial dt if it
130 ! is smaller. Then clamp between max and min dt if provided.
131 time%dt = min(this%cfl_trg / cfl * time%dt, this%init_dt)
132 time%dt = max(min(time%dt, this%max_dt), this%min_dt)
133 this%dt_last_change = 0
134 this%cfl_avg = cfl
135
136 else
137 ! Calculate the average of cfl over the desired interval
138 this%cfl_avg = this%alpha * cfl + (1 - this%alpha) * this%cfl_avg
139
140 if (abs(this%cfl_avg - this%cfl_trg) .ge. this%dev_tol * this%cfl_trg &
141 .and. this%dt_last_change .ge. this%max_update_frequency &
142 .or. this%dt_last_change .ge. this%min_update_frequency) then
143
144 if (this%cfl_trg/cfl .ge. 1) then
145 ! increase of time step
146 scaling_factor = min(this%max_dt_increase_factor, this%cfl_trg/cfl)
147 else
148 ! reduction of time step
149 scaling_factor = max(this%min_dt_decrease_factor, this%cfl_trg/cfl)
150 end if
151
152 dt_old = time%dt
153 time%dt = scaling_factor * dt_old
154 time%dt = max(min(time%dt, this%max_dt), this%min_dt)
155
156 write(log_buf, '(A,E15.7,1x,A,E15.7)') &
157 'Average CFL:', this%cfl_avg, &
158 'Target CFL:', this%cfl_trg
159 call neko_log%message(log_buf)
160
161 write(log_buf, '(A,E15.7,1x,A,E15.7)') 'Old dt:', dt_old, &
162 'New dt:', time%dt
163 call neko_log%message(log_buf)
164
165 this%dt_last_change = 0
166
167 else
168 this%dt_last_change = this%dt_last_change + 1
169 end if
170 end if
171
172 ! If running in mpmd, the new dt is the minimum across simulations
173 if (pe_size .ne. global_pe_size) then
174 global_min_dt = time%dt
175 call mpi_allreduce(mpi_in_place, global_min_dt, 1, mpi_real_precision, &
176 mpi_min, neko_global_comm, ierr)
177
178 ! If my dt is larger that the global min, mark a change
179 if (time%dt .gt. global_min_dt) then
180 time%dt = global_min_dt
181 this%dt_last_change = 0
182 end if
183 end if
184
185 end subroutine time_step_controller_set_dt
186
187
188
189
190end module time_step_controller
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition comm.F90:1
type(mpi_comm), public neko_global_comm
Definition comm.F90:44
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public global_pe_size
Global MPI size of communicator.
Definition comm.F90:68
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
Implements type time_step_controller.
subroutine time_step_controller_init(this, params)
Constructor.
subroutine time_step_controller_set_dt(this, time, cfl)
Set new dt based on cfl if requested.
A struct that contains all info about the time, expand as needed.
#define max(a, b)
Definition tensor.cu:40