Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
blasius.f90
Go to the documentation of this file.
1! Copyright (c) 2021-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!
34module blasius
35 use num_types, only : rp
36 use coefs, only : coef_t
37 use utils, only : nonlinear_index
42 use utils, only : neko_error
43 use, intrinsic :: iso_fortran_env
44 use, intrinsic :: iso_c_binding
45 use bc, only : bc_t
46 use json_module, only : json_file
48 use time_state, only : time_state_t
49 implicit none
50 private
51
55 type, public, extends(bc_t) :: blasius_t
56 real(kind=rp), dimension(3) :: uinf = [0d0, 0d0, 0d0]
57 real(kind=rp) :: delta
58 procedure(blasius_profile), nopass, pointer :: bla => null()
59 type(c_ptr), private :: blax_d = c_null_ptr
60 type(c_ptr), private :: blay_d = c_null_ptr
61 type(c_ptr), private :: blaz_d = c_null_ptr
62 contains
63 procedure, pass(this) :: apply_scalar => blasius_apply_scalar
64 procedure, pass(this) :: apply_vector => blasius_apply_vector
65 procedure, pass(this) :: apply_scalar_dev => blasius_apply_scalar_dev
66 procedure, pass(this) :: apply_vector_dev => blasius_apply_vector_dev
67 procedure, pass(this) :: set_params => blasius_set_params
69 procedure, pass(this) :: init => blasius_init
71 procedure, pass(this) :: init_from_components => &
74 procedure, pass(this) :: free => blasius_free
76 procedure, pass(this) :: finalize => blasius_finalize
77 end type blasius_t
78
79contains
80
84 subroutine blasius_init(this, coef, json)
85 class(blasius_t), intent(inout), target :: this
86 type(coef_t), target, intent(in) :: coef
87 type(json_file), intent(inout) :: json
88 real(kind=rp) :: delta
89 real(kind=rp), allocatable :: uinf(:)
90 character(len=:), allocatable :: approximation
91
92 call this%init_base(coef)
93
94 call json_get_or_lookup(json, 'delta', delta)
95 call json_get(json, 'approximation', approximation)
96 call json_get_or_lookup(json, 'freestream_velocity', uinf)
97
98 if (size(uinf) .ne. 3) then
99 call neko_error("The uinf keyword for the blasius profile should be an &
100 & array of 3 reals")
101 end if
102
103 call this%init_from_components(coef, delta, uinf, approximation)
104
105 end subroutine blasius_init
106
112 subroutine blasius_init_from_components(this, coef, delta, uinf, &
113 approximation)
114 class(blasius_t), intent(inout), target :: this
115 type(coef_t), target, intent(in) :: coef
116 real(kind=rp) :: delta
117 real(kind=rp) :: uinf(3)
118 character(len=*) :: approximation
119
120 call this%init_base(coef)
121
122 this%delta = delta
123 this%uinf = uinf
124
125 select case (trim(approximation))
126 case ('linear')
127 this%bla => blasius_linear
128 case ('quadratic')
129 this%bla => blasius_quadratic
130 case ('cubic')
131 this%bla => blasius_cubic
132 case ('quartic')
133 this%bla => blasius_quartic
134 case ('sin')
135 this%bla => blasius_sin
136 case ('tanh')
137 this%bla => blasius_tanh
138 case default
139 call neko_error('Invalid Blasius approximation')
140 end select
141 end subroutine blasius_init_from_components
142
143 subroutine blasius_free(this)
144 class(blasius_t), target, intent(inout) :: this
145
146 call this%free_base()
147 nullify(this%bla)
148
149 if (c_associated(this%blax_d)) then
150 call device_free(this%blax_d)
151 end if
152
153 if (c_associated(this%blay_d)) then
154 call device_free(this%blay_d)
155 end if
156
157 if (c_associated(this%blaz_d)) then
158 call device_free(this%blaz_d)
159 end if
160
161 end subroutine blasius_free
162
164 subroutine blasius_apply_scalar(this, x, n, time, strong)
165 class(blasius_t), intent(inout) :: this
166 integer, intent(in) :: n
167 real(kind=rp), intent(inout), dimension(n) :: x
168 type(time_state_t), intent(in), optional :: time
169 logical, intent(in), optional :: strong
170 end subroutine blasius_apply_scalar
171
173 subroutine blasius_apply_scalar_dev(this, x_d, time, strong, strm)
174 class(blasius_t), intent(inout), target :: this
175 type(c_ptr), intent(inout) :: x_d
176 type(time_state_t), intent(in), optional :: time
177 logical, intent(in), optional :: strong
178 type(c_ptr), intent(inout) :: strm
179 end subroutine blasius_apply_scalar_dev
180
182 subroutine blasius_apply_vector(this, x, y, z, n, time, strong)
183 class(blasius_t), intent(inout) :: this
184 integer, intent(in) :: n
185 real(kind=rp), intent(inout), dimension(n) :: x
186 real(kind=rp), intent(inout), dimension(n) :: y
187 real(kind=rp), intent(inout), dimension(n) :: z
188 type(time_state_t), intent(in), optional :: time
189 logical, intent(in), optional :: strong
190 integer :: i, m, k, idx(4), facet
191 logical :: strong_
192
193 if (present(strong)) then
194 strong_ = strong
195 else
196 strong_ = .true.
197 end if
198
199 associate(xc => this%coef%dof%x, yc => this%coef%dof%y, &
200 zc => this%coef%dof%z, nx => this%coef%nx, ny => this%coef%ny, &
201 nz => this%coef%nz, lx => this%coef%Xh%lx)
202 m = this%msk(0)
203 if (strong_) then
204 do i = 1, m
205 k = this%msk(i)
206 facet = this%facet(i)
207 idx = nonlinear_index(k, lx, lx, lx)
208 select case (facet)
209 case (1, 2)
210 x(k) = this%bla(zc(idx(1), idx(2), idx(3), idx(4)), &
211 this%delta, this%uinf(1))
212 y(k) = 0.0_rp
213 z(k) = 0.0_rp
214 case (3, 4)
215 x(k) = 0.0_rp
216 y(k) = this%bla(xc(idx(1), idx(2), idx(3), idx(4)), &
217 this%delta, this%uinf(2))
218 z(k) = 0.0_rp
219 case (5, 6)
220 x(k) = 0.0_rp
221 y(k) = 0.0_rp
222 z(k) = this%bla(yc(idx(1), idx(2), idx(3), idx(4)), &
223 this%delta, this%uinf(3))
224 end select
225 end do
226 end if
227 end associate
228 end subroutine blasius_apply_vector
229
231 subroutine blasius_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
232 class(blasius_t), intent(inout), target :: this
233 type(c_ptr), intent(inout) :: x_d
234 type(c_ptr), intent(inout) :: y_d
235 type(c_ptr), intent(inout) :: z_d
236 type(time_state_t), intent(in), optional :: time
237 logical, intent(in), optional :: strong
238 integer :: i, m, k, idx(4), facet
239 integer(c_size_t) :: s
240 real(kind=rp), allocatable :: bla_x(:), bla_y(:), bla_z(:)
241 logical :: strong_
242 type(c_ptr), intent(inout) :: strm
243
244 if (present(strong)) then
245 strong_ = strong
246 else
247 strong_ = .true.
248 end if
249
250 associate(xc => this%coef%dof%x, yc => this%coef%dof%y, &
251 zc => this%coef%dof%z, nx => this%coef%nx, ny => this%coef%ny, &
252 nz => this%coef%nz, lx => this%coef%Xh%lx , &
253 blax_d => this%blax_d, blay_d => this%blay_d, &
254 blaz_d => this%blaz_d)
255
256 m = this%msk(0)
257
258
259 ! Pretabulate values during first call to apply
260 if (.not. c_associated(blax_d) .and. strong_ .and. m .gt. 0) then
261 allocate(bla_x(m), bla_y(m), bla_z(m)) ! Temp arrays
262
263 if (rp .eq. real32) then
264 s = m * 4
265 else if (rp .eq. real64) then
266 s = m * 8
267 end if
268
269 call device_alloc(blax_d, s)
270 call device_alloc(blay_d, s)
271 call device_alloc(blaz_d, s)
272
273 do i = 1, m
274 k = this%msk(i)
275 facet = this%facet(i)
276 idx = nonlinear_index(k, lx, lx, lx)
277 select case (facet)
278 case (1,2)
279 bla_x(i) = this%bla(zc(idx(1), idx(2), idx(3), idx(4)), &
280 this%delta, this%uinf(1))
281 bla_y(i) = 0.0_rp
282 bla_z(i) = 0.0_rp
283 case (3,4)
284 bla_x(i) = 0.0_rp
285 bla_y(i) = this%bla(xc(idx(1), idx(2), idx(3), idx(4)), &
286 this%delta, this%uinf(2))
287 bla_z(i) = 0.0_rp
288 case (5,6)
289 bla_x(i) = 0.0_rp
290 bla_y(i) = 0.0_rp
291 bla_z(i) = this%bla(yc(idx(1), idx(2), idx(3), idx(4)), &
292 this%delta, this%uinf(3))
293 end select
294 end do
295
296 call device_memcpy(bla_x, blax_d, m, host_to_device, sync = .false.)
297 call device_memcpy(bla_y, blay_d, m, host_to_device, sync = .false.)
298 call device_memcpy(bla_z, blaz_d, m, host_to_device, sync = .true.)
299
300 deallocate(bla_x, bla_y, bla_z)
301 end if
302
303 if (strong_ .and. this%msk(0) .gt. 0) then
304 call device_inhom_dirichlet_apply_vector(this%msk_d, x_d, y_d, z_d, &
305 blax_d, blay_d, blaz_d, m, strm)
306 end if
307
308 end associate
309
310 end subroutine blasius_apply_vector_dev
311
313 subroutine blasius_set_params(this, uinf, delta, type)
314 class(blasius_t), intent(inout) :: this
315 real(kind=rp), intent(in) :: uinf(3)
316 real(kind=rp), intent(in) :: delta
317 character(len=*) :: type
318 this%delta = delta
319 this%uinf = uinf
320
321 select case (trim(type))
322 case ('linear')
323 this%bla => blasius_linear
324 case ('quadratic')
325 this%bla => blasius_quadratic
326 case ('cubic')
327 this%bla => blasius_cubic
328 case ('quartic')
329 this%bla => blasius_quartic
330 case ('sin')
331 this%bla => blasius_sin
332 case ('tanh')
333 this%bla => blasius_tanh
334 case default
335 call neko_error('Invalid Blasius approximation')
336 end select
337 end subroutine blasius_set_params
338
340 subroutine blasius_finalize(this, only_facets)
341 class(blasius_t), target, intent(inout) :: this
342 logical, optional, intent(in) :: only_facets
343 logical :: only_facets_
344
345 if (present(only_facets)) then
346 only_facets_ = only_facets
347 else
348 only_facets_ = .false.
349 end if
350
351 call this%finalize_base(only_facets_)
352 end subroutine blasius_finalize
353end module blasius
__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
Abstract interface for computing a Blasius flow profile.
Retrieves a parameter by name or throws an error.
Defines a boundary condition.
Definition bc.f90:34
Defines a Blasius profile dirichlet condition.
Definition blasius.f90:34
subroutine blasius_finalize(this, only_facets)
Finalize.
Definition blasius.f90:341
subroutine blasius_apply_scalar_dev(this, x_d, time, strong, strm)
No-op scalar apply (device version)
Definition blasius.f90:174
subroutine blasius_init(this, coef, json)
Constructor.
Definition blasius.f90:85
subroutine blasius_init_from_components(this, coef, delta, uinf, approximation)
Constructor from components.
Definition blasius.f90:114
subroutine blasius_apply_vector(this, x, y, z, n, time, strong)
Apply blasius conditions (vector valued)
Definition blasius.f90:183
subroutine blasius_apply_scalar(this, x, n, time, strong)
No-op scalar apply.
Definition blasius.f90:165
subroutine blasius_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Apply blasius conditions (vector valued) (device version)
Definition blasius.f90:232
subroutine blasius_free(this)
Definition blasius.f90:144
subroutine blasius_set_params(this, uinf, delta, type)
Set Blasius parameters.
Definition blasius.f90:314
Coefficients.
Definition coef.f90:34
subroutine device_inhom_dirichlet_apply_vector(msk, x, y, z, bla_x, bla_y, bla_z, 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 flow profile.
real(kind=rp) function, public blasius_quadratic(y, delta, u)
Quadratic approximate Blasius Profile .
real(kind=rp) function, public blasius_quartic(y, delta, u)
Quartic approximate Blasius Profile .
real(kind=rp) function, public blasius_sin(y, delta, u)
Sinusoidal approximate Blasius Profile .
real(kind=rp) function, public blasius_cubic(y, delta, u)
Cubic approximate Blasius Profile .
real(kind=rp) function, public blasius_tanh(y, delta, u)
Hyperbolic tangent approximate Blasius Profile from O. Savas (2012) where is the 99 percent thickne...
real(kind=rp) function, public blasius_linear(y, delta, u)
Linear approximate Blasius profile .
Utilities for retrieving parameters from the case files.
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
Base type for a boundary condition.
Definition bc.f90:62
Blasius profile for inlet (vector valued).
Definition blasius.f90:55
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
A struct that contains all info about the time, expand as needed.