Neko 1.99.4
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
entropy_viscosity_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
44
45contains
46
56 subroutine entropy_viscosity_compute_residual_cpu(entropy_residual, &
57 S, S_lag1, S_lag2, S_lag3, bdf_coeffs, dt, n)
58 integer, intent(in) :: n
59 real(kind=rp), dimension(n), intent(out) :: entropy_residual
60 real(kind=rp), dimension(n), intent(in) :: s, s_lag1, s_lag2, s_lag3
61 real(kind=rp), intent(in) :: bdf_coeffs(4)
62 real(kind=rp), intent(in) :: dt
63 integer :: i
64
65 !OCL NORECURRENCE, NOVREC, NOALIAS
66 !DIR$ CONCURRENT
67 !GCC$ ivdep
68 !$omp parallel do simd
69 do i = 1, n
70 entropy_residual(i) = (bdf_coeffs(1) * s(i) &
71 - bdf_coeffs(2) * s_lag1(i) &
72 - bdf_coeffs(3) * s_lag2(i) &
73 - bdf_coeffs(4) * s_lag3(i)) / dt
74 end do
75 !$omp end parallel do simd
76
78
87 entropy_residual, h, c_avisc_entropy, n_S, n)
88 integer, intent(in) :: n
89 real(kind=rp), dimension(n), intent(out) :: reg_coeff
90 real(kind=rp), dimension(n), intent(in) :: entropy_residual, h
91 real(kind=rp), intent(in) :: c_avisc_entropy, n_s
92 integer :: i
93
94 !OCL NORECURRENCE, NOVREC, NOALIAS
95 !DIR$ CONCURRENT
96 !GCC$ ivdep
97 !$omp parallel do simd
98 do i = 1, n
99 reg_coeff(i) = c_avisc_entropy * h(i) * h(i) * &
100 entropy_residual(i) / n_s
101 end do
102 !$omp end parallel do simd
103
105
110 subroutine entropy_viscosity_apply_element_max_cpu(reg_coeff, lx, nelv)
111 integer, intent(in) :: lx, nelv
112 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: reg_coeff
113 integer :: i, j, k, el
114 real(kind=rp) :: max_visc_el
115
116 !$omp parallel do private(i, j, k, max_visc_el)
117 do el = 1, nelv
118 max_visc_el = 0.0_rp
119
120 do k = 1, lx
121 do j = 1, lx
122 !OCL NORECURRENCE, NOVREC, NOALIAS
123 !DIR$ CONCURRENT
124 !GCC$ ivdep
125 !$omp simd reduction(max:max_visc_el)
126 do i = 1, lx
127 max_visc_el = max(max_visc_el, reg_coeff(i, j, k, el))
128 end do
129 end do
130 end do
131
132 do k = 1, lx
133 do j = 1, lx
134 !OCL NORECURRENCE, NOVREC, NOALIAS
135 !DIR$ CONCURRENT
136 !GCC$ ivdep
137 !$omp simd
138 do i = 1, lx
139 reg_coeff(i, j, k, el) = max_visc_el
140 end do
141 end do
142 end do
143 end do
144 !$omp end parallel do
145
147
155 h, max_wave_speed, c_avisc_low, n)
156 integer, intent(in) :: n
157 real(kind=rp), dimension(n), intent(inout) :: reg_coeff
158 real(kind=rp), dimension(n), intent(in) :: h, max_wave_speed
159 real(kind=rp), intent(in) :: c_avisc_low
160 integer :: i
161 real(kind=rp) :: low_order_visc
162
163 !OCL NORECURRENCE, NOVREC, NOALIAS
164 !DIR$ CONCURRENT
165 !GCC$ ivdep
166 !$omp parallel do simd private(low_order_visc)
167 do i = 1, n
168 low_order_visc = c_avisc_low * h(i) * max_wave_speed(i)
169 reg_coeff(i) = min(reg_coeff(i), low_order_visc)
170 end do
171 !$omp end parallel do simd
172
174
181 temp_field, mult_field, n)
182 integer, intent(in) :: n
183 real(kind=rp), dimension(n), intent(out) :: reg_coeff
184 real(kind=rp), dimension(n), intent(in) :: temp_field, mult_field
185 integer :: i
186
187 !OCL NORECURRENCE, NOVREC, NOALIAS
188 !DIR$ CONCURRENT
189 !GCC$ ivdep
190 !$omp parallel do simd
191 do i = 1, n
192 reg_coeff(i) = temp_field(i) / mult_field(i)
193 end do
194 !$omp end parallel do simd
195
197
198end module entropy_viscosity_cpu
CPU backend for entropy viscosity regularization.
subroutine, public entropy_viscosity_compute_viscosity_cpu(reg_coeff, entropy_residual, h, c_avisc_entropy, n_s, n)
Compute viscosity from entropy residual on CPU.
subroutine, public entropy_viscosity_compute_residual_cpu(entropy_residual, s, s_lag1, s_lag2, s_lag3, bdf_coeffs, dt, n)
Compute entropy residual on CPU.
subroutine, public entropy_viscosity_smooth_divide_cpu(reg_coeff, temp_field, mult_field, n)
Divide by multiplicity for smoothing on CPU.
subroutine, public entropy_viscosity_clamp_to_low_order_cpu(reg_coeff, h, max_wave_speed, c_avisc_low, n)
Clamp regularization coefficient to low-order viscosity on CPU.
subroutine, public entropy_viscosity_apply_element_max_cpu(reg_coeff, lx, nelv)
Apply element-wise maximum on CPU.
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
#define max(a, b)
Definition tensor.cu:40