Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
spalding_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
37 implicit none
38 private
39
40 public :: spalding_compute_cpu
41
42contains
46 subroutine spalding_compute_cpu(u, v, w, ind_r, ind_s, ind_t, ind_e, &
47 n_x, n_y, n_z, nu, h, tau_x, tau_y, tau_z, n_nodes, lx, nelv, &
48 kappa, B, tstep)
49 integer, intent(in) :: n_nodes, lx, nelv, tstep
50 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, v, w
51 integer, intent(in), dimension(n_nodes) :: ind_r, ind_s, ind_t, ind_e
52 real(kind=rp), dimension(n_nodes), intent(in) :: n_x, n_y, n_z, h, nu
53 real(kind=rp), dimension(n_nodes), intent(inout) :: tau_x, tau_y, tau_z
54 real(kind=rp), intent(in) :: kappa, b
55 integer :: i
56 real(kind=rp) :: ui, vi, wi, magu, utau, normu, guess
57
58 do i=1, n_nodes
59 ! Sample the velocity
60 ui = u(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
61 vi = v(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
62 wi = w(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
63
64 ! Project on tangential direction
65 normu = ui * n_x(i) + vi * n_y(i) + wi * n_z(i)
66
67 ui = ui - normu * n_x(i)
68 vi = vi - normu * n_y(i)
69 wi = wi - normu * n_z(i)
70
71 magu = sqrt(ui**2 + vi**2 + wi**2)
72
73 ! Get initial guess for Newton solver
74 if (tstep .eq. 1) then
75 guess = sqrt(magu * nu(i) / h(i))
76 else
77 guess = tau_x(i)**2 + tau_y(i)**2 + tau_z(i)**2
78 guess = sqrt(sqrt(guess))
79 end if
80
81 utau = solve_cpu(magu, h(i), guess, nu(i), kappa, b)
82
83 ! Distribute according to the velocity vector
84 tau_x(i) = -utau**2 * ui / magu
85 tau_y(i) = -utau**2 * vi / magu
86 tau_z(i) = -utau**2 * wi / magu
87 end do
88
89 end subroutine spalding_compute_cpu
90
98 function solve_cpu(u, y, guess, nu, kappa, B) result(utau)
99 real(kind=rp), intent(in) :: u
100 real(kind=rp), intent(in) :: y
101 real(kind=rp), intent(in) :: guess
102 real(kind=rp), intent(in) :: nu, kappa, b
103 real(kind=rp) :: yp, up, utau
104 real(kind=rp) :: error, f, df, old
105 integer :: niter, k, maxiter
106 character(len=LOG_SIZE) :: log_msg
107
108 utau = guess
109
110 maxiter = 100
111
112 do k=1, maxiter
113 up = u / utau
114 yp = y * utau / nu
115 niter = k
116 old = utau
117
118 ! Evaluate function and its derivative
119 f = (up + exp(-kappa*b)* &
120 (exp(kappa*up) - 1.0_rp - kappa*up - 0.5_rp*(kappa*up)**2 - &
121 1.0_rp/6*(kappa*up)**3) - yp)
122
123 df = (-y / nu - u/utau**2 - kappa*up/utau*exp(-kappa*b) * &
124 (exp(kappa*up) - 1 - kappa*up - 0.5*(kappa*up)**2))
125
126 ! Update solution
127 utau = utau - f / df
128
129 error = abs((old - utau)/old)
130
131 if (error < 1e-3) then
132 exit
133 endif
134
135 enddo
136
137 if (niter .eq. maxiter) then
138 write(log_msg, *) "Newton not converged", error, f, utau, old, guess
139 call neko_log%message(log_msg, neko_log_debug)
140 end if
141 end function solve_cpu
142end module spalding_cpu
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:78
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
integer, parameter, public log_size
Definition log.f90:46
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements the CPU kernel for the spalding_t type.
real(kind=rp) function solve_cpu(u, y, guess, nu, kappa, b)
Newton solver for the algebraic equation defined by the law on cpu.
subroutine, public spalding_compute_cpu(u, v, w, ind_r, ind_s, ind_t, ind_e, n_x, n_y, n_z, nu, h, tau_x, tau_y, tau_z, n_nodes, lx, nelv, kappa, b, tstep)
Compute the wall shear stress on cpu using Spalding's model.