Neko  0.8.1
A portable framework for high-order spectral element flow simulations
neumann.f90
Go to the documentation of this file.
1 ! Copyright (c) 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 !
34 module neumann
35  use num_types
36  use bc, only : bc_t
37  use, intrinsic :: iso_c_binding, only : c_ptr
38  use utils, only : neko_error, nonlinear_index
39  use coefs, only : coef_t
40  implicit none
41  private
42 
47  type, public, extends(bc_t) :: neumann_t
48  real(kind=rp), private :: flux_
49  contains
50  procedure, pass(this) :: apply_scalar => neumann_apply_scalar
51  procedure, pass(this) :: apply_vector => neumann_apply_vector
52  procedure, pass(this) :: apply_scalar_dev => neumann_apply_scalar_dev
53  procedure, pass(this) :: apply_vector_dev => neumann_apply_vector_dev
54  procedure, pass(this) :: init_neumann => neumann_init_neumann
55  procedure, pass(this) :: flux => neumann_flux
56  end type neumann_t
57 
58 contains
59 
62  subroutine neumann_apply_scalar(this, x, n, t, tstep)
63  class(neumann_t), intent(inout) :: this
64  integer, intent(in) :: n
65  real(kind=rp), intent(inout), dimension(n) :: x
66  real(kind=rp), intent(in), optional :: t
67  integer, intent(in), optional :: tstep
68  integer :: i, m, k, facet
69  ! Store non-linear index
70  integer :: idx(4)
71 
72  m = this%msk(0)
73  do i = 1, m
74  k = this%msk(i)
75  facet = this%facet(i)
76  idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx,&
77  this%coef%Xh%lx)
78  select case(facet)
79  case(1,2)
80  x(k) = x(k) + this%flux_*this%coef%area(idx(2), idx(3), facet, idx(4))
81  case(3,4)
82  x(k) = x(k) + this%flux_*this%coef%area(idx(1), idx(3), facet, idx(4))
83  case(5,6)
84  x(k) = x(k) + this%flux_*this%coef%area(idx(1), idx(2), facet, idx(4))
85  end select
86  end do
87  end subroutine neumann_apply_scalar
88 
91  subroutine neumann_apply_vector(this, x, y, z, n, t, tstep)
92  class(neumann_t), intent(inout) :: this
93  integer, intent(in) :: n
94  real(kind=rp), intent(inout), dimension(n) :: x
95  real(kind=rp), intent(inout), dimension(n) :: y
96  real(kind=rp), intent(inout), dimension(n) :: z
97  real(kind=rp), intent(in), optional :: t
98  integer, intent(in), optional :: tstep
99  integer :: i, m, k
100 
101  call neko_error("Neumann bc not implemented for vectors")
102 
103  end subroutine neumann_apply_vector
104 
107  subroutine neumann_apply_scalar_dev(this, x_d, t, tstep)
108  class(neumann_t), intent(inout), target :: this
109  type(c_ptr) :: x_d
110  real(kind=rp), intent(in), optional :: t
111  integer, intent(in), optional :: tstep
112 
113  call neko_error("Neumann bc not implemented on the device")
114 
115  end subroutine neumann_apply_scalar_dev
116 
119  subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
120  class(neumann_t), intent(inout), target :: this
121  type(c_ptr) :: x_d
122  type(c_ptr) :: y_d
123  type(c_ptr) :: z_d
124  real(kind=rp), intent(in), optional :: t
125  integer, intent(in), optional :: tstep
126 
127  call neko_error("Neumann bc not implemented on the device")
128 
129  end subroutine neumann_apply_vector_dev
130 
134  subroutine neumann_init_neumann(this, flux)
135  class(neumann_t), intent(inout) :: this
136  real(kind=rp), intent(in) :: flux
137 
138  this%flux_ = flux
139  end subroutine neumann_init_neumann
140 
142  pure function neumann_flux(this) result(flux)
143  class(neumann_t), intent(in) :: this
144  real(kind=rp) :: flux
145 
146  flux = this%flux_
147  end function neumann_flux
148 
149 end module neumann
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Defines a boundary condition.
Definition: bc.f90:34
Coefficients.
Definition: coef.f90:34
Defines a Neumann boundary condition.
Definition: neumann.f90:34
subroutine neumann_apply_vector(this, x, y, z, n, t, tstep)
Boundary condition apply for a generic Neumann condition to vectors x, y and z.
Definition: neumann.f90:92
subroutine neumann_apply_scalar(this, x, n, t, tstep)
Boundary condition apply for a generic Neumann condition to a vector x.
Definition: neumann.f90:63
subroutine neumann_init_neumann(this, flux)
Constructor.
Definition: neumann.f90:135
pure real(kind=rp) function neumann_flux(this)
Get the set flux.
Definition: neumann.f90:143
subroutine neumann_apply_scalar_dev(this, x_d, t, tstep)
Boundary condition apply for a generic Neumann condition to a vector x (device version)
Definition: neumann.f90:108
subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
Boundary condition apply for a generic Neumann condition to vectors x, y and z (device version)
Definition: neumann.f90:120
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35
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 Neumann boundary condition. This sets the flux of the field to the chosen value.
Definition: neumann.f90:47