Neko  0.8.99
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  use math, only : cfill
41  implicit none
42  private
43 
48  type, public, extends(bc_t) :: neumann_t
49  real(kind=rp), allocatable, private :: flux_(:)
50  contains
51  procedure, pass(this) :: apply_scalar => neumann_apply_scalar
52  procedure, pass(this) :: apply_vector => neumann_apply_vector
53  procedure, pass(this) :: apply_scalar_dev => neumann_apply_scalar_dev
54  procedure, pass(this) :: apply_vector_dev => neumann_apply_vector_dev
55  procedure, pass(this) :: finalize_neumann => neumann_finalize_neumann
56  procedure, pass(this) :: flux => neumann_flux
58  procedure, pass(this) :: free => neumann_free
59  end type neumann_t
60 
61 contains
62 
65  subroutine neumann_apply_scalar(this, x, n, t, tstep)
66  class(neumann_t), intent(inout) :: this
67  integer, intent(in) :: n
68  real(kind=rp), intent(inout), dimension(n) :: x
69  real(kind=rp), intent(in), optional :: t
70  integer, intent(in), optional :: tstep
71  integer :: i, m, k, facet
72  ! Store non-linear index
73  integer :: idx(4)
74 
75  m = this%msk(0)
76  do i = 1, m
77  k = this%msk(i)
78  facet = this%facet(i)
79  idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx,&
80  this%coef%Xh%lx)
81  select case(facet)
82  case(1,2)
83  x(k) = x(k) + this%flux_(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
84  case(3,4)
85  x(k) = x(k) + this%flux_(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
86  case(5,6)
87  x(k) = x(k) + this%flux_(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
88  end select
89  end do
90  end subroutine neumann_apply_scalar
91 
94  subroutine neumann_apply_vector(this, x, y, z, n, t, tstep)
95  class(neumann_t), intent(inout) :: this
96  integer, intent(in) :: n
97  real(kind=rp), intent(inout), dimension(n) :: x
98  real(kind=rp), intent(inout), dimension(n) :: y
99  real(kind=rp), intent(inout), dimension(n) :: z
100  real(kind=rp), intent(in), optional :: t
101  integer, intent(in), optional :: tstep
102 
103  call neko_error("Neumann bc not implemented for vectors")
104 
105  end subroutine neumann_apply_vector
106 
109  subroutine neumann_apply_scalar_dev(this, x_d, t, tstep)
110  class(neumann_t), intent(inout), target :: this
111  type(c_ptr) :: x_d
112  real(kind=rp), intent(in), optional :: t
113  integer, intent(in), optional :: tstep
114 
115  call neko_error("Neumann bc not implemented on the device")
116 
117  end subroutine neumann_apply_scalar_dev
118 
121  subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
122  class(neumann_t), intent(inout), target :: this
123  type(c_ptr) :: x_d
124  type(c_ptr) :: y_d
125  type(c_ptr) :: z_d
126  real(kind=rp), intent(in), optional :: t
127  integer, intent(in), optional :: tstep
128 
129  call neko_error("Neumann bc not implemented on the device")
130 
131  end subroutine neumann_apply_vector_dev
132 
134  subroutine neumann_free(this)
135  class(neumann_t), target, intent(inout) :: this
136 
137  call this%free_base()
138 
139  end subroutine neumann_free
140 
142  pure function neumann_flux(this) result(flux)
143  class(neumann_t), intent(in) :: this
144  real(kind=rp) :: flux(this%msk(0))
145 
146  flux = this%flux_
147  end function neumann_flux
148 
151  subroutine neumann_finalize_neumann(this, flux)
152  class(neumann_t), intent(inout) :: this
153  real(kind=rp), intent(in) :: flux
154 
155  call this%finalize()
156  allocate(this%flux_(this%msk(0)))
157 
158  call cfill(this%flux_, flux, this%msk(0))
159  end subroutine neumann_finalize_neumann
160 
161 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
Definition: math.f90:60
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition: math.f90:314
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:95
subroutine neumann_apply_scalar(this, x, n, t, tstep)
Boundary condition apply for a generic Neumann condition to a vector x.
Definition: neumann.f90:66
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:110
pure real(kind=rp) function, dimension(this%msk(0)) neumann_flux(this)
Get the flux.
Definition: neumann.f90:143
subroutine neumann_finalize_neumann(this, flux)
Finalize by setting the flux.
Definition: neumann.f90:152
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:122
subroutine neumann_free(this)
Destructor.
Definition: neumann.f90:135
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:55
Generic Neumann boundary condition. This sets the flux of the field to the chosen value.
Definition: neumann.f90:48