Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
dong_outflow.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2025, 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 neko_config
36 use dirichlet, only : dirichlet_t
38 use num_types, only : rp, c_rp
39 use bc, only : bc_t
40 use field, only : field_t
41 use dofmap, only : dofmap_t
42 use coefs, only : coef_t
43 use utils, only : nonlinear_index
45 use registry, only : neko_registry
46 use, intrinsic :: iso_c_binding, only : c_ptr, c_sizeof, c_null_ptr, &
47 c_associated
48 use json_module, only : json_file
50 use utils, only : neko_error
51 use time_state, only : time_state_t
52 implicit none
53 private
54
60 type, public, extends(bc_t) :: dong_outflow_t
61 type(field_t), pointer :: u
62 type(field_t), pointer :: v
63 type(field_t), pointer :: w
64 real(kind=rp) :: delta
65 real(kind=rp) :: uinf
66 type(c_ptr) :: normal_x_d = c_null_ptr
67 type(c_ptr) :: normal_y_d = c_null_ptr
68 type(c_ptr) :: normal_z_d = c_null_ptr
69 contains
70 procedure, pass(this) :: apply_scalar => dong_outflow_apply_scalar
71 procedure, pass(this) :: apply_vector => dong_outflow_apply_vector
72 procedure, pass(this) :: apply_scalar_dev => dong_outflow_apply_scalar_dev
73 procedure, pass(this) :: apply_vector_dev => dong_outflow_apply_vector_dev
75 procedure, pass(this) :: init => dong_outflow_init
77 procedure, pass(this) :: free => dong_outflow_free
79 procedure, pass(this) :: finalize => dong_outflow_finalize
80 end type dong_outflow_t
81
82contains
86 subroutine dong_outflow_init(this, coef, json)
87 class(dong_outflow_t), target, intent(inout) :: this
88 type(coef_t), target, intent(in) :: coef
89 type(json_file), intent(inout) :: json
90 call this%free()
91 call this%init_base(coef)
92
93 call json_get_or_default(json, 'delta', this%delta, 0.01_rp)
94 call json_get_or_default(json, 'velocity_scale', this%uinf, 1.0_rp)
95
96 end subroutine dong_outflow_init
97
100 subroutine dong_outflow_apply_scalar(this, x, n, time, strong)
101 class(dong_outflow_t), intent(inout) :: this
102 integer, intent(in) :: n
103 real(kind=rp), intent(inout), dimension(n) :: x
104 type(time_state_t), intent(in), optional :: time
105 logical, intent(in), optional :: strong
106 integer :: i, m, k, facet, idx(4)
107 real(kind=rp) :: vn, s0, ux, uy, uz, normal_xyz(3)
108 logical :: strong_
109
110 if (present(strong)) then
111 strong_ = strong
112 else
113 strong_ = .true.
114 end if
115
116 !Im actually not sure what to do if one has two dong that share a corner.
117 if (strong_) then
118 m = this%msk(0)
119 do i = 1, m
120 k = this%msk(i)
121 facet = this%facet(i)
122 ux = this%u%x(k,1,1,1)
123 uy = this%v%x(k,1,1,1)
124 uz = this%w%x(k,1,1,1)
125 idx = nonlinear_index(k, this%Xh%lx, this%Xh%lx, this%Xh%lx)
126 normal_xyz = this%coef%get_normal(idx(1), idx(2), idx(3), idx(4), &
127 facet)
128 vn = ux*normal_xyz(1) + uy*normal_xyz(2) + uz*normal_xyz(3)
129 s0 = 0.5_rp*(1.0_rp - tanh(vn / (this%uinf * this%delta)))
130
131 x(k) = -0.5*(ux*ux+uy*uy+uz*uz)*s0
132 end do
133 end if
134 end subroutine dong_outflow_apply_scalar
135
138 subroutine dong_outflow_apply_vector(this, x, y, z, n, time, strong)
139 class(dong_outflow_t), intent(inout) :: this
140 integer, intent(in) :: n
141 real(kind=rp), intent(inout), dimension(n) :: x
142 real(kind=rp), intent(inout), dimension(n) :: y
143 real(kind=rp), intent(inout), dimension(n) :: z
144 type(time_state_t), intent(in), optional :: time
145 logical, intent(in), optional :: strong
146
147 end subroutine dong_outflow_apply_vector
148
151 subroutine dong_outflow_apply_scalar_dev(this, x_d, time, strong, strm)
152 class(dong_outflow_t), intent(inout), target :: this
153 type(c_ptr), intent(inout) :: x_d
154 type(time_state_t), intent(in), optional :: time
155 logical, intent(in), optional :: strong
156 type(c_ptr), intent(inout) :: strm
157 logical :: strong_
158
159 if (present(strong)) then
160 strong_ = strong
161 else
162 strong_ = .true.
163 end if
164
165 if (strong_ .and. this%msk(0) .gt. 0) then
166 call device_dong_outflow_apply_scalar(this%msk_d, x_d, &
167 this%normal_x_d, this%normal_y_d, this%normal_z_d, &
168 this%u%x_d, this%v%x_d, this%w%x_d, &
169 this%uinf, this%delta, &
170 this%msk(0), strm)
171 end if
172
173 end subroutine dong_outflow_apply_scalar_dev
174
177 subroutine dong_outflow_apply_vector_dev(this, x_d, y_d, z_d, time, &
178 strong, strm)
179 class(dong_outflow_t), intent(inout), target :: this
180 type(c_ptr), intent(inout) :: x_d
181 type(c_ptr), intent(inout) :: y_d
182 type(c_ptr), intent(inout) :: z_d
183 type(time_state_t), intent(in), optional :: time
184 logical, intent(in), optional :: strong
185 type(c_ptr), intent(inout) :: strm
186
187 !call device_dong_outflow_apply_vector(this%msk_d, x_d, y_d, z_d, &
188 ! this%g, size(this%msk))
189
190 end subroutine dong_outflow_apply_vector_dev
191
193 subroutine dong_outflow_free(this)
194 class(dong_outflow_t), target, intent(inout) :: this
195
196 call this%free_base
197 nullify(this%u)
198 nullify(this%v)
199 nullify(this%w)
200
201 if (c_associated(this%normal_x_d)) then
202 call device_free(this%normal_x_d)
203 this%normal_x_d = c_null_ptr
204 end if
205
206 if (c_associated(this%normal_y_d)) then
207 call device_free(this%normal_y_d)
208 this%normal_y_d = c_null_ptr
209 end if
210
211 if (c_associated(this%normal_z_d)) then
212 call device_free(this%normal_z_d)
213 this%normal_z_d = c_null_ptr
214 end if
215
216 end subroutine dong_outflow_free
217
219 subroutine dong_outflow_finalize(this, only_facets)
220 class(dong_outflow_t), target, intent(inout) :: this
221 logical, optional, intent(in) :: only_facets
222 real(kind=rp), allocatable :: temp_x(:)
223 real(kind=rp), allocatable :: temp_y(:)
224 real(kind=rp), allocatable :: temp_z(:)
225 real(c_rp) :: dummy
226 integer :: i, m, k, facet, idx(4)
227 real(kind=rp) :: normal_xyz(3)
228
229 if (present(only_facets)) then
230 if (only_facets .eqv. .false.) then
231 call neko_error("For dong_outflow_t, only_facets has to be true.")
232 end if
233 end if
234
235 call this%finalize_base(.true.)
236
237 this%u => neko_registry%get_field("u")
238 this%v => neko_registry%get_field("v")
239 this%w => neko_registry%get_field("w")
240 if ((neko_bcknd_device .eq. 1) .and. (this%msk(0) .gt. 0)) then
241 call device_alloc(this%normal_x_d, c_sizeof(dummy)*this%msk(0))
242 call device_alloc(this%normal_y_d, c_sizeof(dummy)*this%msk(0))
243 call device_alloc(this%normal_z_d, c_sizeof(dummy)*this%msk(0))
244 m = this%msk(0)
245 allocate(temp_x(m))
246 allocate(temp_y(m))
247 allocate(temp_z(m))
248 do i = 1, m
249 k = this%msk(i)
250 facet = this%facet(i)
251 idx = nonlinear_index(k, this%Xh%lx, this%Xh%lx, this%Xh%lx)
252 normal_xyz = &
253 this%coef%get_normal(idx(1), idx(2), idx(3), idx(4), facet)
254 temp_x(i) = normal_xyz(1)
255 temp_y(i) = normal_xyz(2)
256 temp_z(i) = normal_xyz(3)
257 end do
258 call device_memcpy(temp_x, this%normal_x_d, m, host_to_device, &
259 sync = .false.)
260 call device_memcpy(temp_y, this%normal_y_d, m, host_to_device, &
261 sync = .false.)
262 call device_memcpy(temp_z, this%normal_z_d, m, host_to_device, &
263 sync = .true.)
264 deallocate( temp_x, temp_y, temp_z)
265 end if
266 end subroutine dong_outflow_finalize
267
268end module dong_outflow
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Copy data between host and device (or device and device)
Definition device.F90:71
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
subroutine, public device_dong_outflow_apply_scalar(msk, x, normal_x, normal_y, normal_z, u, v, w, uinf, delta, m, strm)
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:219
subroutine, public device_alloc(x_d, s)
Allocate memory on the device.
Definition device.F90:192
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a dong outflow condition.
subroutine dong_outflow_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Boundary condition apply for a generic Dirichlet condition to vectors x, y and z (device version)
subroutine dong_outflow_apply_vector(this, x, y, z, n, time, strong)
Boundary condition apply for a generic Dirichlet condition to vectors x, y and z.
subroutine dong_outflow_apply_scalar_dev(this, x_d, time, strong, strm)
Boundary condition apply for a generic Dirichlet condition to a vector x (device version)
subroutine dong_outflow_free(this)
Destructor.
subroutine dong_outflow_apply_scalar(this, x, n, time, strong)
Boundary condition apply for a generic Dirichlet condition to a vector x.
subroutine dong_outflow_init(this, coef, json)
Constructor.
subroutine dong_outflow_finalize(this, only_facets)
Finalize.
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public c_rp
Definition num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:128
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
Base type for a boundary condition.
Definition bc.f90:61
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:48
Dong outflow condition Follows "A Convective-like Energy-Stable Open Boundary Condition for Simulati...
A struct that contains all info about the time, expand as needed.