Neko  0.8.99
A portable framework for high-order spectral element flow simulations
dong_outflow.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022, 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
37  use device
38  use num_types
39  use bc
40  use field
41  use dofmap
42  use coefs
43  use utils
46  use, intrinsic :: iso_c_binding, only : c_ptr, c_sizeof
47  use json_module, only : json_file
49  implicit none
50  private
51 
57  type, public, extends(bc_t) :: dong_outflow_t
58  type(field_t), pointer :: u
59  type(field_t), pointer :: v
60  type(field_t), pointer :: w
61  real(kind=rp) :: delta
62  real(kind=rp) :: uinf
63  type(c_ptr) :: normal_x_d
64  type(c_ptr) :: normal_y_d
65  type(c_ptr) :: normal_z_d
66  contains
67  procedure, pass(this) :: apply_scalar => dong_outflow_apply_scalar
68  procedure, pass(this) :: apply_vector => dong_outflow_apply_vector
69  procedure, pass(this) :: apply_scalar_dev => dong_outflow_apply_scalar_dev
70  procedure, pass(this) :: apply_vector_dev => dong_outflow_apply_vector_dev
71  procedure, pass(this) :: init => dong_outflow_init
73  procedure, pass(this) :: free => dong_outflow_free
74  end type dong_outflow_t
75 
76 contains
77  subroutine dong_outflow_init(this, coef, json)
78  class(dong_outflow_t), intent(inout) :: this
79  type(coef_t), intent(in) :: coef
80  type(json_file), intent(inout) :: json
81  real(kind=rp), allocatable :: temp_x(:)
82  real(kind=rp), allocatable :: temp_y(:)
83  real(kind=rp), allocatable :: temp_z(:)
84  real(c_rp) :: dummy
85  integer :: i, m, k, facet, idx(4)
86  real(kind=rp) :: normal_xyz(3)
87 
88 ! call this%dirichlet_t%init
89 
90  call json_get_or_default(json, 'case.fluid.outflow_condition.delta', &
91  this%delta, 0.01_rp)
92  call json_get_or_default(json, 'case.fluid.outflow_condition.velocity_scale', &
93  this%uinf, 1.0_rp)
94 
95  this%u => neko_field_registry%get_field("u")
96  this%v => neko_field_registry%get_field("v")
97  this%w => neko_field_registry%get_field("w")
98 
99  if ((neko_bcknd_device .eq. 1) .and. (this%msk(0) .gt. 0)) then
100  call device_alloc(this%normal_x_d, c_sizeof(dummy)*this%msk(0))
101  call device_alloc(this%normal_y_d, c_sizeof(dummy)*this%msk(0))
102  call device_alloc(this%normal_z_d, c_sizeof(dummy)*this%msk(0))
103  m = this%msk(0)
104  allocate(temp_x(m))
105  allocate(temp_y(m))
106  allocate(temp_z(m))
107  do i = 1, m
108  k = this%msk(i)
109  facet = this%facet(i)
110  idx = nonlinear_index(k,this%Xh%lx, this%Xh%lx,this%Xh%lx)
111  normal_xyz = &
112  this%coef%get_normal(idx(1), idx(2), idx(3), idx(4),facet)
113  temp_x(i) = normal_xyz(1)
114  temp_y(i) = normal_xyz(2)
115  temp_z(i) = normal_xyz(3)
116  end do
117  call device_memcpy(temp_x, this%normal_x_d, m, &
118  host_to_device, sync=.false.)
119  call device_memcpy(temp_y, this%normal_y_d, m, &
120  host_to_device, sync=.false.)
121  call device_memcpy(temp_z, this%normal_z_d, m, &
122  host_to_device, sync=.true.)
123  deallocate( temp_x, temp_y, temp_z)
124  end if
125  end subroutine dong_outflow_init
126 
129  subroutine dong_outflow_apply_scalar(this, x, n, t, tstep)
130  class(dong_outflow_t), intent(inout) :: this
131  integer, intent(in) :: n
132  real(kind=rp), intent(inout), dimension(n) :: x
133  real(kind=rp), intent(in), optional :: t
134  integer, intent(in), optional :: tstep
135  integer :: i, m, k, facet, idx(4)
136  real(kind=rp) :: vn, s0, ux, uy, uz, normal_xyz(3)
137 
138  m = this%msk(0)
139  do i = 1, m
140  k = this%msk(i)
141  facet = this%facet(i)
142  ux = this%u%x(k,1,1,1)
143  uy = this%v%x(k,1,1,1)
144  uz = this%w%x(k,1,1,1)
145  idx = nonlinear_index(k,this%Xh%lx, this%Xh%lx,this%Xh%lx)
146  normal_xyz = this%coef%get_normal(idx(1), idx(2), idx(3), idx(4),facet)
147  vn = ux*normal_xyz(1) + uy*normal_xyz(2) + uz*normal_xyz(3)
148  s0 = 0.5_rp*(1.0_rp - tanh(vn / (this%uinf * this%delta)))
149 
150  x(k)=-0.5*(ux*ux+uy*uy+uz*uz)*s0
151  end do
152  end subroutine dong_outflow_apply_scalar
153 
156  subroutine dong_outflow_apply_vector(this, x, y, z, n, t, tstep)
157  class(dong_outflow_t), intent(inout) :: this
158  integer, intent(in) :: n
159  real(kind=rp), intent(inout), dimension(n) :: x
160  real(kind=rp), intent(inout), dimension(n) :: y
161  real(kind=rp), intent(inout), dimension(n) :: z
162  real(kind=rp), intent(in), optional :: t
163  integer, intent(in), optional :: tstep
164 
165  end subroutine dong_outflow_apply_vector
166 
169  subroutine dong_outflow_apply_scalar_dev(this, x_d, t, tstep)
170  class(dong_outflow_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 
175  call device_dong_outflow_apply_scalar(this%msk_d,x_d, this%normal_x_d, &
176  this%normal_y_d, this%normal_z_d,&
177  this%u%x_d, this%v%x_d, this%w%x_d,&
178  this%uinf, this%delta,&
179  this%msk(0))
180 
181  end subroutine dong_outflow_apply_scalar_dev
182 
185  subroutine dong_outflow_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
186  class(dong_outflow_t), intent(inout), target :: this
187  type(c_ptr) :: x_d
188  type(c_ptr) :: y_d
189  type(c_ptr) :: z_d
190  real(kind=rp), intent(in), optional :: t
191  integer, intent(in), optional :: tstep
192 
193  !call device_dong_outflow_apply_vector(this%msk_d, x_d, y_d, z_d, &
194  ! this%g, size(this%msk))
195 
196  end subroutine dong_outflow_apply_vector_dev
197 
199  subroutine dong_outflow_free(this)
200  class(dong_outflow_t), target, intent(inout) :: this
201 
202  call this%free_base
203 
204  end subroutine dong_outflow_free
205 
206 end module dong_outflow
__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:51
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:53
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:44
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)
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public host_to_device
Definition: device.F90:47
subroutine, public device_alloc(x_d, s)
Allocate memory on the device.
Definition: device.F90:164
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(this, x, y, z, n, t, tstep)
Boundary condition apply for a generic Dirichlet condition to vectors x, y and z.
subroutine dong_outflow_apply_scalar(this, x, n, t, tstep)
Boundary condition apply for a generic Dirichlet condition to a vector x.
subroutine dong_outflow_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
Boundary condition apply for a generic Dirichlet condition to vectors x, y and z (device version)
subroutine dong_outflow_apply_scalar_dev(this, x_d, t, tstep)
Boundary condition apply for a generic Dirichlet condition to a vector x (device version)
subroutine dong_outflow_free(this)
Destructor.
subroutine dong_outflow_init(this, coef, json)
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines a field.
Definition: field.f90:34
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
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:51
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
Dong outflow condition Follows "A Convective-like Energy-Stable Open Boundary Condition for Simulati...