Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
PDE_filter.f90
Go to the documentation of this file.
1! Copyright (c) 2023-2026, 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!
35
37 use num_types, only : rp
38 use json_module, only : json_file
39 use registry, only : neko_registry
40 use field, only : field_t
41 use coefs, only : coef_t
42 use ax_product, only : ax_t, ax_helm_factory
43 use krylov, only : ksp_t, ksp_monitor_t, krylov_solver_factory
44 use precon, only : pc_t, precon_factory, precon_destroy
45 use bc_list, only : bc_list_t
46 use neumann, only : neumann_t
48 use gather_scatter, only : gs_t, gs_op_add
51 use registry, only : neko_registry
52 use filter, only : filter_t
55 use coefs, only : coef_t
56 use logger, only : neko_log, log_size
58 use dofmap, only : dofmap_t
59 use jacobi, only : jacobi_t
61 use sx_jacobi, only : sx_jacobi_t
62 use utils, only : neko_error
65 implicit none
66 private
67
72 type, public, extends(filter_t) :: pde_filter_t
73
75 class(ax_t), allocatable :: ax
77 type(ksp_monitor_t) :: ksp_results(1)
79 class(ksp_t), allocatable :: ksp_filt
81 class(pc_t), allocatable :: pc_filt
83 type(bc_list_t) :: bclst_filt
84
85 ! Inputs from the user
87 real(kind=rp) :: r
89 real(kind=rp) :: abstol_filt
91 integer :: ksp_max_iter
93 character(len=:), allocatable :: ksp_solver
94 ! > preconditioner type
95 character(len=:), allocatable :: precon_type_filt
96 integer :: ksp_n, n, i
97
98
99
100 contains
102 procedure, pass(this) :: init => pde_filter_init_from_json
104 procedure, pass(this) :: init_from_components => &
107 procedure, pass(this) :: free => pde_filter_free
109 procedure, pass(this) :: apply => pde_filter_apply
110 end type pde_filter_t
111
112contains
113
115 subroutine pde_filter_init_from_json(this, json, coef)
116 class(pde_filter_t), intent(inout) :: this
117 type(json_file), intent(inout) :: json
118 type(coef_t), target, intent(in) :: coef
119 real(kind=rp) :: r, tol
120 integer :: max_iter
121 character(len=:), allocatable :: ksp_solver, precon_type
122
123 ! user parameters
124 call json_get(json, "radius", r)
125 call json_get_or_default(json, "tolerance", tol, 1e-10_rp)
126 call json_get_or_default(json, "max_iter", max_iter, 200)
127 call json_get_or_default(json, "solver", ksp_solver, "cg")
128 call json_get_or_default(json, "preconditioner", precon_type, "jacobi")
129
130 call this%init_from_components(coef, r, tol, max_iter, ksp_solver, &
131 precon_type)
132
133 end subroutine pde_filter_init_from_json
134
136 subroutine pde_filter_init_from_components(this, coef, r, tol, max_iter, &
137 ksp_solver, precon_type)
138 class(pde_filter_t), intent(inout) :: this
139 type(coef_t), target, intent(in) :: coef
140 real(kind=rp), intent(in) :: r, tol
141 integer, intent(in) :: max_iter
142 character(len=*), intent(in) :: ksp_solver, precon_type
143 integer :: n
144
145 call this%init_base(coef)
146
147 this%r = r
148 this%abstol_filt = tol
149 this%ksp_max_iter = max_iter
150 this%ksp_solver = ksp_solver
151 this%precon_type_filt = precon_type
152
153 ! set the number of dofs
154 n = this%coef%dof%size()
155
156 ! init the bc list (all Neuman BCs, will remain empty)
157 call this%bclst_filt%init()
158
159 ! Setup backend dependent Ax routines
160 call ax_helm_factory(this%Ax, full_formulation = .false.)
161
162 ! set up krylov solver
163 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
164 this%ksp_max_iter, this%abstol_filt)
165
166 ! set up preconditioner
167 call filter_precon_factory(this%pc_filt, this%ksp_filt, &
168 this%coef, this%coef%dof, this%coef%gs_h, this%bclst_filt, &
169 this%precon_type_filt)
170
172
174 subroutine pde_filter_free(this)
175 class(pde_filter_t), intent(inout) :: this
176
177 if (allocated(this%Ax)) then
178 deallocate(this%Ax)
179 end if
180
181 if (allocated(this%ksp_filt)) then
182 call this%ksp_filt%free()
183 deallocate(this%ksp_filt)
184 end if
185
186 if (allocated(this%pc_filt)) then
187 call precon_destroy(this%pc_filt)
188 deallocate(this%pc_filt)
189 end if
190
191 if (allocated(this%ksp_solver)) then
192 deallocate(this%ksp_solver)
193 end if
194
195 if (allocated(this%precon_type_filt)) then
196 deallocate(this%precon_type_filt)
197 end if
198
199 call this%bclst_filt%free()
200
201 call this%free_base()
202
203 end subroutine pde_filter_free
204
209 subroutine pde_filter_apply(this, F_out, F_in)
210 class(pde_filter_t), intent(inout) :: this
211 type(field_t), intent(in) :: F_in
212 type(field_t), intent(inout) :: F_out
213 integer :: n, i
214 type(field_t), pointer :: RHS, d_F_out
215 character(len=LOG_SIZE) :: log_buf
216 integer :: temp_indices(2)
217
218 n = this%coef%dof%size()
219 call neko_scratch_registry%request_field(rhs, temp_indices(1), .false.)
220 call neko_scratch_registry%request_field(d_f_out, temp_indices(2), .false.)
221 ! in a similar fasion to pressure/velocity, we will solve for d_F_out.
222
223 ! to improve convergence, we use F_in as an initial guess for F_out.
224 ! so F_out = F_in + d_F_in.
225
226 ! Defining the operator A = -r^2 \nabla^2 + I
227 ! the system changes from:
228 ! A (F_out) = F_in
229 ! to
230 ! A (d_F_out) = F_in - A(F_in)
231
232 ! set up Helmholtz operators and RHS
233 if (neko_bcknd_device .eq. 1) then
234 call device_cfill(this%coef%h1_d, (this%r / (2.0_rp * sqrt(3.0_rp)))**2, n)
235 call device_cfill(this%coef%h2_d, 1.0_rp, n)
236 else
237 ! h1 is already negative in its definition
238 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
239 ! ax_helm includes the mass matrix in h2
240 this%coef%h2 = 1.0_rp
241 end if
242 this%coef%ifh2 = .true.
243
244 ! compute the A(F_in) component of the RHS
245 ! (note, to be safe with the inout intent we first copy F_in to the
246 ! temporary d_F_out)
247 call field_copy(d_f_out, f_in)
248 call this%Ax%compute(rhs%x, d_f_out%x, this%coef, this%coef%msh, &
249 this%coef%Xh)
250
251 if (neko_bcknd_device .eq. 1) then
252 call device_subcol3(rhs%x_d, f_in%x_d, this%coef%B_d, n)
253 call device_cmult(rhs%x_d, -1.0_rp, n)
254 else
255 do i = 1, n
256 ! mass matrix should be included here
257 rhs%x(i,1,1,1) = f_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
258 - rhs%x(i,1,1,1)
259 end do
260 end if
261
262 ! gather scatter
263 call this%coef%gs_h%op(rhs, gs_op_add)
264
265 ! set BCs
266 call this%bclst_filt%apply_scalar(rhs%x, n)
267
268 ! Solve Helmholtz equation
269 call profiler_start_region("filter solve")
270 this%ksp_results(1) = &
271 this%ksp_filt%solve(this%Ax, d_f_out, rhs%x, n, this%coef, &
272 this%bclst_filt, this%coef%gs_h)
273
275
276 ! add result
277 call field_add3(f_out, f_in, d_f_out)
278 ! update preconditioner (needed?)
279 call this%pc_filt%update()
280
281 ! write it all out
282 call neko_log%section('PDE Filter')
283
284 write(log_buf, '(A,A,A)') 'Iterations: ', 'Start residual: ', &
285 'Final residual:'
286 call neko_log%message(log_buf)
287 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
288 this%ksp_results%res_start, this%ksp_results%res_final
289 call neko_log%message(log_buf)
290 call neko_log%end_section()
291
292 call neko_scratch_registry%relinquish_field(temp_indices)
293
294 end subroutine pde_filter_apply
295
297 subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype)
298 class(pc_t), allocatable, target, intent(inout) :: pc
299 class(ksp_t), target, intent(inout) :: ksp
300 type(coef_t), target, intent(in) :: coef
301 type(dofmap_t), target, intent(in) :: dof
302 type(gs_t), target, intent(inout) :: gs
303 type(bc_list_t), target, intent(inout) :: bclst
304 character(len=*) :: pctype
305
306 call precon_factory(pc, pctype)
307
308 select type (pcp => pc)
309 type is (jacobi_t)
310 call pcp%init(coef, dof, gs)
311 type is (sx_jacobi_t)
312 call pcp%init(coef, dof, gs)
313 type is (device_jacobi_t)
314 call pcp%init(coef, dof, gs)
315 end select
316
317 call ksp%set_pc(pc)
318
319 end subroutine filter_precon_factory
320
321end module pde_filter
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.
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Coefficients.
Definition coef.f90:34
Jacobi preconditioner accelerator backend.
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_subcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
subroutine, public field_copy(a, b, n)
Copy a vector .
subroutine, public field_add3(a, b, c, n)
Vector addition .
Defines a field.
Definition field.f90:34
Filter to be applied to a scalar field.
Definition filter.f90:38
Gather-scatter.
Jacobi preconditioner.
Definition pc_jacobi.f90:34
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
integer, parameter, public ksp_max_iter
Maximum number of iters.
Definition krylov.f90:51
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:64
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:66
Build configurations.
integer, parameter neko_bcknd_device
Defines a Neumann boundary condition.
Definition neumann.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
A PDE based filter.
subroutine pde_filter_free(this)
Destructor.
subroutine pde_filter_init_from_json(this, json, coef)
Constructor from json.
subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype)
Initialize a Krylov preconditioner.
subroutine pde_filter_init_from_components(this, coef, r, tol, max_iter, ksp_solver, precon_type)
Actual constructor.
subroutine pde_filter_apply(this, f_out, f_in)
Apply the filter.
Defines Pressure and velocity residuals in the Pn-Pn formulation.
Definition pnpn_res.f90:34
Krylov preconditioner.
Definition precon.f90:34
Profiling interface.
Definition profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:107
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Jacobi preconditioner SX-Aurora backend.
Utilities.
Definition utils.f90:35
Base type for a matrix-vector product providing .
Definition ax.f90:43
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Defines a jacobi preconditioner.
Base abstract class for filter.
Definition filter.f90:47
Defines a jacobi preconditioner.
Definition pc_jacobi.f90:45
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:73
A Neumann boundary condition. Sets the flux of the field to the chosen values.
Definition neumann.f90:60
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .
Abstract type to compute pressure residual.
Definition pnpn_res.f90:48
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
Defines a jacobi preconditioner for SX-Aurora.