Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
blasius.f90
Go to the documentation of this file.
1! Copyright (c) 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
40 use flow_profile
41 use, intrinsic :: iso_fortran_env
42 use, intrinsic :: iso_c_binding
43 use bc, only : bc_t
44 use json_module, only : json_file
45 use json_utils, only : json_get
46 implicit none
47 private
48
52 type, public, extends(bc_t) :: blasius_t
53 real(kind=rp), dimension(3) :: uinf = [0d0, 0d0, 0d0]
54 real(kind=rp) :: delta
55 procedure(blasius_profile), nopass, pointer :: bla => null()
56 type(c_ptr), private :: blax_d = c_null_ptr
57 type(c_ptr), private :: blay_d = c_null_ptr
58 type(c_ptr), private :: blaz_d = c_null_ptr
59 contains
60 procedure, pass(this) :: apply_scalar => blasius_apply_scalar
61 procedure, pass(this) :: apply_vector => blasius_apply_vector
62 procedure, pass(this) :: apply_scalar_dev => blasius_apply_scalar_dev
63 procedure, pass(this) :: apply_vector_dev => blasius_apply_vector_dev
64 procedure, pass(this) :: set_params => blasius_set_params
66 procedure, pass(this) :: init => blasius_init
68 procedure, pass(this) :: init_from_components => &
71 procedure, pass(this) :: free => blasius_free
73 procedure, pass(this) :: finalize => blasius_finalize
74 end type blasius_t
75
76contains
77
81 subroutine blasius_init(this, coef, json)
82 class(blasius_t), intent(inout), target :: this
83 type(coef_t), intent(in) :: coef
84 type(json_file), intent(inout) :: json
85 real(kind=rp) :: delta
86 real(kind=rp), allocatable :: uinf(:)
87 character(len=:), allocatable :: approximation
88
89 call this%init_base(coef)
90
91 call json_get(json, 'delta', delta)
92 call json_get(json, 'approximation', approximation)
93 call json_get(json, 'freestream_velocity', uinf)
94
95 if (size(uinf) .ne. 3) then
96 call neko_error("The uinf keyword for the blasius profile should be an &
97& array of 3 reals")
98 end if
99
100 call this%init_from_components(coef, delta, uinf, approximation)
101
102 end subroutine blasius_init
103
109 subroutine blasius_init_from_components(this, coef, delta, uinf, &
110 approximation)
111 class(blasius_t), intent(inout), target :: this
112 type(coef_t), intent(in) :: coef
113 real(kind=rp) :: delta
114 real(kind=rp) :: uinf(3)
115 character(len=*) :: approximation
116
117 call this%init_base(coef)
118
119 this%delta = delta
120 this%uinf = uinf
121
122 select case (trim(approximation))
123 case ('linear')
124 this%bla => blasius_linear
125 case ('quadratic')
126 this%bla => blasius_quadratic
127 case ('cubic')
128 this%bla => blasius_cubic
129 case ('quartic')
130 this%bla => blasius_quartic
131 case ('sin')
132 this%bla => blasius_sin
133 case default
134 call neko_error('Invalid Blasius approximation')
135 end select
136 end subroutine blasius_init_from_components
137
138 subroutine blasius_free(this)
139 class(blasius_t), target, intent(inout) :: this
140
141 call this%free_base()
142 nullify(this%bla)
143
144 if (c_associated(this%blax_d)) then
145 call device_free(this%blax_d)
146 end if
147
148 if (c_associated(this%blay_d)) then
149 call device_free(this%blay_d)
150 end if
151
152 if (c_associated(this%blaz_d)) then
153 call device_free(this%blaz_d)
154 end if
155
156 end subroutine blasius_free
157
159 subroutine blasius_apply_scalar(this, x, n, t, tstep, strong)
160 class(blasius_t), intent(inout) :: this
161 integer, intent(in) :: n
162 real(kind=rp), intent(inout), dimension(n) :: x
163 real(kind=rp), intent(in), optional :: t
164 integer, intent(in), optional :: tstep
165 logical, intent(in), optional :: strong
166 end subroutine blasius_apply_scalar
167
169 subroutine blasius_apply_scalar_dev(this, x_d, t, tstep, strong)
170 class(blasius_t), intent(inout), target :: this
171 type(c_ptr) :: x_d
172 real(kind=rp), intent(in), optional :: t
173 integer, intent(in), optional :: tstep
174 logical, intent(in), optional :: strong
175 end subroutine blasius_apply_scalar_dev
176
178 subroutine blasius_apply_vector(this, x, y, z, n, t, tstep, strong)
179 class(blasius_t), intent(inout) :: this
180 integer, intent(in) :: n
181 real(kind=rp), intent(inout), dimension(n) :: x
182 real(kind=rp), intent(inout), dimension(n) :: y
183 real(kind=rp), intent(inout), dimension(n) :: z
184 real(kind=rp), intent(in), optional :: t
185 integer, intent(in), optional :: tstep
186 logical, intent(in), optional :: strong
187 integer :: i, m, k, idx(4), facet
188 logical :: strong_ = .true.
189
190 if (present(strong)) strong_ = strong
191
192 associate(xc => this%coef%dof%x, yc => this%coef%dof%y, &
193 zc => this%coef%dof%z, nx => this%coef%nx, ny => this%coef%ny, &
194 nz => this%coef%nz, lx => this%coef%Xh%lx)
195 m = this%msk(0)
196 if (strong_) then
197 do i = 1, m
198 k = this%msk(i)
199 facet = this%facet(i)
200 idx = nonlinear_index(k, lx, lx, lx)
201 select case (facet)
202 case (1, 2)
203 x(k) = this%bla(zc(idx(1), idx(2), idx(3), idx(4)), &
204 this%delta, this%uinf(1))
205 y(k) = 0.0_rp
206 z(k) = 0.0_rp
207 case (3, 4)
208 x(k) = 0.0_rp
209 y(k) = this%bla(xc(idx(1), idx(2), idx(3), idx(4)), &
210 this%delta, this%uinf(2))
211 z(k) = 0.0_rp
212 case (5, 6)
213 x(k) = 0.0_rp
214 y(k) = 0.0_rp
215 z(k) = this%bla(yc(idx(1), idx(2), idx(3), idx(4)), &
216 this%delta, this%uinf(3))
217 end select
218 end do
219 end if
220 end associate
221 end subroutine blasius_apply_vector
222
224 subroutine blasius_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
225 class(blasius_t), intent(inout), target :: this
226 type(c_ptr) :: x_d
227 type(c_ptr) :: y_d
228 type(c_ptr) :: z_d
229 real(kind=rp), intent(in), optional :: t
230 integer, intent(in), optional :: tstep
231 logical, intent(in), optional :: strong
232 integer :: i, m, k, idx(4), facet
233 integer(c_size_t) :: s
234 real(kind=rp), allocatable :: bla_x(:), bla_y(:), bla_z(:)
235 logical :: strong_ = .true.
236
237 if (present(strong)) strong_ = strong
238
239 associate(xc => this%coef%dof%x, yc => this%coef%dof%y, &
240 zc => this%coef%dof%z, nx => this%coef%nx, ny => this%coef%ny, &
241 nz => this%coef%nz, lx => this%coef%Xh%lx , &
242 blax_d => this%blax_d, blay_d => this%blay_d, &
243 blaz_d => this%blaz_d)
244
245 m = this%msk(0)
246
247
248 ! Pretabulate values during first call to apply
249 if (.not. c_associated(blax_d) .and. strong_ .and. this%msk(0) .gt. 0) then
250 allocate(bla_x(m), bla_y(m), bla_z(m)) ! Temp arrays
251
252 if (rp .eq. real32) then
253 s = m * 4
254 else if (rp .eq. real64) then
255 s = m * 8
256 end if
257
258 call device_alloc(blax_d, s)
259 call device_alloc(blay_d, s)
260 call device_alloc(blaz_d, s)
261
262 do i = 1, m
263 k = this%msk(i)
264 facet = this%facet(i)
265 idx = nonlinear_index(k, lx, lx, lx)
266 select case (facet)
267 case (1,2)
268 bla_x(i) = this%bla(zc(idx(1), idx(2), idx(3), idx(4)), &
269 this%delta, this%uinf(1))
270 bla_y(i) = 0.0_rp
271 bla_z(i) = 0.0_rp
272 case (3,4)
273 bla_x(i) = 0.0_rp
274 bla_y(i) = this%bla(xc(idx(1), idx(2), idx(3), idx(4)), &
275 this%delta, this%uinf(2))
276 bla_z(i) = 0.0_rp
277 case (5,6)
278 bla_x(i) = 0.0_rp
279 bla_y(i) = 0.0_rp
280 bla_z(i) = this%bla(yc(idx(1), idx(2), idx(3), idx(4)), &
281 this%delta, this%uinf(3))
282 end select
283 end do
284
285 call device_memcpy(bla_x, blax_d, m, host_to_device, sync = .false.)
286 call device_memcpy(bla_y, blay_d, m, host_to_device, sync = .false.)
287 call device_memcpy(bla_z, blaz_d, m, host_to_device, sync = .true.)
288
289 deallocate(bla_x, bla_y, bla_z)
290 end if
291
292 if (strong_ .and. this%msk(0) .gt. 0) then
293 call device_inhom_dirichlet_apply_vector(this%msk_d, x_d, y_d, z_d, &
294 blax_d, blay_d, blaz_d, m)
295 end if
296
297 end associate
298
299 end subroutine blasius_apply_vector_dev
300
302 subroutine blasius_set_params(this, uinf, delta, type)
303 class(blasius_t), intent(inout) :: this
304 real(kind=rp), intent(in) :: uinf(3)
305 real(kind=rp), intent(in) :: delta
306 character(len=*) :: type
307 this%delta = delta
308 this%uinf = uinf
309
310 select case (trim(type))
311 case ('linear')
312 this%bla => blasius_linear
313 case ('quadratic')
314 this%bla => blasius_quadratic
315 case ('cubic')
316 this%bla => blasius_cubic
317 case ('quartic')
318 this%bla => blasius_quartic
319 case ('sin')
320 this%bla => blasius_sin
321 case default
322 call neko_error('Invalid Blasius approximation')
323 end select
324 end subroutine blasius_set_params
325
327 subroutine blasius_finalize(this)
328 class(blasius_t), target, intent(inout) :: this
329
330 call this%finalize_base()
331 end subroutine blasius_finalize
332end module blasius
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Copy data between host and device (or device and device)
Definition device.F90:65
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_apply_scalar(this, x, n, t, tstep, strong)
No-op scalar apply.
Definition blasius.f90:160
subroutine blasius_init(this, coef, json)
Constructor.
Definition blasius.f90:82
subroutine blasius_apply_scalar_dev(this, x_d, t, tstep, strong)
No-op scalar apply (device version)
Definition blasius.f90:170
subroutine blasius_init_from_components(this, coef, delta, uinf, approximation)
Constructor from components.
Definition blasius.f90:111
subroutine blasius_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
Apply blasius conditions (vector valued) (device version)
Definition blasius.f90:225
subroutine blasius_free(this)
Definition blasius.f90:139
subroutine blasius_apply_vector(this, x, y, z, n, t, tstep, strong)
Apply blasius conditions (vector valued)
Definition blasius.f90:179
subroutine blasius_finalize(this)
Finalize.
Definition blasius.f90:328
subroutine blasius_set_params(this, uinf, delta, type)
Set Blasius parameters.
Definition blasius.f90:303
Coefficients.
Definition coef.f90:34
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:46
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:200
subroutine, public device_alloc(x_d, s)
Allocate memory on the device.
Definition device.F90:179
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_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
Utilities.
Definition utils.f90:35
Base type for a boundary condition.
Definition bc.f90:57
Blasius profile for inlet (vector valued).
Definition blasius.f90:52
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55