Neko  0.8.99
A portable framework for high-order spectral element flow simulations
dirichlet.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2021, 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 !
34 module dirichlet
36  use num_types, only : rp
37  use bc, only : bc_t
38  use, intrinsic :: iso_c_binding, only : c_ptr
39  implicit none
40  private
41 
44  type, public, extends(bc_t) :: dirichlet_t
45  real(kind=rp), private :: g
46  contains
47  procedure, pass(this) :: apply_scalar => dirichlet_apply_scalar
48  procedure, pass(this) :: apply_vector => dirichlet_apply_vector
49  procedure, pass(this) :: apply_scalar_dev => dirichlet_apply_scalar_dev
50  procedure, pass(this) :: apply_vector_dev => dirichlet_apply_vector_dev
51  procedure, pass(this) :: set_g => dirichlet_set_g
53  procedure, pass(this) :: free => dirichlet_free
54  end type dirichlet_t
55 
56 contains
57 
60  subroutine dirichlet_apply_scalar(this, x, n, t, tstep)
61  class(dirichlet_t), intent(inout) :: this
62  integer, intent(in) :: n
63  real(kind=rp), intent(inout), dimension(n) :: x
64  real(kind=rp), intent(in), optional :: t
65  integer, intent(in), optional :: tstep
66  integer :: i, m, k
67 
68  m = this%msk(0)
69  do i = 1, m
70  k = this%msk(i)
71  x(k) = this%g
72  end do
73  end subroutine dirichlet_apply_scalar
74 
77  subroutine dirichlet_apply_vector(this, x, y, z, n, t, tstep)
78  class(dirichlet_t), intent(inout) :: this
79  integer, intent(in) :: n
80  real(kind=rp), intent(inout), dimension(n) :: x
81  real(kind=rp), intent(inout), dimension(n) :: y
82  real(kind=rp), intent(inout), dimension(n) :: z
83  real(kind=rp), intent(in), optional :: t
84  integer, intent(in), optional :: tstep
85  integer :: i, m, k
86 
87  m = this%msk(0)
88  do i = 1, m
89  k = this%msk(i)
90  x(k) = this%g
91  y(k) = this%g
92  z(k) = this%g
93  end do
94 
95  end subroutine dirichlet_apply_vector
96 
99  subroutine dirichlet_apply_scalar_dev(this, x_d, t, tstep)
100  class(dirichlet_t), intent(inout), target :: this
101  type(c_ptr) :: x_d
102  real(kind=rp), intent(in), optional :: t
103  integer, intent(in), optional :: tstep
104 
105  call device_dirichlet_apply_scalar(this%msk_d, x_d, &
106  this%g, size(this%msk))
107 
108  end subroutine dirichlet_apply_scalar_dev
109 
112  subroutine dirichlet_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
113  class(dirichlet_t), intent(inout), target :: this
114  type(c_ptr) :: x_d
115  type(c_ptr) :: y_d
116  type(c_ptr) :: z_d
117  real(kind=rp), intent(in), optional :: t
118  integer, intent(in), optional :: tstep
119 
120  call device_dirichlet_apply_vector(this%msk_d, x_d, y_d, z_d, &
121  this%g, size(this%msk))
122 
123  end subroutine dirichlet_apply_vector_dev
124 
126  subroutine dirichlet_set_g(this, g)
127  class(dirichlet_t), intent(inout) :: this
128  real(kind=rp), intent(in) :: g
129 
130  this%g = g
131 
132  end subroutine dirichlet_set_g
133 
135  subroutine dirichlet_free(this)
136  class(dirichlet_t), target, intent(inout) :: this
137 
138  call this%free_base
139 
140  end subroutine dirichlet_free
141 
142 end module dirichlet
Defines a boundary condition.
Definition: bc.f90:34
subroutine, public device_dirichlet_apply_scalar(msk, x, g, m)
subroutine, public device_dirichlet_apply_vector(msk, x, y, z, g, m)
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
subroutine dirichlet_free(this)
Destructor.
Definition: dirichlet.f90:136
subroutine dirichlet_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)
Definition: dirichlet.f90:113
subroutine dirichlet_apply_scalar_dev(this, x_d, t, tstep)
Boundary condition apply for a generic Dirichlet condition to a vector x (device version)
Definition: dirichlet.f90:100
subroutine dirichlet_set_g(this, g)
Set value of .
Definition: dirichlet.f90:127
subroutine dirichlet_apply_scalar(this, x, n, t, tstep)
Boundary condition apply for a generic Dirichlet condition to a vector x.
Definition: dirichlet.f90:61
subroutine dirichlet_apply_vector(this, x, y, z, n, t, tstep)
Boundary condition apply for a generic Dirichlet condition to vectors x, y and z.
Definition: dirichlet.f90:78
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Base type for a boundary condition.
Definition: bc.f90:51
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44