Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
compressible_ops_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!
35 use num_types, only : rp
36 implicit none
37 private
38
42
43
44contains
45
48 u, v, w, gamma, p, rho, n)
49 integer, intent(in) :: n
50 real(kind=rp), intent(in) :: gamma
51 real(kind=rp), dimension(n), intent(in) :: u, v, w, p, rho
52 real(kind=rp), dimension(n), intent(inout) :: max_wave_speed
53 integer :: i
54 real(kind=rp) :: vel_mag, sound_speed
55
56 ! Compute maximum wave speed:
57 ! |u| + c = sqrt(u^2 + v^2 + w^2) + sqrt(gamma * p / rho)
58 !OCL NORECURRENCE, NOVREC, NOALIAS
59 !DIR$ CONCURRENT
60 !GCC$ ivdep
61 !$omp parallel do simd private(vel_mag, sound_speed)
62 do i = 1, n
63 vel_mag = sqrt(u(i)*u(i) + v(i)*v(i) + &
64 w(i)*w(i))
65 sound_speed = sqrt(gamma * p(i) / rho(i))
66 max_wave_speed(i) = vel_mag + sound_speed
67 end do
68 !$omp end parallel do simd
69
71
74 subroutine compressible_ops_cpu_compute_entropy(S, p, rho, gamma, n)
75 integer, intent(in) :: n
76 real(kind=rp), intent(in) :: gamma
77 real(kind=rp), dimension(n), intent(in) :: p, rho
78 real(kind=rp), dimension(n), intent(inout) :: s
79 integer :: i
80 real(kind=rp) :: inv_gamma_m1
81
82 inv_gamma_m1 = 1.0_rp / (gamma - 1.0_rp)
83
84 ! Compute entropy: S = 1/(gamma-1) * rho * (log(p) - gamma * log(rho))
85 !OCL NORECURRENCE, NOVREC, NOALIAS
86 !DIR$ CONCURRENT
87 !GCC$ ivdep
88 !$omp parallel do simd
89 do i = 1, n
90 s(i) = inv_gamma_m1 * rho(i) * &
91 (log(p(i)) - gamma * log(rho(i)))
92 end do
93 !$omp end parallel do simd
94
96
98 subroutine compressible_ops_cpu_update_uvw(u, v, w, m_x, m_y, m_z, rho, n)
99 integer, intent(in) :: n
100 real(kind=rp), dimension(n), intent(inout) :: u, v, w
101 real(kind=rp), dimension(n), intent(in) :: m_x, m_y, m_z, rho
102 integer :: i
103
104 !OCL NORECURRENCE, NOVREC, NOALIAS
105 !DIR$ CONCURRENT
106 !GCC$ ivdep
107 !$omp parallel do simd
108 do i = 1, n
109 u(i) = m_x(i) / rho(i)
110 v(i) = m_y(i) / rho(i)
111 w(i) = m_z(i) / rho(i)
112 end do
113 !$omp end parallel do simd
114
116
118 subroutine compressible_ops_cpu_update_mxyz_p_ruvw(m_x, m_y, m_z, p, ruvw, &
119 u, v, w, E, rho, gamma, n)
120 integer, intent(in) :: n
121 real(kind=rp), dimension(n), intent(inout) :: m_x, m_y, m_z, p, ruvw
122 real(kind=rp), dimension(n), intent(in) :: u, v, w, e, rho
123 real(kind=rp), intent(in) :: gamma
124 real(kind=rp) :: tmp
125 integer :: i
126
127 !OCL NORECURRENCE, NOVREC, NOALIAS
128 !DIR$ CONCURRENT
129 !GCC$ ivdep
130 !$omp parallel do simd private(tmp)
131 do i = 1, n
132 m_x(i) = u(i) * rho(i)
133 m_y(i) = v(i) * rho(i)
134 m_z(i) = w(i) * rho(i)
135 tmp = 0.5_rp * rho(i) * (u(i)**2 + v(i)**2 + w(i)**2)
136 p(i) = (gamma - 1.0_rp) * (e(i) - tmp)
137 ruvw(i) = tmp
138 end do
139 !$omp end parallel do simd
140
142
144 subroutine compressible_ops_cpu_update_e(E, p, ruvw, gamma, n)
145 integer, intent(in) :: n
146 real(kind=rp), dimension(n), intent(inout) :: e, p
147 ! ruvw = 0.5 * rho * (u^2 + v^2 + w^2)
148 real(kind=rp), dimension(n), intent(in) :: ruvw
149 real(kind=rp), intent(in) :: gamma
150 integer :: i
151 real(kind=rp) :: inv_gamma_m1
152
153 inv_gamma_m1 = 1.0_rp / (gamma - 1.0_rp)
154
155 !OCL NORECURRENCE, NOVREC, NOALIAS
156 !DIR$ CONCURRENT
157 !GCC$ ivdep
158 !$omp parallel do simd
159 do i = 1, n
160 ! Ensure pressure is positive
161 p(i) = max(p(i), 1.0e-12_rp)
162 ! E = p / (gamma - 1) + 0.5 * rho * (u^2 + v^2 + w^2)
163 e(i) = p(i) * inv_gamma_m1 + ruvw(i)
164 end do
165 !$omp end parallel do simd
166 end subroutine compressible_ops_cpu_update_e
167
168end module compressible_ops_cpu
CPU implementation of compressible flow operations.
subroutine, public compressible_ops_cpu_update_uvw(u, v, w, m_x, m_y, m_z, rho, n)
Update u,v,w fields.
subroutine, public compressible_ops_cpu_update_e(e, p, ruvw, gamma, n)
Update E field.
subroutine, public compressible_ops_cpu_compute_entropy(s, p, rho, gamma, n)
Compute entropy field S = 1/(gamma-1) * rho * (log(p) - gamma * log(rho)) on CPU.
subroutine, public compressible_ops_cpu_compute_max_wave_speed(max_wave_speed, u, v, w, gamma, p, rho, n)
Compute maximum wave speed for compressible flows on CPU.
subroutine, public compressible_ops_cpu_update_mxyz_p_ruvw(m_x, m_y, m_z, p, ruvw, u, v, w, e, rho, gamma, n)
Update m_x, m_y, m_z, p, ruvw, fields.
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
#define max(a, b)
Definition tensor.cu:40