Neko  0.9.0
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
39  use num_types, only: rp
40  use utils, only: neko_error
41 
46 
47  implicit none
48 
49  private
51 
52 contains
53 
69  subroutine smooth_step_field(F, edge0, edge1)
70  type(field_t), intent(inout) :: f
71  real(kind=rp), intent(in) :: edge0, edge1
72  integer :: n
73 
74  n = f%size()
75  if (neko_bcknd_device .eq. 1) then
76  call smooth_step_device(f%x_d, edge0, edge1, n)
77  else
78  call smooth_step_cpu(f%x, edge0, edge1, n)
79  end if
80  end subroutine smooth_step_field
81 
89  subroutine permeability_field(F, k_0, k_1, q)
90  type(field_t), intent(inout) :: f
91  real(kind=rp), intent(in) :: k_0, k_1
92  real(kind=rp), intent(in) :: q
93  integer :: n
94 
95  n = f%size()
96  if (neko_bcknd_device .eq. 1) then
97  call permeability_device(f%x_d, k_0, k_1, q, n)
98  else
99  call permeability_cpu(f%x, k_0, k_1, q, n)
100  end if
101  end subroutine permeability_field
102 
108  subroutine step_function_field(F, x0, value0, value1)
109  type(field_t), intent(inout) :: f
110  real(kind=rp), intent(in) :: x0, value0, value1
111  integer :: n
112 
113  n = f%size()
114  if (neko_bcknd_device .eq. 1) then
115  call step_function_device(f%x_d, x0, value0, value1, n)
116  else
117  call step_function_cpu(f%x, x0, value0, value1, n)
118  end if
119  end subroutine step_function_field
120 
121 end module filters
Defines a field.
Definition: field.f90:34
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
subroutine, public step_function_cpu(x, x_step, value0, value1, n)
Apply a step function to a scalar.
Definition: filters_cpu.f90:54
subroutine, public smooth_step_cpu(x, edge0, edge1, n)
Apply a smooth step function to a scalar.
Definition: filters_cpu.f90:44
Device implementations of the filter functions.
subroutine permeability_device(x, k_0, k_1, q, n)
Apply a permeability function to a scalar.
subroutine smooth_step_device(x, edge0, edge1, n)
Apply a smooth step function to a scalar.
subroutine step_function_device(x, x_step, value0, value1, n)
Apply a step function to a scalar.
A module containing filter functions and subroutines. These functions are used to modify fields in a ...
Definition: filters.f90:36
subroutine, public permeability_field(F, k_0, k_1, q)
Apply a permeability function to a field.
Definition: filters.f90:90
subroutine, public smooth_step_field(F, edge0, edge1)
Apply a smooth step function to a field.
Definition: filters.f90:70
subroutine, public step_function_field(F, x0, value0, value1)
Apply a step function to a field.
Definition: filters.f90:109
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35