Neko  0.9.0
A portable framework for high-order spectral element flow simulations
filters_cpu.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 !
35  use num_types, only: rp
36  implicit none
37  private
39 
40 contains
41 
43  subroutine smooth_step_cpu(x, edge0, edge1, n)
44  integer, intent(in) :: n
45  real(kind=rp), dimension(n), intent(inout) :: x
46  real(kind=rp), intent(in) :: edge0, edge1
47 
48  x = smooth_step_kernel(x, edge0, edge1)
49 
50  end subroutine smooth_step_cpu
51 
53  subroutine step_function_cpu(x, x_step, value0, value1, n)
54  integer, intent(in) :: n
55  real(kind=rp), dimension(n), intent(inout) :: x
56  real(kind=rp), intent(in) :: x_step, value0, value1
57 
58  x = step_function_kernel(x, x_step, value0, value1)
59 
60  end subroutine step_function_cpu
61 
63  subroutine permeability_cpu(x, k_0, k_1, q, n)
64  integer, intent(in) :: n
65  real(kind=rp), dimension(n), intent(inout) :: x
66  real(kind=rp), intent(in) :: k_0, k_1, q
67 
68  x = permeability_kernel(x, k_0, k_1, q)
69 
70  end subroutine permeability_cpu
71 
72 
73  ! ========================================================================== !
74  ! Internal functions and subroutines
75  ! ========================================================================== !
76 
78  elemental function smooth_step_kernel(x, edge0, edge1) result(res)
79  real(kind=rp), intent(in) :: x
80  real(kind=rp), intent(in) :: edge0, edge1
81  real(kind=rp) :: res, t
82 
83  t = clamp_kernel((x - edge0) / (edge1 - edge0), 0.0_rp, 1.0_rp)
84 
85  res = t**3 * (t * (6.0_rp * t - 15.0_rp) + 10.0_rp)
86 
87  end function smooth_step_kernel
88 
90  elemental function clamp_kernel(x, lowerlimit, upperlimit) result(res)
91  real(kind=rp), intent(in) :: x
92  real(kind=rp), intent(in) :: lowerlimit, upperlimit
93  real(kind=rp) :: res
94 
95  res = max(lowerlimit, min(upperlimit, x))
96  end function clamp_kernel
97 
99  elemental function step_function_kernel(x, x_step, value0, value1) result(res)
100  real(kind=rp), intent(in) :: x, x_step, value0, value1
101  real(kind=rp) :: res
102 
103  res = merge(value0, value1, x > x_step)
104 
105  end function step_function_kernel
106 
108  elemental function permeability_kernel(x, k_0, k_1, q) result(perm)
109  real(kind=rp), intent(in) :: x, k_0, k_1, q
110  real(kind=rp) :: perm
111 
112  perm = k_0 + (k_1 - k_0) * x * (q + 1.0_rp) / (q + x)
113 
114  end function permeability_kernel
115 
116 
117 end module filters_cpu
CPU implementations of the filter functions.
Definition: filters_cpu.f90:34
subroutine, public permeability_cpu(x, k_0, k_1, q, n)
Apply a permeability function to a scalar.
Definition: filters_cpu.f90:64
elemental real(kind=rp) function smooth_step_kernel(x, edge0, edge1)
Apply a smooth step function to a scalar.
Definition: filters_cpu.f90:79
elemental real(kind=rp) function step_function_kernel(x, x_step, value0, value1)
Apply a step function to a scalar.
subroutine, public step_function_cpu(x, x_step, value0, value1, n)
Apply a step function to a scalar.
Definition: filters_cpu.f90:54
elemental real(kind=rp) function clamp_kernel(x, lowerlimit, upperlimit)
Clamp a value between two limits.
Definition: filters_cpu.f90:91
elemental real(kind=rp) function permeability_kernel(x, k_0, k_1, q)
Apply a permeability function to a scalar.
subroutine, public smooth_step_cpu(x, edge0, edge1, n)
Apply a smooth step function to a scalar.
Definition: filters_cpu.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
#define max(a, b)
Definition: tensor.cu:40