Neko  0.9.99
A portable framework for high-order spectral element flow simulations
field_dirichlet_vector.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, bc_list_free
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
47  use utils, only: neko_error
48  use field_list, only : field_list_t
49  implicit none
50  private
51 
53  ! for the application on a vector field.
54  type, public, extends(bc_t) :: field_dirichlet_vector_t
55  ! The bc for the first compoent.
56  type(field_dirichlet_t) :: bc_u
57  ! The bc for the second compoent.
58  type(field_dirichlet_t) :: bc_v
59  ! The bc for the third compoent.
60  type(field_dirichlet_t) :: bc_w
64  type(bc_list_t) :: bc_list
67  procedure(field_dirichlet_update), nopass, pointer :: update => null()
68  contains
70  procedure, pass(this) :: init_field => field_dirichlet_vector_init
72  procedure, pass(this) :: free => field_dirichlet_vector_free
74  procedure, pass(this) :: apply_scalar => field_dirichlet_vector_apply_scalar
76  procedure, pass(this) :: apply_vector => field_dirichlet_vector_apply_vector
78  procedure, pass(this) :: apply_vector_dev => &
81  procedure, pass(this) :: apply_scalar_dev => &
84 
85 contains
86 
88  subroutine field_dirichlet_vector_init(this, bc_name)
89  class(field_dirichlet_vector_t), intent(inout) :: this
90  character(len=*), intent(in) :: bc_name
91 
92  call neko_error("Fields must be initialized individually!")
93 
94  end subroutine field_dirichlet_vector_init
95 
98  subroutine field_dirichlet_vector_free(this)
99  class(field_dirichlet_vector_t), target, intent(inout) :: this
100 
101  call this%bc_u%free()
102  call this%bc_v%free()
103  call this%bc_w%free()
104 
105  call this%field_list%free()
106  call bc_list_free(this%bc_list)
107 
108  if (associated(this%update)) then
109  nullify(this%update)
110  end if
111 
112 
113  end subroutine field_dirichlet_vector_free
114 
120  subroutine field_dirichlet_vector_apply_scalar(this, x, n, t, tstep)
121  class(field_dirichlet_vector_t), intent(inout) :: this
122  integer, intent(in) :: n
123  real(kind=rp), intent(inout), dimension(n) :: x
124  real(kind=rp), intent(in), optional :: t
125  integer, intent(in), optional :: tstep
126 
127  call neko_error("field_dirichlet_vector cannot apply scalar BCs.&
128 &Use field_dirichlet instead!")
129 
131 
136  subroutine field_dirichlet_vector_apply_scalar_dev(this, x_d, t, tstep)
137  class(field_dirichlet_vector_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  call neko_error("field_dirichlet_vector cannot apply scalar BCs.&
143 &Use field_dirichlet instead!")
144 
146 
154  subroutine field_dirichlet_vector_apply_vector(this, x, y, z, n, t, tstep)
155  class(field_dirichlet_vector_t), intent(inout) :: this
156  integer, intent(in) :: n
157  real(kind=rp), intent(inout), dimension(n) :: x
158  real(kind=rp), intent(inout), dimension(n) :: y
159  real(kind=rp), intent(inout), dimension(n) :: z
160  real(kind=rp), intent(in), optional :: t
161  integer, intent(in), optional :: tstep
162 
163  if (present(t) .and. present(tstep)) then
164  call this%bc_u%apply_scalar(x, n, t, tstep)
165  call this%bc_v%apply_scalar(y, n, t, tstep)
166  call this%bc_w%apply_scalar(z, n, t, tstep)
167  else if (present(t)) then
168  call this%bc_u%apply_scalar(x, n, t=t)
169  call this%bc_v%apply_scalar(y, n, t=t)
170  call this%bc_w%apply_scalar(z, n, t=t)
171  else if (present(tstep)) then
172  call this%bc_u%apply_scalar(x, n, tstep=tstep)
173  call this%bc_v%apply_scalar(y, n, tstep=tstep)
174  call this%bc_w%apply_scalar(z, n, tstep=tstep)
175  end if
176 
178 
185  subroutine field_dirichlet_vector_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
186  class(field_dirichlet_vector_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  if (present(t) .and. present(tstep)) then
194  call this%bc_u%apply_scalar_dev(x_d, t, tstep)
195  call this%bc_v%apply_scalar_dev(y_d, t, tstep)
196  call this%bc_w%apply_scalar_dev(z_d, t, tstep)
197  else if (present(t)) then
198  call this%bc_u%apply_scalar_dev(x_d, t=t)
199  call this%bc_v%apply_scalar_dev(y_d, t=t)
200  call this%bc_w%apply_scalar_dev(z_d, t=t)
201  else if (present(tstep)) then
202  call this%bc_u%apply_scalar_dev(x_d, tstep=tstep)
203  call this%bc_v%apply_scalar_dev(y_d, tstep=tstep)
204  call this%bc_w%apply_scalar_dev(z_d, tstep=tstep)
205  end if
206 
208 
209 end module field_dirichlet_vector
Abstract interface defining a dirichlet condition on a list of fields.
Defines a boundary condition.
Definition: bc.f90:34
subroutine, public bc_list_free(bclst)
Destructor for a list of boundary conditions.
Definition: bc.f90:491
Coefficients.
Definition: coef.f90:34
subroutine, public device_masked_copy(a_d, b_d, mask_d, n, m)
Copy a masked vector .
Definition: device_math.F90:91
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_vector_apply_scalar(this, x, n, t, tstep)
Apply scalar by performing a masked copy.
subroutine field_dirichlet_vector_init(this, bc_name)
Initializes thisfield_bc.
subroutine field_dirichlet_vector_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
(No-op) Apply vector (device).
subroutine field_dirichlet_vector_apply_scalar_dev(this, x_d, t, tstep)
Apply scalar (device).
subroutine field_dirichlet_vector_apply_vector(this, x, y, z, n, t, tstep)
(No-op) Apply vector.
subroutine field_dirichlet_vector_free(this)
Destructor. Currently unused as is, all field_dirichlet attributes are freed in fluid_scheme::free.
Defines inflow dirichlet conditions.
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:258
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:136
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....
Extension of the user defined dirichlet condition field_dirichlet
field_list_t, To be able to group fields together
Definition: field_list.f90:13