43 use,
intrinsic :: iso_fortran_env
44 use,
intrinsic :: iso_c_binding
46 use json_module,
only : json_file
56 real(kind=
rp),
dimension(3) :: uinf = [0d0, 0d0, 0d0]
57 real(kind=
rp) :: delta
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
71 procedure, pass(this) :: init_from_components => &
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
92 call this%init_base(coef)
95 call json_get(json,
'approximation', approximation)
98 if (
size(uinf) .ne. 3)
then
99 call neko_error(
"The uinf keyword for the blasius profile should be an &
103 call this%init_from_components(coef, delta, uinf, 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
120 call this%init_base(coef)
125 select case (trim(approximation))
139 call neko_error(
'Invalid Blasius approximation')
144 class(
blasius_t),
target,
intent(inout) :: this
146 call this%free_base()
149 if (c_associated(this%blax_d))
then
153 if (c_associated(this%blay_d))
then
157 if (c_associated(this%blaz_d))
then
166 integer,
intent(in) :: n
167 real(kind=
rp),
intent(inout),
dimension(n) :: x
169 logical,
intent(in),
optional :: strong
174 class(
blasius_t),
intent(inout),
target :: this
175 type(c_ptr),
intent(inout) :: x_d
177 logical,
intent(in),
optional :: strong
178 type(c_ptr),
intent(inout) :: strm
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
189 logical,
intent(in),
optional :: strong
190 integer :: i, m, k, idx(4), facet
193 if (
present(strong))
then
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)
206 facet = this%facet(i)
210 x(k) = this%bla(zc(idx(1), idx(2), idx(3), idx(4)), &
211 this%delta, this%uinf(1))
216 y(k) = this%bla(xc(idx(1), idx(2), idx(3), idx(4)), &
217 this%delta, this%uinf(2))
222 z(k) = this%bla(yc(idx(1), idx(2), idx(3), idx(4)), &
223 this%delta, this%uinf(3))
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(:)
242 type(c_ptr),
intent(inout) :: strm
244 if (
present(strong))
then
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)
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))
263 if (rp .eq. real32)
then
265 else if (rp .eq. real64)
then
269 call device_alloc(blax_d, s)
270 call device_alloc(blay_d, s)
271 call device_alloc(blaz_d, s)
275 facet = this%facet(i)
279 bla_x(i) = this%bla(zc(idx(1), idx(2), idx(3), idx(4)), &
280 this%delta, this%uinf(1))
285 bla_y(i) = this%bla(xc(idx(1), idx(2), idx(3), idx(4)), &
286 this%delta, this%uinf(2))
291 bla_z(i) = this%bla(yc(idx(1), idx(2), idx(3), idx(4)), &
292 this%delta, this%uinf(3))
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.)
300 deallocate(bla_x, bla_y, bla_z)
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)
315 real(kind=rp),
intent(in) :: uinf(3)
316 real(kind=rp),
intent(in) :: delta
317 character(len=*) :: type
321 select case (trim(type))
323 this%bla => blasius_linear
325 this%bla => blasius_quadratic
327 this%bla => blasius_cubic
329 this%bla => blasius_quartic
331 this%bla => blasius_sin
333 this%bla => blasius_tanh
335 call neko_error(
'Invalid Blasius approximation')
341 class(
blasius_t),
target,
intent(inout) :: this
342 logical,
optional,
intent(in) :: only_facets
343 logical :: only_facets_
345 if (
present(only_facets))
then
346 only_facets_ = only_facets
348 only_facets_ = .false.
351 call this%finalize_base(only_facets_)
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Copy data between host and device (or device and device)
Abstract interface for computing a Blasius flow profile.
Retrieves a parameter by name or throws an error.
Defines a boundary condition.
Defines a Blasius profile dirichlet condition.
subroutine blasius_finalize(this, only_facets)
Finalize.
subroutine blasius_apply_scalar_dev(this, x_d, time, strong, strm)
No-op scalar apply (device version)
subroutine blasius_init(this, coef, json)
Constructor.
subroutine blasius_init_from_components(this, coef, delta, uinf, approximation)
Constructor from components.
subroutine blasius_apply_vector(this, x, y, z, n, time, strong)
Apply blasius conditions (vector valued)
subroutine blasius_apply_scalar(this, x, n, time, strong)
No-op scalar apply.
subroutine blasius_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Apply blasius conditions (vector valued) (device version)
subroutine blasius_free(this)
subroutine blasius_set_params(this, uinf, delta, type)
Set Blasius parameters.
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.
integer, parameter, public host_to_device
subroutine, public device_free(x_d)
Deallocate memory on the device.
subroutine, public device_alloc(x_d, s)
Allocate memory on the device.
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.
Module with things related to the simulation time.
Base type for a boundary condition.
Blasius profile for inlet (vector valued).
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
A struct that contains all info about the time, expand as needed.