Neko  0.8.1
A portable framework for high-order spectral element flow simulations
filters.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 !
36 module filters
37  use field, only: field_t
38  use num_types, only: rp
39  implicit none
40 
41  private
43 
44 contains
45 
61  subroutine smooth_step_field(F, edge0, edge1)
62  use filters_cpu, only: smooth_step_cpu
63 
64  type(field_t), intent(inout) :: f
65  real(kind=rp), intent(in) :: edge0, edge1
66 
67  f%x = smooth_step_cpu(f%x, edge0, edge1)
68  end subroutine smooth_step_field
69 
77  subroutine permeability_field(F_out, x, k_0, k_1, q)
78  use filters_cpu, only: permeability_cpu
79 
80  type(field_t), intent(inout) :: f_out
81  type(field_t), intent(in) :: x
82  real(kind=rp), intent(in) :: k_0, k_1
83  real(kind=rp), intent(in) :: q
84 
85  f_out%x = permeability_cpu(x%x, k_0, k_1, q)
86  end subroutine permeability_field
87 
93  subroutine step_function_field(F, x0, value0, value1)
95 
96  type(field_t), intent(inout) :: f
97  real(kind=rp), intent(in) :: x0, value0, value1
98 
99  f%x = step_function_cpu(f%x, x0, value0, value1)
100  end subroutine step_function_field
101 
102 end module filters
Defines a field.
Definition: field.f90:34
CPU implementations of the filter functions.
Definition: filters_cpu.f90:34
elemental real(kind=rp) function smooth_step_cpu(x, edge0, edge1)
Apply a smooth step function to a scalar.
Definition: filters_cpu.f90:46
elemental real(kind=rp) function permeability_cpu(x, k_0, k_1, q)
Apply a permeability function to a scalar.
Definition: filters_cpu.f90:76
elemental real(kind=rp) function step_function_cpu(x, x_step, value0, value1)
Apply a step function to a scalar.
Definition: filters_cpu.f90:67
A module containing filter functions and subroutines. These functions are used to modify fields in a ...
Definition: filters.f90:36
subroutine, public smooth_step_field(F, edge0, edge1)
Apply a smooth step function to a field.
Definition: filters.f90:62
subroutine, public step_function_field(F, x0, value0, value1)
Apply a step function to a field.
Definition: filters.f90:94
subroutine, public permeability_field(F_out, x, k_0, k_1, q)
Apply a permeability function to a field.
Definition: filters.f90:78
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12