Neko  0.8.99
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  use field_list, only : field_list_t
48  implicit none
49  private
50 
61  type, public, extends(bc_t) :: field_dirichlet_t
64  type(field_t) :: field_bc
69  procedure(field_dirichlet_update), nopass, pointer :: update => null()
70  contains
72  procedure, pass(this) :: init_field => field_dirichlet_init
74  procedure, pass(this) :: apply_scalar => field_dirichlet_apply_scalar
76  procedure, pass(this) :: apply_vector => field_dirichlet_apply_vector
78  procedure, pass(this) :: apply_vector_dev => field_dirichlet_apply_vector_dev
80  procedure, pass(this) :: apply_scalar_dev => field_dirichlet_apply_scalar_dev
82  procedure, pass(this) :: free => field_dirichlet_free
83 
84  end type field_dirichlet_t
85 
95  abstract interface
96  subroutine field_dirichlet_update(dirichlet_field_list, dirichlet_bc_list, coef, t, tstep, which_solver)
97  import rp
98  import field_list_t
99  import bc_list_t
100  import coef_t
101  type(field_list_t), intent(inout) :: dirichlet_field_list
102  type(bc_list_t), intent(inout) :: dirichlet_bc_list
103  type(coef_t), intent(inout) :: coef
104  real(kind=rp), intent(in) :: t
105  integer, intent(in) :: tstep
106  character(len=*), intent(in) :: which_solver
107  end subroutine field_dirichlet_update
108  end interface
109 
110  public :: field_dirichlet_update
111 
112 contains
113 
116  subroutine field_dirichlet_init(this, bc_name)
117  class(field_dirichlet_t), target, intent(inout) :: this
118  character(len=*), intent(in) :: bc_name
119 
120  call this%field_bc%init(this%dof, bc_name)
121  call this%field_list%init(1)
122  call this%field_list%assign_to_field(1, this%field_bc)
123 
124  end subroutine field_dirichlet_init
125 
127  subroutine field_dirichlet_free(this)
128  class(field_dirichlet_t), target, intent(inout) :: this
129 
130  call this%free_base
131  call this%field_bc%free()
132  call this%field_list%free()
133 
134  if (associated(this%update)) then
135  this%update => null()
136  end if
137 
138  end subroutine field_dirichlet_free
139 
145  subroutine field_dirichlet_apply_scalar(this, x, n, t, tstep)
146  class(field_dirichlet_t), intent(inout) :: this
147  integer, intent(in) :: n
148  real(kind=rp), intent(inout), dimension(n) :: x
149  real(kind=rp), intent(in), optional :: t
150  integer, intent(in), optional :: tstep
151 
152  if (this%msk(0) .gt. 0) then
153  call masked_copy(x, this%field_bc%x, this%msk, n, this%msk(0))
154  end if
155 
156  end subroutine field_dirichlet_apply_scalar
157 
162  subroutine field_dirichlet_apply_scalar_dev(this, x_d, t, tstep)
163  class(field_dirichlet_t), intent(inout), target :: this
164  type(c_ptr) :: x_d
165  real(kind=rp), intent(in), optional :: t
166  integer, intent(in), optional :: tstep
167 
168  if (this%msk(0) .gt. 0) then
169  call device_masked_copy(x_d, this%field_bc%x_d, this%msk_d, &
170  this%field_bc%dof%size(), this%msk(0))
171  end if
172 
173  end subroutine field_dirichlet_apply_scalar_dev
174 
182  subroutine field_dirichlet_apply_vector(this, x, y, z, n, t, tstep)
183  class(field_dirichlet_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  real(kind=rp), intent(in), optional :: t
189  integer, intent(in), optional :: tstep
190 
191  call neko_error("field_dirichlet cannot apply vector BCs.&
192 &Use field_dirichlet_vector instead!")
193 
194  end subroutine field_dirichlet_apply_vector
195 
202  subroutine field_dirichlet_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
203  class(field_dirichlet_t), intent(inout), target :: this
204  type(c_ptr) :: x_d
205  type(c_ptr) :: y_d
206  type(c_ptr) :: z_d
207  real(kind=rp), intent(in), optional :: t
208  integer, intent(in), optional :: tstep
209 
210  call neko_error("field_dirichlet cannot apply vector BCs.&
211 &Use field_dirichlet_vector instead!")
212 
213  end subroutine field_dirichlet_apply_vector_dev
214 
215 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)
Copy a masked vector .
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:247
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, public split_string(string, delimiter)
Split a string based on delimiter (tokenizer) OBS: very hacky, this should really be improved,...
Definition: utils.f90:89
A list of boundary conditions.
Definition: bc.f90:104
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
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:13