43 use json_module,
only : json_file
51 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
59 character(len=:),
allocatable :: filter_type
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(:)
74 procedure, pass(this) :: init_from_components => &
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
94 filter_type,
"nonBoyd")
96 if (json%valid_path(
'transfer_function'))
then
97 call json_get(json,
'transfer_function', transfer)
100 if (
allocated(transfer))
then
101 call this%init_from_components(coef, filter_type, transfer)
103 call this%init_from_components(coef, filter_type)
107 if (
allocated(transfer))
then
111 if (
allocated(filter_type))
then
112 deallocate(filter_type)
124 type(
coef_t),
intent(in),
target :: coef
125 character(len=*),
intent(in) :: filter_type
126 real(kind=
rp),
intent(in),
optional :: transfer(:)
132 call this%init_base(coef)
136 this%filter_type = filter_type
138 allocate(this%fh(nx, nx))
139 allocate(this%fht(nx, nx))
140 allocate(this%transfer(nx))
142 call rzero(this%fh, nx*nx)
143 call rzero(this%fht, nx*nx)
144 call rone(this%transfer, nx)
146 if (
present(transfer))
then
147 if (
size(transfer) .eq. nx)
then
148 this%transfer = transfer
150 call neko_error(
"The transfer function of the elementwise " // &
151 "filter must correspond the order of the polynomial")
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)
159 call device_cfill(this%fht_d, 0.0_rp, this%nx * this%nx)
170 if (
allocated(this%filter_type))
then
171 deallocate(this%filter_type)
174 if (
allocated(this%fh))
then
178 if (
allocated(this%fht))
then
182 if (
allocated(this%transfer))
then
183 deallocate(this%transfer)
186 if (c_associated(this%fh_d))
then
190 if (c_associated(this%fht_d))
then
194 this%filter_type =
""
198 call this%free_base()
207 this%nx, this%filter_type)
220 type(
field_t),
intent(inout) :: F_out
221 type(
field_t),
intent(in) :: F_in
224 call tnsr3d(f_out%x, this%nx, f_in%x, this%nx, this%fh, this%fht, &
225 this%fht, this%coef%msh%nelv)
237 integer,
intent(in) :: nx
238 real(kind=
rp),
intent(inout) :: fh(nx, nx), fht(nx, nx)
239 real(kind=
rp),
intent(in) :: transfer(nx)
240 real(kind=
rp) :: diag(nx, nx), rmult(nx), lj(nx), zpts(nx)
242 integer :: n, i, j, k
244 character(len=*),
intent(in) :: filter_type
246 call phi%init(nx, nx)
247 call pht%init(nx, nx)
249 call zwgll(zpts, rmult, nx)
255 select case (filter_type)
260 pht%x(k, j) = lj(k) - lj(k - 2)
267 call trsp(phi%x, nx, pht%x, nx)
275 diag(i, i) = transfer(i)
278 call mxm (diag, nx, pht%x, nx, fh, nx)
279 call mxm (phi%x, nx, fh, nx, pht%x, nx)
281 call copy (fh, pht%x, nx*nx)
282 call trsp (fht, nx, fh, nx)
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
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.
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
subroutine, public device_free(x_d)
Deallocate memory on the device.
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.
Filter to be applied to a scalar field.
Utilities for retrieving parameters from the case files.
subroutine, public rone(a, n)
Set all elements to one.
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public rzero(a, n)
Zero a real vector.
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.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
LIBRARY ROUTINES FOR SPECTRAL METHODS.
subroutine zwgll(z, w, np)
Generate NP Gauss-Lobatto Legendre points (Z) and weights (W) associated with Jacobi polynomial P(N)(...
subroutine legendre_poly(l, x, n)
Evaluate Legendre polynomials of degrees 0-N at point x and store in array L.
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Implements the elementwise filter for SEM.
Base abstract class for filter.