Neko  0.9.99
A portable framework for high-order spectral element flow simulations
smagorinsky_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
39  use field, only : field_t
40  use operators, only : strain_rate
41  use coefs, only : coef_t
42  use gs_ops, only : gs_op_add
43  implicit none
44  private
45 
46  public :: smagorinsky_compute_cpu
47 
48 contains
49 
57  subroutine smagorinsky_compute_cpu(t, tstep, coef, nut, delta, c_s)
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_s
64  type(field_t), pointer :: u, v, w
65  ! double of the strain rate tensor
66  type(field_t), pointer :: s11, s22, s33, s12, s13, s23
67  real(kind=rp) :: s_abs
68  integer :: temp_indices(6)
69  integer :: e, i
70 
71  u => neko_field_registry%get_field_by_name("u")
72  v => neko_field_registry%get_field_by_name("v")
73  w => neko_field_registry%get_field_by_name("w")
74 
75  call neko_scratch_registry%request_field(s11, temp_indices(1))
76  call neko_scratch_registry%request_field(s22, temp_indices(2))
77  call neko_scratch_registry%request_field(s33, temp_indices(3))
78  call neko_scratch_registry%request_field(s12, temp_indices(4))
79  call neko_scratch_registry%request_field(s13, temp_indices(5))
80  call neko_scratch_registry%request_field(s23, temp_indices(6))
81 
82  ! Compute the strain rate tensor
83  call strain_rate(s11%x, s22%x, s33%x, s12%x, s13%x, s23%x, u, v, w, coef)
84 
85  call coef%gs_h%op(s11, gs_op_add)
86  call coef%gs_h%op(s22, gs_op_add)
87  call coef%gs_h%op(s33, gs_op_add)
88  call coef%gs_h%op(s12, gs_op_add)
89  call coef%gs_h%op(s13, gs_op_add)
90  call coef%gs_h%op(s23, gs_op_add)
91 
92  do concurrent(i = 1:s11%dof%size())
93  s11%x(i,1,1,1) = s11%x(i,1,1,1) * coef%mult(i,1,1,1)
94  s22%x(i,1,1,1) = s22%x(i,1,1,1) * coef%mult(i,1,1,1)
95  s33%x(i,1,1,1) = s33%x(i,1,1,1) * coef%mult(i,1,1,1)
96  s12%x(i,1,1,1) = s12%x(i,1,1,1) * coef%mult(i,1,1,1)
97  s13%x(i,1,1,1) = s13%x(i,1,1,1) * coef%mult(i,1,1,1)
98  s23%x(i,1,1,1) = s23%x(i,1,1,1) * coef%mult(i,1,1,1)
99  end do
100 
101  do concurrent(e = 1:coef%msh%nelv)
102  do concurrent(i = 1:coef%Xh%lxyz)
103  s_abs = sqrt(2.0_rp * (s11%x(i,1,1,e)*s11%x(i,1,1,e) + &
104  s22%x(i,1,1,e)*s22%x(i,1,1,e) + &
105  s33%x(i,1,1,e)*s33%x(i,1,1,e)) + &
106  4.0_rp * (s12%x(i,1,1,e)*s12%x(i,1,1,e) + &
107  s13%x(i,1,1,e)*s13%x(i,1,1,e) + &
108  s23%x(i,1,1,e)*s23%x(i,1,1,e)))
109 
110  nut%x(i,1,1,e) = c_s**2 * delta%x(i,1,1,e)**2 * s_abs
111  end do
112  end do
113 
114  call neko_scratch_registry%relinquish_field(temp_indices)
115  end subroutine smagorinsky_compute_cpu
116 
117 end module smagorinsky_cpu
118 
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 strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
Definition: operators.f90:430
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 smagorinsky_t type.
subroutine, public smagorinsky_compute_cpu(t, tstep, coef, nut, delta, c_s)
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
Definition: field_list.f90:13