Neko  0.8.1
A portable framework for high-order spectral element flow simulations
field_dirichlet.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2024, 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 num_types, only: rp
36  use coefs, only: coef_t
37  use dirichlet, only: dirichlet_t
38  use bc, only: bc_list_t, bc_t
39  use device, only: c_ptr, c_size_t
40  use utils, only: split_string
41  use field, only : field_t
42  use field_list, only : field_list_t
43  use math, only: masked_copy
45  use dofmap, only : dofmap_t
46  use utils, only: neko_error
47  implicit none
48  private
49 
55  type, public, extends(bc_t) :: field_dirichlet_t
56  type(field_t) :: field_bc
57  contains
59  procedure, pass(this) :: init_field => field_dirichlet_init
61  procedure, pass(this) :: apply_scalar => field_dirichlet_apply_scalar
63  procedure, pass(this) :: apply_vector => field_dirichlet_apply_vector
65  procedure, pass(this) :: apply_vector_dev => field_dirichlet_apply_vector_dev
67  procedure, pass(this) :: apply_scalar_dev => field_dirichlet_apply_scalar_dev
68  end type field_dirichlet_t
69 
77  abstract interface
78  subroutine field_dirichlet_update(dirichlet_field_list, dirichlet_bc_list, coef, t, tstep, which_solver)
79  import rp
80  import field_list_t
81  import bc_list_t
82  import coef_t
83  type(field_list_t), intent(inout) :: dirichlet_field_list
84  type(bc_list_t), intent(inout) :: dirichlet_bc_list
85  type(coef_t), intent(inout) :: coef
86  real(kind=rp), intent(in) :: t
87  integer, intent(in) :: tstep
88  character(len=*), intent(in) :: which_solver
89  end subroutine field_dirichlet_update
90  end interface
91 
92  public :: field_dirichlet_update
93 
94 contains
95 
98  subroutine field_dirichlet_init(this, bc_name)
99  class(field_dirichlet_t), intent(inout) :: this
100  character(len=*), intent(in) :: bc_name
101 
102  call this%field_bc%init(this%dof, bc_name)
103 
104  end subroutine field_dirichlet_init
105 
107  subroutine field_dirichlet_free(this)
108  type(field_dirichlet_t), intent(inout) :: this
109 
110  call this%field_bc%free()
111 
112  end subroutine field_dirichlet_free
113 
119  subroutine field_dirichlet_apply_scalar(this, x, n, t, tstep)
120  class(field_dirichlet_t), intent(inout) :: this
121  integer, intent(in) :: n
122  real(kind=rp), intent(inout), dimension(n) :: x
123  real(kind=rp), intent(in), optional :: t
124  integer, intent(in), optional :: tstep
125 
126  if (this%msk(0) .gt. 0) then
127  call masked_copy(x, this%field_bc%x, this%msk, n, this%msk(0))
128  end if
129 
130  end subroutine field_dirichlet_apply_scalar
131 
136  subroutine field_dirichlet_apply_scalar_dev(this, x_d, t, tstep)
137  class(field_dirichlet_t), intent(inout), target :: this
138  type(c_ptr) :: x_d
139  real(kind=rp), intent(in), optional :: t
140  integer, intent(in), optional :: tstep
141 
142  if (this%msk(0) .gt. 0) then
143  call device_masked_copy(x_d, this%field_bc%x_d, this%msk_d, &
144  this%field_bc%dof%size(), this%msk(0))
145  end if
146 
147  end subroutine field_dirichlet_apply_scalar_dev
148 
156  subroutine field_dirichlet_apply_vector(this, x, y, z, n, t, tstep)
157  class(field_dirichlet_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  call neko_error("field_dirichlet cannot apply vector BCs.&
166 &Use field_dirichlet_vector instead!")
167 
168  end subroutine field_dirichlet_apply_vector
169 
176  subroutine field_dirichlet_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
177  class(field_dirichlet_t), intent(inout), target :: this
178  type(c_ptr) :: x_d
179  type(c_ptr) :: y_d
180  type(c_ptr) :: z_d
181  real(kind=rp), intent(in), optional :: t
182  integer, intent(in), optional :: tstep
183 
184  call neko_error("field_dirichlet cannot apply vector BCs.&
185 &Use field_dirichlet_vector instead!")
186 
187  end subroutine field_dirichlet_apply_vector_dev
188 
189 end module field_dirichlet
Abstract interface defining a dirichlet condition on a list of fields.
Defines a boundary condition.
Definition: bc.f90:34
Coefficients.
Definition: coef.f90:34
subroutine, public device_masked_copy(a_d, b_d, mask_d, n, m)
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Defines inflow dirichlet conditions.
subroutine field_dirichlet_apply_vector(this, x, y, z, n, t, tstep)
(No-op) Apply vector.
subroutine field_dirichlet_apply_scalar_dev(this, x_d, t, tstep)
Apply scalar (device).
subroutine field_dirichlet_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
(No-op) Apply vector (device).
subroutine field_dirichlet_apply_scalar(this, x, n, t, tstep)
Apply scalar by performing a masked copy.
subroutine field_dirichlet_init(this, bc_name)
Initializes thisfield_bc.
subroutine field_dirichlet_free(this)
Destructor. Currently thisfield_bc is being freed in fluid_scheme::free
Defines a field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public masked_copy(a, b, mask, n, m)
Copy a masked vector .
Definition: math.f90:230
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35
character(len=100) function, dimension(:), allocatable split_string(string, delimiter)
Split a string based on delimiter (tokenizer) OBS: very hacky, this should really be improved,...
Definition: utils.f90:82
A list of boundary conditions.
Definition: bc.f90:102
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:54
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44
User defined dirichlet condition, for which the user can work with an entire field....
field_list_t, To be able to group fields together
Definition: field_list.f90:7