Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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 filter, only : filter_t
38 use math, only : rzero, rone, copy
39 use field, only : field_t
40 use coefs, only : coef_t
41 use utils, only : neko_error
43 use json_module, only : json_file
45 use speclib, only : zwgll, legendre_poly
46 use matrix, only : matrix_t
47 use mxm_wrapper, only : mxm
48 use tensor, only : tnsr3d, trsp
50 use device_math, only : device_cfill
51 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
52 implicit none
53 private
54
56 type, public, extends(filter_t) :: elementwise_filter_t
59 character(len=:), allocatable :: filter_type
61 integer :: nx
63 integer :: nt
65 real(kind=rp), allocatable :: fh(:,:), fht(:,:)
66 type(c_ptr) :: fh_d = c_null_ptr
67 type(c_ptr) :: fht_d = c_null_ptr
69 real(kind=rp), allocatable :: transfer(:)
70 contains
72 procedure, pass(this) :: init => elementwise_filter_init_from_json
74 procedure, pass(this) :: init_from_components => &
77 procedure, pass(this) :: free => elementwise_filter_free
79 procedure, pass(this) :: build_1d
81 procedure, pass(this) :: apply => elementwise_field_filter_3d
83
84contains
86 subroutine elementwise_filter_init_from_json(this, json, coef)
87 class(elementwise_filter_t), intent(inout) :: this
88 type(json_file), intent(inout) :: json
89 type(coef_t), intent(in), target :: coef
90 real(kind=rp), allocatable :: transfer(:)
91 character(len=:), allocatable :: filter_type
92
93 call json_get_or_default(json, "elementwise_filter_type", &
94 filter_type, "nonBoyd")
95
96 if (json%valid_path('transfer_function')) then
97 call json_get(json, 'transfer_function', transfer)
98 end if
99
100 if (allocated(transfer)) then
101 call this%init_from_components(coef, filter_type, transfer)
102 else
103 call this%init_from_components(coef, filter_type)
104 end if
105
106
107 if (allocated(transfer)) then
108 deallocate(transfer)
109 end if
110
111 if (allocated(filter_type)) then
112 deallocate(filter_type)
113 end if
114
116
121 subroutine elementwise_filter_init_from_components(this, coef, filter_type, &
122 transfer)
123 class(elementwise_filter_t), intent(inout) :: this
124 type(coef_t), intent(in), target :: coef
125 character(len=*), intent(in) :: filter_type
126 real(kind=rp), intent(in), optional :: transfer(:)
127 integer :: nx
128
129 call this%free()
130
131 ! Filter assumes lx = ly = lz
132 call this%init_base(coef)
133 nx = coef%dof%xh%lx
134 this%nx = nx
135 this%nt = nx ! initialize as if nothing is filtered yet
136 this%filter_type = filter_type
137
138 allocate(this%fh(nx, nx))
139 allocate(this%fht(nx, nx))
140 allocate(this%transfer(nx))
141
142 call rzero(this%fh, nx*nx)
143 call rzero(this%fht, nx*nx)
144 call rone(this%transfer, nx) ! initialize as if nothing is filtered yet
145
146 if (present(transfer)) then
147 if (size(transfer) .eq. nx) then
148 this%transfer = transfer
149 else
150 call neko_error("The transfer function of the elementwise " // &
151 "filter must correspond the order of the polynomial")
152 end if
153 end if
154
155 if (neko_bcknd_device .eq. 1) then
156 call device_map(this%fh, this%fh_d, this%nx * this%nx)
157 call device_map(this%fht, this%fht_d, this%nx * this%nx)
158 call device_cfill(this%fh_d, 0.0_rp, this%nx * this%nx)
159 call device_cfill(this%fht_d, 0.0_rp, this%nx * this%nx)
160 end if
161
162 call this%build_1d()
163
165
167 subroutine elementwise_filter_free(this)
168 class(elementwise_filter_t), intent(inout) :: this
169
170 if (allocated(this%filter_type)) then
171 deallocate(this%filter_type)
172 end if
173
174 if (allocated(this%fh)) then
175 if (neko_bcknd_device .eq. 1) then
176 call device_unmap(this%fh, this%fh_d)
177 end if
178 deallocate(this%fh)
179 end if
180
181 if (allocated(this%fht)) then
182 if (neko_bcknd_device .eq. 1) then
183 call device_unmap(this%fht, this%fht_d)
184 end if
185 deallocate(this%fht)
186 end if
187
188 if (allocated(this%transfer)) then
189 deallocate(this%transfer)
190 end if
191
192 this%filter_type = ""
193 this%nx = 0
194 this%nt = 0
195
196 call this%free_base()
197
198 end subroutine elementwise_filter_free
199
201 subroutine build_1d(this)
202 class(elementwise_filter_t), intent(inout) :: this
203
204 call build_1d_cpu(this%fh, this%fht, this%transfer, &
205 this%nx, this%filter_type)
206 if (neko_bcknd_device .eq. 1) then
207 call device_memcpy(this%fh, this%fh_d, &
208 this%nx * this%nx, host_to_device, sync = .false.)
209 call device_memcpy(this%fht, this%fht_d, &
210 this%nx * this%nx, host_to_device, sync = .false.)
211 end if
212
213 end subroutine build_1d
214
216 subroutine elementwise_field_filter_3d(this, F_out, F_in)
217 class(elementwise_filter_t), intent(inout) :: this
218 type(field_t), intent(inout) :: F_out
219 type(field_t), intent(in) :: F_in
220
221 ! F_out = fh x fh x fh x F_in
222 call tnsr3d(f_out%x, this%nx, f_in%x, this%nx, this%fh, this%fht, &
223 this%fht, this%coef%msh%nelv)
224
225 end subroutine elementwise_field_filter_3d
226
234 subroutine build_1d_cpu(fh, fht, transfer, nx, filter_type)
235 integer, intent(in) :: nx
236 real(kind=rp), intent(inout) :: fh(nx, nx), fht(nx, nx)
237 real(kind=rp), intent(in) :: transfer(nx)
238 real(kind=rp) :: diag(nx, nx), rmult(nx), lj(nx), zpts(nx)
239 type(matrix_t) :: phi, pht
240 integer :: n, i, j, k
241 real(kind=rp) :: z
242 character(len=*), intent(in) :: filter_type
243
244 call phi%init(nx, nx)
245 call pht%init(nx, nx)
246
247 call zwgll(zpts, rmult, nx)
248
249 n = nx-1
250 do j = 1, nx
251 z = zpts(j)
252 call legendre_poly(lj, z, n)
253 select case (filter_type)
254 case ("Boyd")
255 pht%x(1,j) = lj(1)
256 pht%x(2,j) = lj(2)
257 do k = 3, nx
258 pht%x(k, j) = lj(k) - lj(k - 2)
259 end do
260 case ("nonBoyd")
261 pht%x(:,j) = lj
262 end select
263 end do
264
265 call trsp(phi%x, nx, pht%x, nx)
266 pht%x = phi%x
267
268 call pht%inverse(0) ! "0" for cpu implementation
269
270 diag = 0.0_rp
271
272 do i = 1, nx
273 diag(i, i) = transfer(i)
274 end do
275
276 call mxm (diag, nx, pht%x, nx, fh, nx) ! -1
277 call mxm (phi%x, nx, fh, nx, pht%x, nx) ! V D V
278
279 call copy (fh, pht%x, nx*nx)
280 call trsp (fht, nx, fh, nx)
281
282 end subroutine build_1d_cpu
283
284end module elementwise_filter
Map a Fortran array to a device (allocate and associate)
Definition device.F90:78
Copy data between host and device (or device and device)
Definition device.F90:72
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:84
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Coefficients.
Definition coef.f90:34
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:48
Implements elementwise_filter_t.
subroutine build_1d_cpu(fh, fht, transfer, 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...
subroutine elementwise_filter_init_from_components(this, coef, filter_type, transfer)
Actual Constructor.
subroutine elementwise_filter_init_from_json(this, json, coef)
Constructor.
subroutine elementwise_filter_free(this)
Destructor.
subroutine build_1d(this)
Build the 1d filter for an element.
subroutine elementwise_field_filter_3d(this, f_out, f_in)
Filter a 3D field.
Defines a field.
Definition field.f90:34
Filter to be applied to a scalar field.
Definition filter.f90:38
Utilities for retrieving parameters from the case files.
Definition math.f90:60
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:275
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:289
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:233
Defines a matrix.
Definition matrix.f90:34
Wrapper for all matrix-matrix product implementations.
subroutine, public mxm(a, n1, b, n2, c, n3)
Compute matrix-matrix product for contiguously packed matrices A,B, and C.
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
LIBRARY ROUTINES FOR SPECTRAL METHODS.
Definition speclib.f90:148
subroutine zwgll(z, w, np)
Generate NP Gauss-Lobatto Legendre points (Z) and weights (W) associated with Jacobi polynomial P(N)(...
Definition speclib.f90:179
subroutine legendre_poly(l, x, n)
Evaluate Legendre polynomials of degrees 0-N at point x and store in array L.
Definition speclib.f90:1003
Tensor operations.
Definition tensor.f90:61
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
Definition tensor.f90:124
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
Definition tensor.f90:234
Utilities.
Definition utils.f90:35
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Implements the elementwise filter for SEM.
Base abstract class for filter.
Definition filter.f90:47