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