Neko  0.8.99
A portable framework for high-order spectral element flow simulations
elementwise_filter.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 !
33 !
36  use num_types, only : rp
37  use math
38  use field, only : field_t
39  use utils, only : neko_error
40  use neko_config, only : neko_bcknd_device
42  use tensor, only : tnsr3d
43  implicit none
44  private
45 
47  type, public :: elementwise_filter_t
50  character(len=64) :: filter_type
52  integer :: nx
54  integer :: nt
56  real(kind=rp), allocatable :: fh(:,:), fht(:,:)
58  real(kind=rp), allocatable :: trnsfr(:)
59  contains
61  procedure, pass(this) :: init => elementwise_filter_init
63  procedure, pass(this) :: free => elementwise_filter_free
65  procedure, pass(this) :: build_1d
67  procedure, pass(this) :: filter_3d => elementwise_field_filter_3d
68  end type elementwise_filter_t
69 
70 contains
74  subroutine elementwise_filter_init(this, nx, filter_type)
75  class(elementwise_filter_t), intent(inout) :: this
76  character(len=*) :: filter_type
77  integer :: nx
78 
79  this%nx = nx
80  this%nt = nx ! initialize as if nothing is filtered yet
81  this%filter_type = filter_type
82 
83  allocate(this%fh(nx, nx))
84  allocate(this%fht(nx, nx))
85  allocate(this%trnsfr(nx))
86 
87  call rzero(this%fh, nx*nx)
88  call rzero(this%fht, nx*nx)
89  call rone(this%trnsfr, nx) ! initialize as if nothing is filtered yet
90 
91  end subroutine elementwise_filter_init
92 
94  subroutine elementwise_filter_free(this)
95  class(elementwise_filter_t), intent(inout) :: this
96 
97  if (allocated(this%fh)) then
98  deallocate(this%fh)
99  end if
100 
101  if (allocated(this%fht)) then
102  deallocate(this%fht)
103  end if
104 
105  if (allocated(this%trnsfr)) then
106  deallocate(this%trnsfr)
107  end if
108 
109  this%filter_type = ""
110  this%nx = 0
111  this%nt = 0
112 
113  end subroutine elementwise_filter_free
114 
116  subroutine build_1d(this)
117  class(elementwise_filter_t), intent(inout) :: this
118 
119  if (neko_bcknd_device .eq. 1) then
120  call neko_error("build_1d not implemented on accelarators.")
121  else
122  call build_1d_cpu(this%fh, this%fht, this%trnsfr, &
123  this%nx, this%filter_type)
124  end if
125 
126  end subroutine build_1d
127 
129  subroutine elementwise_field_filter_3d(this, v, u, nelv)
130  class(elementwise_filter_t), intent(inout) :: this
131  integer, intent(inout) :: nelv
132  real(kind=rp), intent(inout), dimension(this%nx, this%nx, this%nx, nelv) :: u,v
133 
134  ! v = fh x fh x fh x u
135  call tnsr3d(v, this%nx, u, this%nx, this%fh, this%fht, this%fht, nelv)
136 
137  end subroutine elementwise_field_filter_3d
138 
139 end module elementwise_filter
Implements the CPU kernel for the elementwise_filter_t type.
subroutine, public build_1d_cpu(fh, fht, trnsfr, nx, filter_type)
Build the 1d filter for an element on the CPU. Suppose field x is filtered into x_hat by x_hat = fh*x...
Implements explicit_filter_t.
subroutine elementwise_field_filter_3d(this, v, u, nelv)
Filter a 3D field.
subroutine elementwise_filter_free(this)
Destructor.
subroutine build_1d(this)
Build the 1d filter for an element.
subroutine elementwise_filter_init(this, nx, filter_type)
Constructor.
Defines a field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public rone(a, n)
Set all elements to one.
Definition: math.f90:217
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:184
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
Tensor operations.
Definition: tensor.f90:61
subroutine, public tnsr3d(v, nv, u, nu, A, Bt, Ct, nelv)
Tensor product performed on nelv elements.
Definition: tensor.f90:223
Utilities.
Definition: utils.f90:35
Implements the explicit filter for SEM.