Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
euler_res_cpu.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!
34!? It handles the time advancement of primitive variables using Runge-Kutta methods
35!? and evaluates the right-hand side terms of the Euler equations including artificial viscosity.
37 use euler_residual, only : euler_rhs_t
38 use field, only : field_t
39 use ax_product, only : ax_t
40 use coefs, only : coef_t
41 use gather_scatter, only : gs_t
42 use num_types, only : rp
43 use operators, only: div, rotate_cyc
44 use gs_ops, only : gs_op_add
47 use field_list, only : field_list_t
48 implicit none
49 private
50
51 type, public, extends(euler_rhs_t) :: euler_res_cpu_t
52 contains
53 procedure, nopass :: step => advance_primitive_variables_cpu
54 procedure, nopass :: evaluate_rhs_cpu
55 end type euler_res_cpu_t
56
57contains
72 subroutine advance_primitive_variables_cpu(rho_field, m_x, m_y, m_z, &
73 E, p, u, v, w, Ax, &
74 coef, gs, h, effective_visc, rk_scheme, dt)
75 type(field_t), intent(inout) :: rho_field, m_x, m_y, m_z, E
76 type(field_t), intent(in) :: p, u, v, w, h, effective_visc
77 class(ax_t), intent(inout) :: Ax
78 type(coef_t), intent(inout) :: coef
79 type(gs_t), intent(inout) :: gs
80 class(runge_kutta_time_scheme_t), intent(in) :: rk_scheme
81 real(kind=rp), intent(in) :: dt
82 integer :: n, s, i, j, k
83 type(field_t), pointer :: k_rho_1, k_rho_2, k_rho_3, k_rho_4, &
84 k_m_x_1, k_m_x_2, k_m_x_3, k_m_x_4, &
85 k_m_y_1, k_m_y_2, k_m_y_3, k_m_y_4, &
86 k_m_z_1, k_m_z_2, k_m_z_3, k_m_z_4, &
87 k_E_1, k_E_2, k_E_3, k_E_4, &
88 temp_rho, temp_m_x, temp_m_y, temp_m_z, temp_E
89 integer :: tmp_indices(25)
90 type(field_list_t) :: k_rho, k_m_x, k_m_y, k_m_z, k_E
91
92 n = p%dof%size()
93 s = rk_scheme%order
94 call neko_scratch_registry%request_field(k_rho_1, tmp_indices(1), .true.)
95 call neko_scratch_registry%request_field(k_rho_2, tmp_indices(2), .true.)
96 call neko_scratch_registry%request_field(k_rho_3, tmp_indices(3), .true.)
97 call neko_scratch_registry%request_field(k_rho_4, tmp_indices(4), .true.)
98 call neko_scratch_registry%request_field(k_m_x_1, tmp_indices(5), .true.)
99 call neko_scratch_registry%request_field(k_m_x_2, tmp_indices(6), .true.)
100 call neko_scratch_registry%request_field(k_m_x_3, tmp_indices(7), .true.)
101 call neko_scratch_registry%request_field(k_m_x_4, tmp_indices(8), .true.)
102 call neko_scratch_registry%request_field(k_m_y_1, tmp_indices(9), .true.)
103 call neko_scratch_registry%request_field(k_m_y_2, tmp_indices(10), .true.)
104 call neko_scratch_registry%request_field(k_m_y_3, tmp_indices(11), .true.)
105 call neko_scratch_registry%request_field(k_m_y_4, tmp_indices(12), .true.)
106 call neko_scratch_registry%request_field(k_m_z_1, tmp_indices(13), .true.)
107 call neko_scratch_registry%request_field(k_m_z_2, tmp_indices(14), .true.)
108 call neko_scratch_registry%request_field(k_m_z_3, tmp_indices(15), .true.)
109 call neko_scratch_registry%request_field(k_m_z_4, tmp_indices(16), .true.)
110 call neko_scratch_registry%request_field(k_e_1, tmp_indices(17), .true.)
111 call neko_scratch_registry%request_field(k_e_2, tmp_indices(18), .true.)
112 call neko_scratch_registry%request_field(k_e_3, tmp_indices(19), .true.)
113 call neko_scratch_registry%request_field(k_e_4, tmp_indices(20), .true.)
114 call neko_scratch_registry%request_field(temp_rho, tmp_indices(21), .false.)
115 call neko_scratch_registry%request_field(temp_m_x, tmp_indices(22), .false.)
116 call neko_scratch_registry%request_field(temp_m_y, tmp_indices(23), .false.)
117 call neko_scratch_registry%request_field(temp_m_z, tmp_indices(24), .false.)
118 call neko_scratch_registry%request_field(temp_e, tmp_indices(25), .false.)
119
120 ! Initialize Runge-Kutta stage variables for each conserved quantity
121 call k_rho%init(4)
122 call k_rho%assign(1, k_rho_1)
123 call k_rho%assign(2, k_rho_2)
124 call k_rho%assign(3, k_rho_3)
125 call k_rho%assign(4, k_rho_4)
126 call k_m_x%init(4)
127 call k_m_x%assign(1, k_m_x_1)
128 call k_m_x%assign(2, k_m_x_2)
129 call k_m_x%assign(3, k_m_x_3)
130 call k_m_x%assign(4, k_m_x_4)
131 call k_m_y%init(4)
132 call k_m_y%assign(1, k_m_y_1)
133 call k_m_y%assign(2, k_m_y_2)
134 call k_m_y%assign(3, k_m_y_3)
135 call k_m_y%assign(4, k_m_y_4)
136 call k_m_z%init(4)
137 call k_m_z%assign(1, k_m_z_1)
138 call k_m_z%assign(2, k_m_z_2)
139 call k_m_z%assign(3, k_m_z_3)
140 call k_m_z%assign(4, k_m_z_4)
141 call k_e%init(4)
142 call k_e%assign(1, k_e_1)
143 call k_e%assign(2, k_e_2)
144 call k_e%assign(3, k_e_3)
145 call k_e%assign(4, k_e_4)
146
147 ! Loop over Runge-Kutta stages
148 do i = 1, s
149 ! Copy current solution state to temporary arrays for this RK stage
150 do concurrent(k = 1:n)
151 temp_rho%x(k,1,1,1) = rho_field%x(k,1,1,1)
152 temp_m_x%x(k,1,1,1) = m_x%x(k,1,1,1)
153 temp_m_y%x(k,1,1,1) = m_y%x(k,1,1,1)
154 temp_m_z%x(k,1,1,1) = m_z%x(k,1,1,1)
155 temp_e%x(k,1,1,1) = e%x(k,1,1,1)
156 end do
157
158 ! Accumulate previous stage contributions using RK coefficients
159 do j = 1, i-1
160 do concurrent(k = 1:n)
161 temp_rho%x(k,1,1,1) = temp_rho%x(k,1,1,1) &
162 + dt * rk_scheme%coeffs_A(i, j) * k_rho%items(j)%ptr%x(k,1,1,1)
163 temp_m_x%x(k,1,1,1) = temp_m_x%x(k,1,1,1) &
164 + dt * rk_scheme%coeffs_A(i, j) * k_m_x%items(j)%ptr%x(k,1,1,1)
165 temp_m_y%x(k,1,1,1) = temp_m_y%x(k,1,1,1) &
166 + dt * rk_scheme%coeffs_A(i, j) * k_m_y%items(j)%ptr%x(k,1,1,1)
167 temp_m_z%x(k,1,1,1) = temp_m_z%x(k,1,1,1) &
168 + dt * rk_scheme%coeffs_A(i, j) * k_m_z%items(j)%ptr%x(k,1,1,1)
169 temp_e%x(k,1,1,1) = temp_e%x(k,1,1,1) &
170 + dt * rk_scheme%coeffs_A(i, j) * k_e%items(j)%ptr%x(k,1,1,1)
171 end do
172 end do
173
174 ! Evaluate RHS terms for current stage using intermediate solution values
175 call evaluate_rhs_cpu(k_rho%items(i)%ptr, k_m_x%items(i)%ptr, &
176 k_m_y%items(i)%ptr, k_m_z%items(i)%ptr, &
177 k_e%items(i)%ptr, &
178 temp_rho, temp_m_x, temp_m_y, temp_m_z, temp_e, &
179 p, u, v, w, ax, &
180 coef, gs, h, effective_visc)
181 end do
182
183 ! Update the solution
184 do i = 1, s
185 do concurrent(k = 1:n)
186 rho_field%x(k,1,1,1) = rho_field%x(k,1,1,1) &
187 + dt * rk_scheme%coeffs_b(i) * k_rho%items(i)%ptr%x(k,1,1,1)
188 m_x%x(k,1,1,1) = m_x%x(k,1,1,1) &
189 + dt * rk_scheme%coeffs_b(i) * k_m_x%items(i)%ptr%x(k,1,1,1)
190 m_y%x(k,1,1,1) = m_y%x(k,1,1,1) &
191 + dt * rk_scheme%coeffs_b(i) * k_m_y%items(i)%ptr%x(k,1,1,1)
192 m_z%x(k,1,1,1) = m_z%x(k,1,1,1) &
193 + dt * rk_scheme%coeffs_b(i) * k_m_z%items(i)%ptr%x(k,1,1,1)
194 e%x(k,1,1,1) = e%x(k,1,1,1) &
195 + dt * rk_scheme%coeffs_b(i) * k_e%items(i)%ptr%x(k,1,1,1)
196 end do
197 end do
198
199 call neko_scratch_registry%relinquish_field(tmp_indices)
200
202
224 subroutine evaluate_rhs_cpu(rhs_rho_field, rhs_m_x, rhs_m_y, rhs_m_z, rhs_E, &
225 rho_field, m_x, m_y, m_z, E, p, u, v, w, Ax, &
226 coef, gs, h, effective_visc)
227 type(field_t), intent(inout) :: rhs_rho_field, &
228 rhs_m_x, rhs_m_y, rhs_m_z, rhs_E
229 type(field_t), intent(inout) :: rho_field, m_x, m_y, m_z, E
230 type(field_t), intent(in) :: p, u, v, w, h, effective_visc
231 class(ax_t), intent(inout) :: Ax
232 type(coef_t), intent(inout) :: coef
233 type(gs_t), intent(inout) :: gs
234 integer :: i, n
235 type(field_t), pointer :: f_x, f_y, f_z, &
236 visc_rho, visc_m_x, visc_m_y, visc_m_z, visc_E
237 integer :: tmp_indices(8)
238
239 n = coef%dof%size()
240 call neko_scratch_registry%request_field(f_x, tmp_indices(1), .false.)
241 call neko_scratch_registry%request_field(f_y, tmp_indices(2), .false.)
242 call neko_scratch_registry%request_field(f_z, tmp_indices(3), .false.)
243
245 ! Compute density flux divergence
246 call div(rhs_rho_field%x, m_x%x, m_y%x, m_z%x, coef)
247
249 ! Compute momentum flux divergences
250 ! m_x
251 do concurrent(i = 1:n)
252 f_x%x(i,1,1,1) = m_x%x(i,1,1,1) * m_x%x(i,1,1,1) / rho_field%x(i, 1, 1, 1) &
253 + p%x(i,1,1,1)
254 f_y%x(i,1,1,1) = m_x%x(i,1,1,1) * m_y%x(i,1,1,1) / rho_field%x(i, 1, 1, 1)
255 f_z%x(i,1,1,1) = m_x%x(i,1,1,1) * m_z%x(i,1,1,1) / rho_field%x(i, 1, 1, 1)
256 end do
257 call div(rhs_m_x%x, f_x%x, f_y%x, f_z%x, coef)
258 ! m_y
259 do concurrent(i = 1:n)
260 f_x%x(i,1,1,1) = m_y%x(i,1,1,1) * m_x%x(i,1,1,1) / rho_field%x(i, 1, 1, 1)
261 f_y%x(i,1,1,1) = m_y%x(i,1,1,1) * m_y%x(i,1,1,1) / rho_field%x(i, 1, 1, 1) &
262 + p%x(i,1,1,1)
263 f_z%x(i,1,1,1) = m_y%x(i,1,1,1) * m_z%x(i,1,1,1) / rho_field%x(i, 1, 1, 1)
264 end do
265 call div(rhs_m_y%x, f_x%x, f_y%x, f_z%x, coef)
266 ! m_z
267 do concurrent(i = 1:n)
268 f_x%x(i,1,1,1) = m_z%x(i,1,1,1) * m_x%x(i,1,1,1) / rho_field%x(i, 1, 1, 1)
269 f_y%x(i,1,1,1) = m_z%x(i,1,1,1) * m_y%x(i,1,1,1) / rho_field%x(i, 1, 1, 1)
270 f_z%x(i,1,1,1) = m_z%x(i,1,1,1) * m_z%x(i,1,1,1) / rho_field%x(i, 1, 1, 1) &
271 + p%x(i,1,1,1)
272 end do
273 call div(rhs_m_z%x, f_x%x, f_y%x, f_z%x, coef)
274
276 ! Compute energy flux divergence
277 do concurrent(i = 1:n)
278 f_x%x(i,1,1,1) = (e%x(i,1,1,1) + p%x(i,1,1,1)) * u%x(i,1,1,1)
279 f_y%x(i,1,1,1) = (e%x(i,1,1,1) + p%x(i,1,1,1)) * v%x(i,1,1,1)
280 f_z%x(i,1,1,1) = (e%x(i,1,1,1) + p%x(i,1,1,1)) * w%x(i,1,1,1)
281 end do
282 call div(rhs_e%x, f_x%x, f_y%x, f_z%x, coef)
283
284 ! gs
285 call gs%op(rhs_rho_field, gs_op_add)
286 call rotate_cyc(rhs_m_x%x, rhs_m_y%x, rhs_m_z%x, 1, coef)
287 call gs%op(rhs_m_x, gs_op_add)
288 call gs%op(rhs_m_y, gs_op_add)
289 call gs%op(rhs_m_z, gs_op_add)
290 call rotate_cyc(rhs_m_x%x, rhs_m_y%x, rhs_m_z%x, 0, coef)
291 call gs%op(rhs_e, gs_op_add)
292 do concurrent(i = 1:rhs_e%dof%size())
293 rhs_rho_field%x(i,1,1,1) = rhs_rho_field%x(i,1,1,1) * coef%mult(i,1,1,1)
294 rhs_m_x%x(i,1,1,1) = rhs_m_x%x(i,1,1,1) * coef%mult(i,1,1,1)
295 rhs_m_y%x(i,1,1,1) = rhs_m_y%x(i,1,1,1) * coef%mult(i,1,1,1)
296 rhs_m_z%x(i,1,1,1) = rhs_m_z%x(i,1,1,1) * coef%mult(i,1,1,1)
297 rhs_e%x(i,1,1,1) = rhs_e%x(i,1,1,1) * coef%mult(i,1,1,1)
298 end do
299
300 call neko_scratch_registry%request_field(visc_rho, tmp_indices(4), .false.)
301 call neko_scratch_registry%request_field(visc_m_x, tmp_indices(5), .false.)
302 call neko_scratch_registry%request_field(visc_m_y, tmp_indices(6), .false.)
303 call neko_scratch_registry%request_field(visc_m_z, tmp_indices(7), .false.)
304 call neko_scratch_registry%request_field(visc_e, tmp_indices(8), .false.)
305
306 ! Set h1 coefficient to the effective viscosity for the Laplacian operator
307 do concurrent(i = 1:n)
308 coef%h1(i,1,1,1) = effective_visc%x(i,1,1,1)
309 end do
310
311 ! Calculate artificial diffusion with variable viscosity
312 call ax%compute(visc_rho%x, rho_field%x, coef, p%msh, p%Xh)
313 call ax%compute(visc_m_x%x, m_x%x, coef, p%msh, p%Xh)
314 call ax%compute(visc_m_y%x, m_y%x, coef, p%msh, p%Xh)
315 call ax%compute(visc_m_z%x, m_z%x, coef, p%msh, p%Xh)
316 call ax%compute(visc_e%x, e%x, coef, p%msh, p%Xh)
317
318 ! Reset h1 coefficient back to 1.0 for other operations
319 do concurrent(i = 1:n)
320 coef%h1(i,1,1,1) = 1.0_rp
321 end do
322
323 ! gs
324 call gs%op(visc_rho, gs_op_add)
325 call rotate_cyc(visc_m_x%x, visc_m_y%x, visc_m_z%x, 1, coef)
326 call gs%op(visc_m_x, gs_op_add)
327 call gs%op(visc_m_y, gs_op_add)
328 call gs%op(visc_m_z, gs_op_add)
329 call rotate_cyc(visc_m_x%x, visc_m_y%x, visc_m_z%x, 0, coef)
330 call gs%op(visc_e, gs_op_add)
331
332 ! Move div to the rhs and apply artificial viscosity
333 ! The viscosity coefficient is already included in the Laplacian operator
334 do concurrent(i = 1:n)
335 rhs_rho_field%x(i,1,1,1) = -rhs_rho_field%x(i,1,1,1) &
336 - coef%Binv(i,1,1,1) * visc_rho%x(i,1,1,1)
337 rhs_m_x%x(i,1,1,1) = -rhs_m_x%x(i,1,1,1) &
338 - coef%Binv(i,1,1,1) * visc_m_x%x(i,1,1,1)
339 rhs_m_y%x(i,1,1,1) = -rhs_m_y%x(i,1,1,1) &
340 - coef%Binv(i,1,1,1) * visc_m_y%x(i,1,1,1)
341 rhs_m_z%x(i,1,1,1) = -rhs_m_z%x(i,1,1,1) &
342 - coef%Binv(i,1,1,1) * visc_m_z%x(i,1,1,1)
343 rhs_e%x(i,1,1,1) = -rhs_e%x(i,1,1,1) &
344 - coef%Binv(i,1,1,1) * visc_e%x(i,1,1,1)
345 end do
346
347 call neko_scratch_registry%relinquish_field(tmp_indices)
348 end subroutine evaluate_rhs_cpu
349
350end module euler_res_cpu
Defines a Matrix-vector product.
Definition ax.f90:34
Coefficients.
Definition coef.f90:34
This module implements CPU-based residual calculations for the Euler equations.
subroutine evaluate_rhs_cpu(rhs_rho_field, rhs_m_x, rhs_m_y, rhs_m_z, rhs_e, rho_field, m_x, m_y, m_z, e, p, u, v, w, ax, coef, gs, h, effective_visc)
Evaluates the right-hand side of the Euler equations Computes fluxes, divergence terms,...
subroutine advance_primitive_variables_cpu(rho_field, m_x, m_y, m_z, e, p, u, v, w, ax, coef, gs, h, effective_visc, rk_scheme, dt)
Advances the primitive variables (density, momentum, energy) in time using a Runge-Kutta scheme.
Defines a field.
Definition field.f90:34
Gather-scatter.
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Base type for a matrix-vector product providing .
Definition ax.f90:43
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
Abstract type to compute rhs.
Definition euler_res.f90:44
field_list_t, To be able to group fields together