Neko  0.8.1
A portable framework for high-order spectral element flow simulations
time_scheme_controller.f90
Go to the documentation of this file.
1 ! Copyright (c) 2023, 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 !
36  use neko_config
37  use num_types, only : rp
42  use, intrinsic :: iso_c_binding
43  implicit none
44  private
71  type, public :: time_scheme_controller_t
72  type(ext_time_scheme_t) :: ext
73  type(ab_time_scheme_t) :: ab
74  type(bdf_time_scheme_t) :: bdf
75 
77  real(kind=rp) :: advection_coeffs(4) = 0
79  real(kind=rp) :: diffusion_coeffs(4) = 0
82  integer :: ndiff = 0
85  integer :: nadv = 0
87  integer :: advection_time_order = 3
89  integer :: diffusion_time_order
91  type(c_ptr) :: advection_coeffs_d = c_null_ptr
93  type(c_ptr) :: diffusion_coeffs_d = c_null_ptr
94 
95  contains
97  procedure, pass(this) :: init => time_scheme_controller_init
99  procedure, pass(this) :: free => time_scheme_controller_free
101  procedure, pass(this) :: set_coeffs => &
103  end type time_scheme_controller_t
104 
105 contains
106 
110  subroutine time_scheme_controller_init(this, torder)
111  implicit none
112  class(time_scheme_controller_t) :: this
113  integer :: torder
114 
115  this%diffusion_time_order = torder
116 
117  ! Force 1st order advection when diffusion is 1st order
118  if (torder .eq. 1) then
119  this%advection_time_order = 1
120  end if
121 
122  if (neko_bcknd_device .eq. 1) then
123  call device_map(this%advection_coeffs, this%advection_coeffs_d, 4)
124  call device_map(this%diffusion_coeffs, this%diffusion_coeffs_d, 4)
125  end if
126  end subroutine time_scheme_controller_init
127 
130  implicit none
131  class(time_scheme_controller_t) :: this
132 
133  if (c_associated(this%advection_coeffs_d)) then
134  call device_free(this%advection_coeffs_d)
135  end if
136  if (c_associated(this%diffusion_coeffs_d)) then
137  call device_free(this%diffusion_coeffs_d)
138  end if
139  end subroutine time_scheme_controller_free
140 
145  implicit none
146  class(time_scheme_controller_t) :: this
147  real(kind=rp), intent(inout), dimension(10) :: dt
148  real(kind=rp), dimension(4) :: adv_coeffs_old
149  real(kind=rp), dimension(4) :: diff_coeffs_old
150 
151  associate( &
152  nadv => this%nadv, &
153  ndiff => this%ndiff, &
154  adv_coeffs => this%advection_coeffs, &
155  adv_coeffs_d => this%advection_coeffs_d, &
156  diff_coeffs => this%diffusion_coeffs, &
157  diff_coeffs_d => this%diffusion_coeffs_d)
158 
159  adv_coeffs_old = adv_coeffs
160  diff_coeffs_old = diff_coeffs
161 
162  ! Increment the order of the scheme if below time_order
163  ndiff = ndiff + 1
164  ndiff = min(ndiff, this%diffusion_time_order)
165  nadv = nadv + 1
166  nadv = min(nadv, this%advection_time_order)
167 
168  call this%bdf%compute_coeffs(diff_coeffs, dt, ndiff)
169 
170  if (nadv .eq. 1) then
171  ! Forward euler
172  call this%ext%compute_coeffs(adv_coeffs, dt, nadv)
173  else if (nadv .eq. 2) then
174  if (ndiff .eq. 1) then
175  ! 2nd order Adam-Bashforth, currently never used
176  call this%ab%compute_coeffs(adv_coeffs, dt, nadv)
177  else
178  ! Linear extrapolation
179  call this%ext%compute_coeffs(adv_coeffs, dt, nadv)
180  end if
181  else if (nadv .eq. 3) then
182  if (ndiff .eq. 1) then
183  ! 3rd order Adam-Bashforth, currently never used
184  call this%ab%compute_coeffs(adv_coeffs, dt, nadv)
185  else if (ndiff .eq. 2) then
186  ! The modified EXT scheme
187  call this%ext%compute_modified_coeffs(adv_coeffs, dt)
188  else
189  ! Quadratic extrapolation
190  call this%ext%compute_coeffs(adv_coeffs, dt, nadv)
191  end if
192  end if
193 
194 
195  if (c_associated(adv_coeffs_d)) then
196  if (maxval(abs(adv_coeffs - adv_coeffs_old)) .gt. 1e-10_rp) then
197  call device_memcpy(adv_coeffs, adv_coeffs_d, 4, &
198  host_to_device, sync=.false.)
199  end if
200  end if
201 
202  if (c_associated(diff_coeffs_d)) then
203  if (maxval(abs(diff_coeffs - diff_coeffs_old)) .gt. 1e-10_rp) then
204  call device_memcpy(diff_coeffs, diff_coeffs_d, 4, &
205  host_to_device, sync=.false.)
206  end if
207  end if
208  end associate
209 
210  end subroutine time_scheme_controller_set_coeffs
211 end module time_scheme_controller
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
Copy data between host and device (or device and device)
Definition: device.F90:51
Adam-Bashforth scheme for time integration.
Backward-differencing scheme for time integration.
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public host_to_device
Definition: device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition: device.F90:172
Explicit extrapolation scheme for time integration.
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Compound scheme for the advection and diffusion operators in a transport equation.
subroutine time_scheme_controller_free(this)
Destructor.
subroutine time_scheme_controller_init(this, torder)
Contructor.
subroutine time_scheme_controller_set_coeffs(this, dt)
Set the time coefficients.
Explicit Adam-Bashforth scheme for time integration.
Implicit backward-differencing scheme for time integration.
Explicit extrapolation scheme for time integration.
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...