Neko  0.8.1
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
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  implicit none
49  private
50 
52  ! for the application on a vector field.
53  type, public, extends(bc_t) :: field_dirichlet_vector_t
54  type(field_dirichlet_t) :: field_dirichlet_u
55  type(field_dirichlet_t) :: field_dirichlet_v
56  type(field_dirichlet_t) :: field_dirichlet_w
57  contains
59  procedure, pass(this) :: init_field => field_dirichlet_vector_init
61  procedure, pass(this) :: apply_scalar => field_dirichlet_vector_apply_scalar
63  procedure, pass(this) :: apply_vector => field_dirichlet_vector_apply_vector
65  procedure, pass(this) :: apply_vector_dev => &
68  procedure, pass(this) :: apply_scalar_dev => &
71 
72 contains
73 
75  subroutine field_dirichlet_vector_init(this, bc_name)
76  class(field_dirichlet_vector_t), intent(inout) :: this
77  character(len=*), intent(in) :: bc_name
78 
79  call neko_error("Fields must be initialized individually!")
80 
81  end subroutine field_dirichlet_vector_init
82 
85  subroutine field_dirichlet_vector_free(this)
86  type(field_dirichlet_vector_t), intent(inout) :: this
87 
88  call this%field_dirichlet_u%free()
89  call this%field_dirichlet_v%free()
90  call this%field_dirichlet_w%free()
91 
92  end subroutine field_dirichlet_vector_free
93 
99  subroutine field_dirichlet_vector_apply_scalar(this, x, n, t, tstep)
100  class(field_dirichlet_vector_t), intent(inout) :: this
101  integer, intent(in) :: n
102  real(kind=rp), intent(inout), dimension(n) :: x
103  real(kind=rp), intent(in), optional :: t
104  integer, intent(in), optional :: tstep
105 
106  call neko_error("field_dirichlet_vector cannot apply scalar BCs.&
107 &Use field_dirichlet instead!")
108 
110 
115  subroutine field_dirichlet_vector_apply_scalar_dev(this, x_d, t, tstep)
116  class(field_dirichlet_vector_t), intent(inout), target :: this
117  type(c_ptr) :: x_d
118  real(kind=rp), intent(in), optional :: t
119  integer, intent(in), optional :: tstep
120 
121  call neko_error("field_dirichlet_vector cannot apply scalar BCs.&
122 &Use field_dirichlet instead!")
123 
125 
133  subroutine field_dirichlet_vector_apply_vector(this, x, y, z, n, t, tstep)
134  class(field_dirichlet_vector_t), intent(inout) :: this
135  integer, intent(in) :: n
136  real(kind=rp), intent(inout), dimension(n) :: x
137  real(kind=rp), intent(inout), dimension(n) :: y
138  real(kind=rp), intent(inout), dimension(n) :: z
139  real(kind=rp), intent(in), optional :: t
140  integer, intent(in), optional :: tstep
141 
142  if (present(t) .and. present(tstep)) then
143  call this%field_dirichlet_u%apply_scalar(x, n, t, tstep)
144  call this%field_dirichlet_v%apply_scalar(y, n, t, tstep)
145  call this%field_dirichlet_w%apply_scalar(z, n, t, tstep)
146  else if (present(t)) then
147  call this%field_dirichlet_u%apply_scalar(x, n, t=t)
148  call this%field_dirichlet_v%apply_scalar(y, n, t=t)
149  call this%field_dirichlet_w%apply_scalar(z, n, t=t)
150  else if (present(tstep)) then
151  call this%field_dirichlet_u%apply_scalar(x, n, tstep=tstep)
152  call this%field_dirichlet_v%apply_scalar(y, n, tstep=tstep)
153  call this%field_dirichlet_w%apply_scalar(z, n, tstep=tstep)
154  end if
155 
157 
164  subroutine field_dirichlet_vector_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
165  class(field_dirichlet_vector_t), intent(inout), target :: this
166  type(c_ptr) :: x_d
167  type(c_ptr) :: y_d
168  type(c_ptr) :: z_d
169  real(kind=rp), intent(in), optional :: t
170  integer, intent(in), optional :: tstep
171 
172  if (present(t) .and. present(tstep)) then
173  call this%field_dirichlet_u%apply_scalar_dev(x_d, t, tstep)
174  call this%field_dirichlet_v%apply_scalar_dev(y_d, t, tstep)
175  call this%field_dirichlet_w%apply_scalar_dev(z_d, t, tstep)
176  else if (present(t)) then
177  call this%field_dirichlet_u%apply_scalar_dev(x_d, t=t)
178  call this%field_dirichlet_v%apply_scalar_dev(y_d, t=t)
179  call this%field_dirichlet_w%apply_scalar_dev(z_d, t=t)
180  else if (present(tstep)) then
181  call this%field_dirichlet_u%apply_scalar_dev(x_d, tstep=tstep)
182  call this%field_dirichlet_v%apply_scalar_dev(y_d, tstep=tstep)
183  call this%field_dirichlet_w%apply_scalar_dev(z_d, tstep=tstep)
184  end if
185 
187 
188 end module field_dirichlet_vector
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_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: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....
Extension of the user defined dirichlet condition field_dirichlet
field_list_t, To be able to group fields together
Definition: field_list.f90:7