Neko  0.8.1
A portable framework for high-order spectral element flow simulations
vreman_cpu.f90
Go to the documentation of this file.
1 ! Copyright (c) 2023, 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 !
34 module vreman_cpu
35  use num_types, only : rp
36  use field_list, only : field_list_t
37  use math, only : cadd, neko_eps
40  use field, only : field_t
41  use operators, only : dudxyz
42  use coefs, only : coef_t
43  implicit none
44  private
45 
46  public :: vreman_compute_cpu
47 
48 contains
49 
57  subroutine vreman_compute_cpu(t, tstep, coef, nut, delta, c)
58  real(kind=rp), intent(in) :: t
59  integer, intent(in) :: tstep
60  type(coef_t), intent(in) :: coef
61  type(field_t), intent(inout) :: nut
62  type(field_t), intent(in) :: delta
63  real(kind=rp), intent(in) :: c
64  ! This is the alpha tensor in the paper
65  type(field_t), pointer :: a11, a12, a13, a21, a22, a23, a31, a32, a33
66  type(field_t), pointer :: u, v, w
67 
68  real(kind=rp) :: beta11
69  real(kind=rp) :: beta12
70  real(kind=rp) :: beta13
71  real(kind=rp) :: beta22
72  real(kind=rp) :: beta23
73  real(kind=rp) :: beta33
74  real(kind=rp) :: b_beta
75  real(kind=rp) :: aijaij
76  integer :: temp_indices(9)
77  integer :: e, i
78 
79  u => neko_field_registry%get_field_by_name("u")
80  v => neko_field_registry%get_field_by_name("v")
81  w => neko_field_registry%get_field_by_name("u")
82 
83  call neko_scratch_registry%request_field(a11, temp_indices(1))
84  call neko_scratch_registry%request_field(a12, temp_indices(2))
85  call neko_scratch_registry%request_field(a13, temp_indices(3))
86  call neko_scratch_registry%request_field(a21, temp_indices(4))
87  call neko_scratch_registry%request_field(a22, temp_indices(5))
88  call neko_scratch_registry%request_field(a23, temp_indices(6))
89  call neko_scratch_registry%request_field(a31, temp_indices(7))
90  call neko_scratch_registry%request_field(a32, temp_indices(8))
91  call neko_scratch_registry%request_field(a33, temp_indices(9))
92 
93 
94  ! Compute the derivatives of the velocity (the alpha tensor)
95  call dudxyz (a11%x, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
96  call dudxyz (a12%x, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
97  call dudxyz (a13%x, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
98 
99  call dudxyz (a21%x, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
100  call dudxyz (a22%x, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
101  call dudxyz (a23%x, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
102 
103  call dudxyz (a31%x, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
104  call dudxyz (a32%x, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
105  call dudxyz (a33%x, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
106 
107  do e=1, coef%msh%nelv
108  do i=1, coef%Xh%lxyz
109  ! beta_ij = alpha_mi alpha_mj
110  beta11 = a11%x(i,1,1,e)**2 + a21%x(i,1,1,e)**2 + a31%x(i,1,1,e)**2
111  beta22 = a12%x(i,1,1,e)**2 + a22%x(i,1,1,e)**2 + a32%x(i,1,1,e)**2
112  beta33 = a13%x(i,1,1,e)**2 + a23%x(i,1,1,e)**2 + a33%x(i,1,1,e)**2
113  beta12 = a11%x(i,1,1,e)*a12%x(i,1,1,e) + &
114  a21%x(i,1,1,e)*a22%x(i,1,1,e) + &
115  a31%x(i,1,1,e)*a32%x(i,1,1,e)
116  beta13 = a11%x(i,1,1,e)*a13%x(i,1,1,e) + &
117  a21%x(i,1,1,e)*a23%x(i,1,1,e) + &
118  a31%x(i,1,1,e)*a33%x(i,1,1,e)
119  beta23 = a12%x(i,1,1,e)*a13%x(i,1,1,e) + &
120  a22%x(i,1,1,e)*a23%x(i,1,1,e) + &
121  a32%x(i,1,1,e)*a33%x(i,1,1,e)
122 
123  b_beta = beta11*beta22 - beta12*beta12 + beta11*beta33 - beta13*beta13 &
124  + beta22*beta33 - beta23*beta23
125 
126  b_beta = max(0.0_rp, b_beta)
127 
128  ! alpha_ij alpha_ij
129  aijaij = beta11 + beta22 + beta33
130 
131  nut%x(i,1,1,e) = c*delta%x(i,1,1,e) * sqrt(b_beta/(aijaij + neko_eps))
132  end do
133  end do
134 
135  call neko_scratch_registry%relinquish_field(temp_indices)
136  end subroutine vreman_compute_cpu
137 
138 end module vreman_cpu
139 
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
Definition: math.f90:60
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition: math.f90:258
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition: math.f90:67
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:70
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.
Definition: vreman_cpu.f90:34
subroutine, public vreman_compute_cpu(t, tstep, coef, nut, delta, c)
Compute eddy viscosity on the CPU.
Definition: vreman_cpu.f90:58
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
field_list_t, To be able to group fields together
Definition: field_list.f90:7
#define max(a, b)
Definition: tensor.cu:40