Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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!
37 use num_types, only : rp
42 use vector, only : vector_t
43 use, intrinsic :: iso_c_binding
44 implicit none
45 private
73 type(ext_time_scheme_t) :: ext
74 type(ab_time_scheme_t) :: ab
75 type(bdf_time_scheme_t) :: bdf
76
78 type(vector_t) :: advection_coeffs
80 type(vector_t) :: diffusion_coeffs
83 integer :: ndiff = 0
86 integer :: nadv = 0
88 integer :: advection_time_order = 3
90 integer :: diffusion_time_order
92 type(c_ptr) :: advection_coeffs_d = c_null_ptr
94 type(c_ptr) :: diffusion_coeffs_d = c_null_ptr
95
96 contains
98 procedure, pass(this) :: init => time_scheme_controller_init
100 procedure, pass(this) :: free => time_scheme_controller_free
102 procedure, pass(this) :: set_coeffs => &
105
106contains
107
111 subroutine time_scheme_controller_init(this, torder)
112 implicit none
113 class(time_scheme_controller_t) :: this
114 integer :: torder
115
116 call this%free()
117
118 this%diffusion_time_order = torder
119
120 ! Force 1st order advection when diffusion is 1st order
121 if (torder .eq. 1) then
122 this%advection_time_order = 1
123 end if
124
125 call this%advection_coeffs%init(4)
126 call this%diffusion_coeffs%init(4)
127 end subroutine time_scheme_controller_init
128
131 implicit none
132 class(time_scheme_controller_t) :: this
133
134 call this%advection_coeffs%free()
135 call this%diffusion_coeffs%free()
136 end subroutine time_scheme_controller_free
137
142 implicit none
143 class(time_scheme_controller_t) :: this
144 real(kind=rp), intent(inout), dimension(10) :: dt
145 real(kind=rp), dimension(4) :: adv_coeffs_old
146 real(kind=rp), dimension(4) :: diff_coeffs_old
147
148 associate( &
149 nadv => this%nadv, &
150 ndiff => this%ndiff, &
151 adv_coeffs => this%advection_coeffs, &
152 adv_coeffs_d => this%advection_coeffs_d, &
153 diff_coeffs => this%diffusion_coeffs, &
154 diff_coeffs_d => this%diffusion_coeffs_d)
155
156 adv_coeffs_old = adv_coeffs%x
157 diff_coeffs_old = diff_coeffs%x
158
159 ! Increment the order of the scheme if below time_order
160 ndiff = ndiff + 1
161 ndiff = min(ndiff, this%diffusion_time_order)
162 nadv = nadv + 1
163 nadv = min(nadv, this%advection_time_order)
164
165 call this%bdf%compute_coeffs(diff_coeffs%x, dt, ndiff)
166
167 if (nadv .eq. 1) then
168 ! Forward euler
169 call this%ext%compute_coeffs(adv_coeffs%x, dt, nadv)
170 else if (nadv .eq. 2) then
171 if (ndiff .eq. 1) then
172 ! 2nd order Adam-Bashforth, currently never used
173 call this%ab%compute_coeffs(adv_coeffs%x, dt, nadv)
174 else
175 ! Linear extrapolation
176 call this%ext%compute_coeffs(adv_coeffs%x, dt, nadv)
177 end if
178 else if (nadv .eq. 3) then
179 if (ndiff .eq. 1) then
180 ! 3rd order Adam-Bashforth, currently never used
181 call this%ab%compute_coeffs(adv_coeffs%x, dt, nadv)
182 else if (ndiff .eq. 2) then
183 ! The modified EXT scheme
184 call this%ext%compute_modified_coeffs(adv_coeffs%x, dt)
185 else
186 ! Quadratic extrapolation
187 call this%ext%compute_coeffs(adv_coeffs%x, dt, nadv)
188 end if
189 end if
190
191 if (neko_bcknd_device .eq. 1) then
192 if (maxval(abs(adv_coeffs%x - adv_coeffs_old)) .gt. 1e-10_rp) then
193 call adv_coeffs%copy_from(host_to_device, sync = .false.)
194 end if
195
196 if (maxval(abs(diff_coeffs%x - diff_coeffs_old)) .gt. 1e-10_rp) then
197 call diff_coeffs%copy_from(host_to_device, sync = .false.)
198 end if
199 end if
200 end associate
201
203end module time_scheme_controller
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
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:225
Explicit extrapolation scheme for time integration.
Build configurations.
integer, parameter neko_bcknd_device
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.
Defines a vector.
Definition vector.f90:34
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 ...