Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
vreman_device.f90
Go to the documentation of this file.
1! Copyright (c) 2023-2026, 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 math, only : neko_eps
38 use registry, only : neko_registry
39 use field, only : field_t
40 use operators, only : dudxyz
41 use coefs, only : coef_t
42 use gs_ops, only : gs_op_add
43 use device_math, only : device_col2
46
47 implicit none
48 private
49
50 public :: vreman_compute_device
51
52contains
53
67 subroutine vreman_compute_device(if_ext, t, tstep, coef, nut, delta, c,&
68 if_corr, scalar_name, ri_c, ref_temp, g)
69
70 logical, intent(in) :: if_ext
71 real(kind=rp), intent(in) :: t
72 integer, intent(in) :: tstep
73 type(coef_t), intent(in) :: coef
74 type(field_t), intent(inout) :: nut
75 type(field_t), intent(in) :: delta
76 real(kind=rp), intent(in) :: c
77 logical, intent(in) :: if_corr
78 character(len=*), intent(in) :: scalar_name
79 real(kind=rp), intent(in) :: ri_c, ref_temp
80 real(kind=rp), intent(in) :: g(3)
81
82 ! This is the alpha tensor in the paper
83 type(field_t), pointer :: a11, a12, a13, a21, a22, a23, a31, a32, a33
84 type(field_t), pointer :: u, v, w
85 type(field_t), pointer :: temperature, dtdx, dtdy, dtdz
86
87 type(field_t), pointer :: beta11
88 type(field_t), pointer :: beta12
89 type(field_t), pointer :: beta13
90 type(field_t), pointer :: beta22
91 type(field_t), pointer :: beta23
92 type(field_t), pointer :: beta33
93 type(field_t), pointer :: b_beta
94 type(field_t), pointer :: aijaij
95 integer :: temp_indices(9), temp_indices_buoy(3)
96 integer :: e, i
97
98 if (if_ext .eqv. .true.) then
99 u => neko_registry%get_field_by_name("u_e")
100 v => neko_registry%get_field_by_name("v_e")
101 w => neko_registry%get_field_by_name("w_e")
102 else
103 u => neko_registry%get_field_by_name("u")
104 v => neko_registry%get_field_by_name("v")
105 w => neko_registry%get_field_by_name("w")
106 end if
107
108 call neko_scratch_registry%request_field(a11, temp_indices(1), .false.)
109 call neko_scratch_registry%request_field(a12, temp_indices(2), .false.)
110 call neko_scratch_registry%request_field(a13, temp_indices(3), .false.)
111 call neko_scratch_registry%request_field(a21, temp_indices(4), .false.)
112 call neko_scratch_registry%request_field(a22, temp_indices(5), .false.)
113 call neko_scratch_registry%request_field(a23, temp_indices(6), .false.)
114 call neko_scratch_registry%request_field(a31, temp_indices(7), .false.)
115 call neko_scratch_registry%request_field(a32, temp_indices(8), .false.)
116 call neko_scratch_registry%request_field(a33, temp_indices(9), .false.)
117
118 ! Compute the derivatives of the velocity (the alpha tensor)
119 call dudxyz (a11%x, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
120 call dudxyz (a12%x, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
121 call dudxyz (a13%x, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
122
123 call dudxyz (a21%x, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
124 call dudxyz (a22%x, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
125 call dudxyz (a23%x, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
126
127 call dudxyz (a31%x, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
128 call dudxyz (a32%x, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
129 call dudxyz (a33%x, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
130
131 call coef%gs_h%op(a11, gs_op_add)
132 call coef%gs_h%op(a12, gs_op_add)
133 call coef%gs_h%op(a13, gs_op_add)
134 call coef%gs_h%op(a21, gs_op_add)
135 call coef%gs_h%op(a22, gs_op_add)
136 call coef%gs_h%op(a23, gs_op_add)
137 call coef%gs_h%op(a31, gs_op_add)
138 call coef%gs_h%op(a32, gs_op_add)
139 call coef%gs_h%op(a33, gs_op_add)
140
141 if (if_corr) then
142 temperature => neko_registry%get_field_by_name(scalar_name)
143
144 call neko_scratch_registry%request_field(dtdx, temp_indices_buoy(1), .false.)
145 call neko_scratch_registry%request_field(dtdy, temp_indices_buoy(2), .false.)
146 call neko_scratch_registry%request_field(dtdz, temp_indices_buoy(3), .false.)
147
148 call dudxyz(dtdx%x, temperature%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
149 call dudxyz(dtdy%x, temperature%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
150 call dudxyz(dtdz%x, temperature%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
151
152 call device_vreman_nut_compute_buoy(a11%x_d, a12%x_d, a13%x_d, &
153 a21%x_d, a22%x_d, a23%x_d, &
154 a31%x_d, a32%x_d, a33%x_d, &
155 delta%x_d, nut%x_d, coef%mult_d, &
156 c, neko_eps, a11%dof%size(), &
157 dtdx%x_d, dtdy%x_d, dtdz%x_d, &
158 g, ri_c, ref_temp)
159 call neko_scratch_registry%relinquish_field(temp_indices_buoy)
160 else
161 call device_vreman_nut_compute(a11%x_d, a12%x_d, a13%x_d, &
162 a21%x_d, a22%x_d, a23%x_d, &
163 a31%x_d, a32%x_d, a33%x_d, &
164 delta%x_d, nut%x_d, coef%mult_d, &
165 c, neko_eps, a11%dof%size())
166 end if
167
168 call coef%gs_h%op(nut, gs_op_add)
169 call device_col2(nut%x_d, coef%mult_d, nut%dof%size())
170
171 call neko_scratch_registry%relinquish_field(temp_indices)
172 end subroutine vreman_compute_device
173
174end module vreman_device
Coefficients.
Definition coef.f90:34
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_vreman_nut_compute(a11_d, a12_d, a13_d, a21_d, a22_d, a23_d, a31_d, a32_d, a33_d, delta_d, nut_d, mult_d, c, eps, n)
Compute the eddy viscosity field for the Sigma model indevice.
subroutine, public device_vreman_nut_compute_buoy(a11_d, a12_d, a13_d, a21_d, a22_d, a23_d, a31_d, a32_d, a33_d, delta_d, nut_d, mult_d, c, eps, n, dtdx_d, dtdy_d, dtdz_d, g, ri_c, ref_temp)
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
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:92
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Implements the device kernel for the vreman_t type.
subroutine, public vreman_compute_device(if_ext, t, tstep, coef, nut, delta, c, if_corr, scalar_name, ri_c, ref_temp, g)
Compute eddy viscosity on the device.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:58