Neko  0.9.0
A portable framework for high-order spectral element flow simulations
sigma_cpu.f90
Go to the documentation of this file.
1 ! Copyright (c) 2023-2024, 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 !
37 
38 module sigma_cpu
39  use num_types, only : rp
40  use field_list, only : field_list_t
43  use field, only : field_t
44  use operators, only : dudxyz
45  use coefs, only : coef_t
46  use gs_ops, only : gs_op_add
47  implicit none
48  private
49 
50  public :: sigma_compute_cpu
51 
52 contains
53 
61  subroutine sigma_compute_cpu(t, tstep, coef, nut, delta, c)
62  real(kind=rp), intent(in) :: t
63  integer, intent(in) :: tstep
64  type(coef_t), intent(in) :: coef
65  type(field_t), intent(inout) :: nut
66  type(field_t), intent(in) :: delta
67  real(kind=rp), intent(in) :: c
68  ! This is the velocity gradient tensor
69  type(field_t), pointer :: g11, g12, g13, g21, g22, g23, g31, g32, g33
70  type(field_t), pointer :: u, v, w
71 
72  real(kind=rp) :: sigg11, sigg12, sigg13, sigg22, sigg23, sigg33
73  real(kind=rp) :: sigma1, sigma2, sigma3
74  real(kind=rp) :: invariant1, invariant2, invariant3
75  real(kind=rp) :: alpha1, alpha2, alpha3
76  real(kind=rp) :: dsigma
77  real(kind=rp) :: pi_3 = 4.0_rp/3.0_rp*atan(1.0_rp)
78  real(kind=rp) :: tmp1
79  real(kind=rp) :: eps
80 
81  integer :: temp_indices(9)
82  integer :: e, i
83 
84  ! some constant
85  eps = 1.d-14
86 
87 
88  ! get fields from registry
89  u => neko_field_registry%get_field_by_name("u")
90  v => neko_field_registry%get_field_by_name("v")
91  w => neko_field_registry%get_field_by_name("w")
92 
93  call neko_scratch_registry%request_field(g11, temp_indices(1))
94  call neko_scratch_registry%request_field(g12, temp_indices(2))
95  call neko_scratch_registry%request_field(g13, temp_indices(3))
96  call neko_scratch_registry%request_field(g21, temp_indices(4))
97  call neko_scratch_registry%request_field(g22, temp_indices(5))
98  call neko_scratch_registry%request_field(g23, temp_indices(6))
99  call neko_scratch_registry%request_field(g31, temp_indices(7))
100  call neko_scratch_registry%request_field(g32, temp_indices(8))
101  call neko_scratch_registry%request_field(g33, temp_indices(9))
102 
103 
104  ! Compute the derivatives of the velocity (the components of the g tensor)
105  call dudxyz (g11%x, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
106  call dudxyz (g12%x, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
107  call dudxyz (g13%x, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
108 
109  call dudxyz (g21%x, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
110  call dudxyz (g22%x, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
111  call dudxyz (g23%x, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
112 
113  call dudxyz (g31%x, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
114  call dudxyz (g32%x, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
115  call dudxyz (g33%x, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
116 
117  call coef%gs_h%op(g11, gs_op_add)
118  call coef%gs_h%op(g12, gs_op_add)
119  call coef%gs_h%op(g13, gs_op_add)
120  call coef%gs_h%op(g21, gs_op_add)
121  call coef%gs_h%op(g22, gs_op_add)
122  call coef%gs_h%op(g23, gs_op_add)
123  call coef%gs_h%op(g31, gs_op_add)
124  call coef%gs_h%op(g32, gs_op_add)
125  call coef%gs_h%op(g33, gs_op_add)
126 
127  do concurrent(i = 1:g11%dof%size())
128  g11%x(i,1,1,1) = g11%x(i,1,1,1) * coef%mult(i,1,1,1)
129  g12%x(i,1,1,1) = g12%x(i,1,1,1) * coef%mult(i,1,1,1)
130  g13%x(i,1,1,1) = g13%x(i,1,1,1) * coef%mult(i,1,1,1)
131  g21%x(i,1,1,1) = g21%x(i,1,1,1) * coef%mult(i,1,1,1)
132  g22%x(i,1,1,1) = g22%x(i,1,1,1) * coef%mult(i,1,1,1)
133  g23%x(i,1,1,1) = g23%x(i,1,1,1) * coef%mult(i,1,1,1)
134  g31%x(i,1,1,1) = g31%x(i,1,1,1) * coef%mult(i,1,1,1)
135  g32%x(i,1,1,1) = g32%x(i,1,1,1) * coef%mult(i,1,1,1)
136  g33%x(i,1,1,1) = g33%x(i,1,1,1) * coef%mult(i,1,1,1)
137  end do
138 
139  do concurrent(e = 1:coef%msh%nelv)
140  do concurrent(i = 1:coef%Xh%lxyz)
141  ! G_ij = g^t g = g_mi g_mj
142  sigg11 = g11%x(i,1,1,e)**2 + g21%x(i,1,1,e)**2 + g31%x(i,1,1,e)**2
143  sigg22 = g12%x(i,1,1,e)**2 + g22%x(i,1,1,e)**2 + g32%x(i,1,1,e)**2
144  sigg33 = g13%x(i,1,1,e)**2 + g23%x(i,1,1,e)**2 + g33%x(i,1,1,e)**2
145  sigg12 = g11%x(i,1,1,e)*g12%x(i,1,1,e) + &
146  g21%x(i,1,1,e)*g22%x(i,1,1,e) + &
147  g31%x(i,1,1,e)*g32%x(i,1,1,e)
148  sigg13 = g11%x(i,1,1,e)*g13%x(i,1,1,e) + &
149  g21%x(i,1,1,e)*g23%x(i,1,1,e) + &
150  g31%x(i,1,1,e)*g33%x(i,1,1,e)
151  sigg23 = g12%x(i,1,1,e)*g13%x(i,1,1,e) + &
152  g22%x(i,1,1,e)*g23%x(i,1,1,e) + &
153  g32%x(i,1,1,e)*g33%x(i,1,1,e)
154 
155  ! If LAPACK compute eigenvalues of the semi-definite positive matrix G
156  ! ..........to be done later on......
157  ! ELSE use the analytical method as done in the following
158 
159  ! eigenvalues with the analytical method of Hasan et al. (2001)
160  ! doi:10.1006/jmre.2001.2400
161  if (abs(sigg11) .lt. eps) then
162  sigg11 = 0.0_rp
163  end if
164  if (abs(sigg12) .lt. eps) then
165  sigg12 = 0.0_rp
166  end if
167  if (abs(sigg13) .lt. eps) then
168  sigg13 = 0.0_rp
169  end if
170  if (abs(sigg22) .lt. eps) then
171  sigg22 = 0.0_rp
172  end if
173  if (abs(sigg23) .lt. eps) then
174  sigg23 = 0.0_rp
175  end if
176  if (abs(sigg33) .lt. eps) then
177  sigg33 = 0.0_rp
178  end if
179 
180  if (abs(sigg12*sigg12 + &
181  sigg13*sigg13 + sigg23*sigg23) .lt. eps) then
182  ! G is diagonal
183  ! estimate the singular values according to:
184  sigma1 = sqrt(max(max(max(sigg11, sigg22), sigg33), 0.0_rp))
185  sigma3 = sqrt(max(min(min(sigg11, sigg22), sigg33), 0.0_rp))
186  invariant1 = sigg11 + sigg22 + sigg33
187  sigma2 = sqrt(abs(invariant1 - sigma1*sigma1 - sigma3*sigma3))
188  else
189 
190  ! estimation of invariants
191  invariant1 = sigg11 + sigg22 + sigg33
192  invariant2 = sigg11*sigg22 + sigg11*sigg33 + sigg22*sigg33 - &
193  (sigg12*sigg12 + sigg13*sigg13 + sigg23*sigg23)
194  invariant3 = sigg11*sigg22*sigg33 + &
195  2.0_rp*sigg12*sigg13*sigg23 - &
196  (sigg11*sigg23*sigg23 + sigg22*sigg13*sigg13 + &
197  sigg33*sigg12*sigg12)
198 
199  ! G is symmetric semi-definite positive matrix:
200  ! the invariants have to be larger-equal zero
201  ! which is obtained via forcing
202  invariant1 = max(invariant1, 0.0_rp)
203  invariant2 = max(invariant2, 0.0_rp)
204  invariant3 = max(invariant3, 0.0_rp)
205 
206  ! compute the following angles from the invariants
207  alpha1 = invariant1*invariant1/9.0_rp - invariant2/3.0_rp
208 
209  ! since alpha1 is always positive (see Hasan et al. (2001))
210  ! forcing is applied
211  alpha1 = max(alpha1, 0.0_rp)
212 
213  alpha2 = invariant1*invariant1*invariant1/27.0_rp - &
214  invariant1*invariant2/6.0_rp + invariant3/2.0_rp
215 
216  ! since acos(alpha2/(alpha1^(3/2)))/3.0_rp only valid for
217  ! alpha2^2 < alpha1^3.0_rp and arccos(x) only valid for -1<=x<=1
218  ! alpha3 is between 0 and pi/3
219  tmp1 = alpha2/sqrt(alpha1 * alpha1 * alpha1)
220 
221  if (tmp1 .le. -1.0_rp) then
222  ! alpha3=pi/3 -> cos(alpha3)=0.5
223  ! compute the singular values
224  sigma1 = sqrt(max(invariant1/3.0_rp + sqrt(alpha1), 0.0_rp))
225  sigma2 = sigma1
226  sigma3 = sqrt(invariant1/3.0_rp - 2.0_rp*sqrt(alpha1))
227 
228  elseif (tmp1 .ge. 1.0_rp) then
229  ! alpha3=0.0_rp -> cos(alpha3)=1.0
230  sigma1 = sqrt(max(invariant1/3.0_rp + 2.0_rp*sqrt(alpha1), &
231  0.0_rp))
232  sigma2 = sqrt(invariant1/3.0_rp - sqrt(alpha1))
233  sigma3 = sigma2
234  else
235  alpha3 = acos(tmp1)/3.0_rp
236 
237  if (abs(invariant3) .lt. eps) then
238  ! In case of Invariant3=0, one or more eigenvalues are equal to zero
239  ! Therefore force sigma3 to 0 and compute sigma1 and sigma2
240  sigma1 = sqrt(max(invariant1/3.0_rp + &
241  2.0_rp*sqrt(alpha1)*cos(alpha3), 0.0_rp))
242  sigma2 = sqrt(abs(invariant1 - sigma1*sigma1))
243  sigma3 = 0.0_rp
244  else
245  sigma1 = sqrt(max(invariant1/3.0_rp + &
246  2.0_rp*sqrt(alpha1)*cos(alpha3), 0.0_rp))
247  sigma2 = sqrt(invariant1/3.0_rp - &
248  2.0_rp*sqrt(alpha1)*cos(pi_3 + alpha3))
249  sigma3 = sqrt(abs(invariant1 - &
250  sigma1*sigma1-sigma2*sigma2))
251  end if ! Invariant3=0 ?
252  end if ! tmp1
253  end if ! G diagonal ?
254 
255  ! Estimate Dsigma
256  if (sigma1 .gt. 0.0_rp) then
257  dsigma = &
258  sigma3*(sigma1 - sigma2)*(sigma2 - sigma3)/(sigma1*sigma1)
259  else
260  dsigma = 0.0_rp
261  end if
262 
263  !clipping to avoid negative values
264  dsigma = max(dsigma, 0.0_rp)
265 
266  ! estimate turbulent viscosity
267 
268  nut%x(i,1,1,e) = (c*delta%x(i,1,1,e))**2 * dsigma
269 
270 
271  end do
272  end do
273 
274  call neko_scratch_registry%relinquish_field(temp_indices)
275  end subroutine sigma_compute_cpu
276 
277 end module sigma_cpu
278 
Coefficients.
Definition: coef.f90:34
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines a field.
Definition: field.f90:34
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 dudxyz(du, u, dr, ds, dt, coef)
Compute derivative of a scalar field along a single direction.
Definition: operators.f90:76
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Implements the CPU kernel for the sigma_t type. Following Nicoud et al. "Using singular values to bui...
Definition: sigma_cpu.f90:38
subroutine, public sigma_compute_cpu(t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the CPU.
Definition: sigma_cpu.f90:62
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
field_list_t, To be able to group fields together
Definition: field_list.f90:13
#define max(a, b)
Definition: tensor.cu:40