Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
runge_kutta_scheme.f90
Go to the documentation of this file.
1! Copyright (c) 2025, 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 time_scheme, only: time_scheme_t
37 use math, only : rzero
38 use utils, only : neko_error
40 use, intrinsic :: iso_c_binding, only : c_ptr, c_associated, c_null_ptr
41 implicit none
42 private
43
45 real(kind=rp), allocatable :: coeffs_a(:, :)
46 real(kind=rp), allocatable :: coeffs_b(:)
47 real(kind=rp), allocatable :: coeffs_c(:)
48 type(c_ptr) :: coeffs_a_d = c_null_ptr
49 type(c_ptr) :: coeffs_b_d = c_null_ptr
50 type(c_ptr) :: coeffs_c_d = c_null_ptr
51 integer :: order
52 contains
54 procedure, pass(this) :: init => runge_kutta_scheme_coeffs_init
56 procedure, pass(this) :: free => runge_kutta_scheme_coeffs_free
58
59contains
67 subroutine runge_kutta_scheme_coeffs_init(this, order)
68 class(runge_kutta_time_scheme_t), intent(inout) :: this
69 integer, intent(in) :: order
70 integer :: s
71
72 ! Assume number of stages is equal to the order
73 s = order
74 this%order = order
75 allocate(this%coeffs_A(order, order))
76 allocate(this%coeffs_b(order))
77 allocate(this%coeffs_c(order))
78
79 associate( &
80 coeffs_a => this%coeffs_A, &
81 coeffs_b => this%coeffs_b, &
82 coeffs_c => this%coeffs_c, &
83 coeffs_a_d => this%coeffs_A_d, &
84 coeffs_b_d => this%coeffs_b_d, &
85 coeffs_c_d => this%coeffs_c_d)
86
87 if (neko_bcknd_device .eq. 1) then
88 call device_map(coeffs_a, coeffs_a_d, s)
89 call device_map(coeffs_b, coeffs_b_d, s)
90 call device_map(coeffs_c, coeffs_c_d, s)
91 end if
92
93 select case (order)
94 case (1)
96 coeffs_a(1, 1) = 0.0_rp
97 coeffs_b(1) = 1.0_rp
98 coeffs_c(1) = 0.0_rp
99 case (2)
101 coeffs_a = 0.0_rp
102 coeffs_a(2, 1) = 1.0_rp
103
104 coeffs_b(1) = 0.5_rp
105 coeffs_b(2) = 0.5_rp
106
107 coeffs_c(1) = 0.0_rp
108 coeffs_c(2) = 1.0_rp
109 case (3)
111 coeffs_a = 0.0_rp
112 coeffs_a(2, 1) = 1.0_rp
113 coeffs_a(3, 1) = 0.25_rp
114 coeffs_a(3, 2) = 0.25_rp
115
116 coeffs_b(1) = 1.0_rp / 6.0_rp
117 coeffs_b(2) = 1.0_rp / 6.0_rp
118 coeffs_b(3) = 2.0_rp / 3.0_rp
119
120 coeffs_c(1) = 0.0_rp
121 coeffs_c(2) = 1.0_rp
122 coeffs_c(3) = 0.5_rp
123 case (4)
125 coeffs_a = 0.0_rp
126 coeffs_a(2, 1) = 0.5_rp
127 coeffs_a(3, 2) = 0.5_rp
128 coeffs_a(4, 3) = 1.0_rp
129
130 coeffs_b(1) = 1.0_rp / 6.0_rp
131 coeffs_b(2) = 1.0_rp / 3.0_rp
132 coeffs_b(3) = 1.0_rp / 3.0_rp
133 coeffs_b(4) = 1.0_rp / 6.0_rp
134
135 coeffs_c(1) = 0.0_rp
136 coeffs_c(2) = 0.5_rp
137 coeffs_c(3) = 0.5_rp
138 coeffs_c(4) = 1.0_rp
139 case default
140 call neko_error("The time order must be 1 to 4.")
141 end select
142
143 if (c_associated(coeffs_a_d)) then
144 call device_memcpy(coeffs_a, coeffs_a_d, s, &
145 host_to_device, sync = .false.)
146 end if
147
148 if (c_associated(coeffs_b_d)) then
149 call device_memcpy(coeffs_b, coeffs_b_d, s, &
150 host_to_device, sync = .false.)
151 end if
152
153 if (c_associated(coeffs_c_d)) then
154 call device_memcpy(coeffs_c, coeffs_c_d, s, &
155 host_to_device, sync = .false.)
156 end if
157 end associate
158
159 end subroutine runge_kutta_scheme_coeffs_init
160
164 implicit none
165 class(runge_kutta_time_scheme_t) :: this
166
167 if (c_associated(this%coeffs_A_d)) then
168 call device_free(this%coeffs_A_d)
169 end if
170 if (c_associated(this%coeffs_b_d)) then
171 call device_free(this%coeffs_b_d)
172 end if
173 if (c_associated(this%coeffs_c_d)) then
174 call device_free(this%coeffs_c_d)
175 end if
176 end subroutine runge_kutta_scheme_coeffs_free
177
Map a Fortran array to a device (allocate and associate)
Definition device.F90:71
Copy data between host and device (or device and device)
Definition device.F90:65
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:46
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:200
Definition math.f90:60
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:194
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
subroutine runge_kutta_scheme_coeffs_free(this)
Destructor - deallocates device memory.
subroutine runge_kutta_scheme_coeffs_init(this, order)
Constructor for Runge-Kutta time scheme.
Base class for time integration schemes.
Utilities.
Definition utils.f90:35
Base abstract class for time integration schemes.