Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
vreman_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!
35 use num_types, only : rp
36 use field_list, only : field_list_t
37 use math, only : cadd, neko_eps, col2
40 use field, only : field_t
41 use operators, only : dudxyz
42 use coefs, only : coef_t
43 use gs_ops, only : gs_op_add
44 implicit none
45 private
46
47 public :: vreman_compute_cpu
48
49contains
50
59 subroutine vreman_compute_cpu(if_ext, t, tstep, coef, nut, delta, c)
60 logical, intent(in) :: if_ext
61 real(kind=rp), intent(in) :: t
62 integer, intent(in) :: tstep
63 type(coef_t), intent(in) :: coef
64 type(field_t), intent(inout) :: nut
65 type(field_t), intent(in) :: delta
66 real(kind=rp), intent(in) :: c
67 ! This is the alpha tensor in the paper
68 type(field_t), pointer :: a11, a12, a13, a21, a22, a23, a31, a32, a33
69 type(field_t), pointer :: u, v, w
70
71 real(kind=rp) :: beta11
72 real(kind=rp) :: beta12
73 real(kind=rp) :: beta13
74 real(kind=rp) :: beta22
75 real(kind=rp) :: beta23
76 real(kind=rp) :: beta33
77 real(kind=rp) :: b_beta
78 real(kind=rp) :: aijaij
79 integer :: temp_indices(9)
80 integer :: e, i
81
82 if (if_ext .eqv. .true.) then
83 u => neko_field_registry%get_field_by_name("u_e")
84 v => neko_field_registry%get_field_by_name("v_e")
85 w => neko_field_registry%get_field_by_name("w_e")
86 else
87 u => neko_field_registry%get_field_by_name("u")
88 v => neko_field_registry%get_field_by_name("v")
89 w => neko_field_registry%get_field_by_name("w")
90 end if
91
92 call neko_scratch_registry%request_field(a11, temp_indices(1))
93 call neko_scratch_registry%request_field(a12, temp_indices(2))
94 call neko_scratch_registry%request_field(a13, temp_indices(3))
95 call neko_scratch_registry%request_field(a21, temp_indices(4))
96 call neko_scratch_registry%request_field(a22, temp_indices(5))
97 call neko_scratch_registry%request_field(a23, temp_indices(6))
98 call neko_scratch_registry%request_field(a31, temp_indices(7))
99 call neko_scratch_registry%request_field(a32, temp_indices(8))
100 call neko_scratch_registry%request_field(a33, temp_indices(9))
101
102
103 ! Compute the derivatives of the velocity (the alpha tensor)
104 call dudxyz (a11%x, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
105 call dudxyz (a12%x, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
106 call dudxyz (a13%x, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
107
108 call dudxyz (a21%x, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
109 call dudxyz (a22%x, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
110 call dudxyz (a23%x, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
111
112 call dudxyz (a31%x, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
113 call dudxyz (a32%x, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
114 call dudxyz (a33%x, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
115
116 call coef%gs_h%op(a11, gs_op_add)
117 call coef%gs_h%op(a12, gs_op_add)
118 call coef%gs_h%op(a13, gs_op_add)
119 call coef%gs_h%op(a21, gs_op_add)
120 call coef%gs_h%op(a22, gs_op_add)
121 call coef%gs_h%op(a23, gs_op_add)
122 call coef%gs_h%op(a31, gs_op_add)
123 call coef%gs_h%op(a32, gs_op_add)
124 call coef%gs_h%op(a33, gs_op_add)
125
126 do concurrent(e = 1:coef%msh%nelv)
127 do concurrent(i = 1:coef%Xh%lxyz)
128 ! beta_ij = alpha_mi alpha_mj
129 beta11 = a11%x(i,1,1,e)**2 + a21%x(i,1,1,e)**2 + a31%x(i,1,1,e)**2
130 beta22 = a12%x(i,1,1,e)**2 + a22%x(i,1,1,e)**2 + a32%x(i,1,1,e)**2
131 beta33 = a13%x(i,1,1,e)**2 + a23%x(i,1,1,e)**2 + a33%x(i,1,1,e)**2
132 beta12 = a11%x(i,1,1,e)*a12%x(i,1,1,e) + &
133 a21%x(i,1,1,e)*a22%x(i,1,1,e) + &
134 a31%x(i,1,1,e)*a32%x(i,1,1,e)
135 beta13 = a11%x(i,1,1,e)*a13%x(i,1,1,e) + &
136 a21%x(i,1,1,e)*a23%x(i,1,1,e) + &
137 a31%x(i,1,1,e)*a33%x(i,1,1,e)
138 beta23 = a12%x(i,1,1,e)*a13%x(i,1,1,e) + &
139 a22%x(i,1,1,e)*a23%x(i,1,1,e) + &
140 a32%x(i,1,1,e)*a33%x(i,1,1,e)
141
142 b_beta = beta11*beta22 - beta12*beta12 + beta11*beta33 &
143 - beta13*beta13 + beta22*beta33 - beta23*beta23
144
145 b_beta = max(0.0_rp, b_beta)
146
147 ! alpha_ij alpha_ij
148 aijaij = beta11 + beta22 + beta33
149
150 nut%x(i,1,1,e) = c*delta%x(i,1,1,e)*delta%x(i,1,1,e) &
151 * sqrt(b_beta/(aijaij + neko_eps)) &
152 * coef%mult(i,1,1,e)
153 end do
154 end do
155
156 call coef%gs_h%op(nut, gs_op_add)
157 call col2(nut%x, coef%mult, nut%dof%size())
158
159 call neko_scratch_registry%relinquish_field(temp_indices)
160 end subroutine vreman_compute_cpu
161
162end module vreman_cpu
163
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
Definition math.f90:60
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:322
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:728
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
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:78
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 vreman_t type.
subroutine, public vreman_compute_cpu(if_ext, t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the CPU.
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
#define max(a, b)
Definition tensor.cu:40