Loading [MathJax]/jax/input/TeX/config.js
Neko 0.9.1
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
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!
34module neumann
35 use num_types, only : rp
36 use bc, only : bc_t
37 use, intrinsic :: iso_c_binding, only : c_ptr
39 use coefs, only : coef_t
40 use math, only : cfill, copy
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
57 procedure, pass(this) :: flux => neumann_flux
59 procedure, pass(this) :: set_flux_scalar => neumann_set_flux_scalar
61 procedure, pass(this) :: set_flux_array => neumann_set_flux_array
63 generic :: set_flux => set_flux_scalar, set_flux_array
65 procedure, pass(this) :: free => neumann_free
66 end type neumann_t
67
68contains
69
72 subroutine neumann_apply_scalar(this, x, n, t, tstep)
73 class(neumann_t), intent(inout) :: this
74 integer, intent(in) :: n
75 real(kind=rp), intent(inout), dimension(n) :: x
76 real(kind=rp), intent(in), optional :: t
77 integer, intent(in), optional :: tstep
78 integer :: i, m, k, facet
79 ! Store non-linear index
80 integer :: idx(4)
81
82 m = this%msk(0)
83 do i = 1, m
84 k = this%msk(i)
85 facet = this%facet(i)
86 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx,&
87 this%coef%Xh%lx)
88 select case (facet)
89 case (1,2)
90 x(k) = x(k) + this%flux_(i)*this%coef%area(idx(2), idx(3), facet, &
91 idx(4))
92 case (3,4)
93 x(k) = x(k) + this%flux_(i)*this%coef%area(idx(1), idx(3), facet, &
94 idx(4))
95 case (5,6)
96 x(k) = x(k) + this%flux_(i)*this%coef%area(idx(1), idx(2), facet, &
97 idx(4))
98 end select
99 end do
100 end subroutine neumann_apply_scalar
101
104 subroutine neumann_apply_vector(this, x, y, z, n, t, tstep)
105 class(neumann_t), intent(inout) :: this
106 integer, intent(in) :: n
107 real(kind=rp), intent(inout), dimension(n) :: x
108 real(kind=rp), intent(inout), dimension(n) :: y
109 real(kind=rp), intent(inout), dimension(n) :: z
110 real(kind=rp), intent(in), optional :: t
111 integer, intent(in), optional :: tstep
112
113 call neko_error("Neumann bc not implemented for vectors")
114
115 end subroutine neumann_apply_vector
116
119 subroutine neumann_apply_scalar_dev(this, x_d, t, tstep)
120 class(neumann_t), intent(inout), target :: this
121 type(c_ptr) :: x_d
122 real(kind=rp), intent(in), optional :: t
123 integer, intent(in), optional :: tstep
124
125 call neko_error("Neumann bc not implemented on the device")
126
127 end subroutine neumann_apply_scalar_dev
128
131 subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
132 class(neumann_t), intent(inout), target :: this
133 type(c_ptr) :: x_d
134 type(c_ptr) :: y_d
135 type(c_ptr) :: z_d
136 real(kind=rp), intent(in), optional :: t
137 integer, intent(in), optional :: tstep
138
139 call neko_error("Neumann bc not implemented on the device")
140
141 end subroutine neumann_apply_vector_dev
142
144 subroutine neumann_free(this)
145 class(neumann_t), target, intent(inout) :: this
146
147 call this%free_base()
148
149 end subroutine neumann_free
150
152 pure function neumann_flux(this) result(flux)
153 class(neumann_t), intent(in) :: this
154 real(kind=rp) :: flux(this%msk(0))
155
156 flux = this%flux_
157 end function neumann_flux
158
160 subroutine neumann_set_flux_scalar(this, flux)
161 class(neumann_t), intent(inout) :: this
162 real(kind=rp), intent(in) :: flux
163
164 this%flux_ = flux
165 end subroutine neumann_set_flux_scalar
166
169 subroutine neumann_set_flux_array(this, flux)
170 class(neumann_t), intent(inout) :: this
171 real(kind=rp), intent(in) :: flux(this%msk(0))
172
173 call copy(this%flux_, flux, this%msk(0))
174 end subroutine neumann_set_flux_array
175
178 subroutine neumann_finalize_neumann(this, flux)
179 class(neumann_t), intent(inout) :: this
180 real(kind=rp), intent(in) :: flux
181
182 call this%finalize()
183 allocate(this%flux_(this%msk(0)))
184
185 call cfill(this%flux_, flux, this%msk(0))
186 end subroutine neumann_finalize_neumann
187
188end 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:348
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:239
Defines a Neumann boundary condition.
Definition neumann.f90:34
subroutine neumann_set_flux_array(this, flux)
Set the flux using an array of values.
Definition neumann.f90:170
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:105
subroutine neumann_apply_scalar(this, x, n, t, tstep)
Boundary condition apply for a generic Neumann condition to a vector x.
Definition neumann.f90:73
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:120
subroutine neumann_set_flux_scalar(this, flux)
Set the flux using a scalar.
Definition neumann.f90:161
pure real(kind=rp) function, dimension(this%msk(0)) neumann_flux(this)
Get the flux.
Definition neumann.f90:153
subroutine neumann_finalize_neumann(this, flux)
Finalize by setting the flux.
Definition neumann.f90:179
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:132
subroutine neumann_free(this)
Destructor.
Definition neumann.f90:145
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
A Neumann boundary condition for scalar fields. Sets the flux of the field to the chosen value.
Definition neumann.f90:48