Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
PDE_filter.f90
Go to the documentation of this file.
1! Copyright (c) 2023, 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
41 use field, only: field_t
42 use coefs, only: coef_t
43 use ax_product, only: ax_t, ax_helm_factory
44 use krylov, only: ksp_t, ksp_monitor_t, krylov_solver_factory, &
45 krylov_solver_destroy
46 use precon, only: pc_t, precon_factory, precon_destroy
47 use bc_list, only : bc_list_t
48 use neumann, only: neumann_t
50 use gather_scatter, only: gs_t, gs_op_add
54 use filter, only: filter_t
57 use coefs, only: coef_t
58 use logger, only: neko_log, log_size
60 use dofmap, only: dofmap_t
61 use jacobi, only: jacobi_t
63 use sx_jacobi, only: sx_jacobi_t
64 use hsmg, only: hsmg_t
65 use utils, only: neko_error
67 implicit none
68 private
69
74 type, public, extends(filter_t) :: pde_filter_t
75
77 class(ax_t), allocatable :: ax
79 type(ksp_monitor_t) :: ksp_results(1)
81 class(ksp_t), allocatable :: ksp_filt
83 class(pc_t), allocatable :: pc_filt
85 type(bc_list_t) :: bclst_filt
86
87 ! Inputs from the user
89 real(kind=rp) :: r
91 real(kind=rp) :: abstol_filt
93 integer :: ksp_max_iter
95 character(len=:), allocatable :: ksp_solver
96 ! > preconditioner type
97 character(len=:), allocatable :: precon_type_filt
98 integer :: ksp_n, n, i
99
100
101
102 contains
104 procedure, pass(this) :: init => pde_filter_init_from_json
106 procedure, pass(this) :: init_from_attributes => &
109 procedure, pass(this) :: free => pde_filter_free
111 procedure, pass(this) :: apply => pde_filter_apply
112 end type pde_filter_t
113
114contains
115
117 subroutine pde_filter_init_from_json(this, json, coef)
118 class(pde_filter_t), intent(inout) :: this
119 type(json_file), intent(inout) :: json
120 type(coef_t), intent(in) :: coef
121
122 ! user parameters
123 call json_get(json, "filter.radius", this%r)
124
125 call json_get_or_default(json, "filter.tolerance", this%abstol_filt, &
126 1.0e-10_rp)
127
128 call json_get_or_default(json, "filter.max_iter", this%ksp_max_iter, 200)
129
130 call json_get_or_default(json, "filter.solver", this%ksp_solver, 'cg')
131
132 call json_get_or_default(json, "filter.preconditioner", &
133 this%precon_type_filt, 'jacobi')
134
135 call this%init_base(json, coef)
136 call pde_filter_init_from_attributes(this, coef)
137
138 end subroutine pde_filter_init_from_json
139
141 subroutine pde_filter_init_from_attributes(this, coef)
142 class(pde_filter_t), intent(inout) :: this
143 type(coef_t), intent(in) :: coef
144 integer :: n
145
146 n = this%coef%dof%size()
147
148 ! init the bc list (all Neuman BCs, will remain empty)
149 call this%bclst_filt%init()
150
151 ! Setup backend dependent Ax routines
152 call ax_helm_factory(this%Ax, full_formulation = .false.)
153
154 ! set up krylov solver
155 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
156 this%ksp_max_iter, this%abstol_filt)
157
158 ! set up preconditioner
159 call filter_precon_factory(this%pc_filt, this%ksp_filt, &
160 this%coef, this%coef%dof, &
161 this%coef%gs_h, &
162 this%bclst_filt, this%precon_type_filt)
163
165
167 subroutine pde_filter_free(this)
168 class(pde_filter_t), intent(inout) :: this
169
170 if (allocated(this%Ax)) then
171 deallocate(this%Ax)
172 end if
173
174 if (allocated(this%ksp_filt)) then
175 call krylov_solver_destroy(this%ksp_filt)
176 deallocate(this%ksp_filt)
177 end if
178
179 if (allocated(this%pc_filt)) then
180 call precon_destroy(this%pc_filt)
181 deallocate(this%pc_filt)
182 end if
183
184 call this%bclst_filt%free()
185
186 call this%free_base()
187
188 end subroutine pde_filter_free
189
193 subroutine pde_filter_apply(this, F_out, F_in)
194 class(pde_filter_t), intent(inout) :: this
195 type(field_t), intent(in) :: F_in
196 type(field_t), intent(inout) :: F_out
197 integer :: n, i
198 ! type(field_t), pointer :: RHS
199 type(field_t) :: RHS, d_F_out
200 character(len=LOG_SIZE) :: log_buf
201 ! integer :: temp_indices(1)
202
203 n = this%coef%dof%size()
204 ! TODO
205 ! This is a bit awkward, because the init for the source terms occurs
206 ! before the init of the scratch registry.
207 ! So we can't use the scratch registry here.
208 ! call neko_scratch_registry%request_field(RHS, temp_indices(1))
209 call rhs%init(this%coef%dof)
210 call d_f_out%init(this%coef%dof)
211
212 ! in a similar fasion to pressure/velocity, we will solve for d_F_out.
213
214 ! to improve convergence, we use F_in as an initial guess for F_out.
215 ! so F_out = F_in + d_F_in.
216
217 ! Defining the operator A = -r^2 \nabla^2 + I
218 ! the system changes from:
219 ! A (F_out) = F_in
220 ! to
221 ! A (d_F_out) = F_in - A(F_in)
222
223 ! set up Helmholtz operators and RHS
224 if (neko_bcknd_device .eq. 1) then
225 ! TODO
226 ! I think this is correct but I've never tested it
227 call device_cfill(this%coef%h1_d, this%r**2, n)
228 call device_cfill(this%coef%h2_d, 1.0_rp, n)
229 else
230 do i = 1, n
231 ! h1 is already negative in its definition
232 this%coef%h1(i,1,1,1) = this%r**2
233 ! ax_helm includes the mass matrix in h2
234 this%coef%h2(i,1,1,1) = 1.0_rp
235 end do
236 end if
237 this%coef%ifh2 = .true.
238
239 ! compute the A(F_in) component of the RHS
240 ! (note, to be safe with the inout intent we first copy F_in to the
241 ! temporary d_F_out)
242 call field_copy(d_f_out, f_in)
243 call this%Ax%compute(rhs%x, d_f_out%x, this%coef, this%coef%msh, &
244 this%coef%Xh)
245
246 if (neko_bcknd_device .eq. 1) then
247 call device_subcol3(rhs%x_d, f_in%x_d, this%coef%B_d, n)
248 call device_cmult(rhs%x_d, -1.0_rp, n)
249 else
250 do i = 1, n
251 ! mass matrix should be included here
252 rhs%x(i,1,1,1) = f_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
253 - rhs%x(i,1,1,1)
254 end do
255 end if
256
257 ! gather scatter
258 call this%coef%gs_h%op(rhs, gs_op_add)
259
260 ! set BCs
261 call this%bclst_filt%apply_scalar(rhs%x, n)
262
263 ! Solve Helmholtz equation
264 call profiler_start_region('filter solve')
265 this%ksp_results(1) = &
266 this%ksp_filt%solve(this%Ax, d_f_out, rhs%x, n, this%coef, &
267 this%bclst_filt, this%coef%gs_h)
268
270
271 ! add result
272 call field_add3(f_out, f_in, d_f_out)
273 ! update preconditioner (needed?)
274 call this%pc_filt%update()
275
276 ! write it all out
277 call neko_log%message('Filter')
278
279 write(log_buf, '(A,A,A)') 'Iterations: ',&
280 'Start residual: ', 'Final residual:'
281 call neko_log%message(log_buf)
282 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
283 this%ksp_results%res_start, this%ksp_results%res_final
284 call neko_log%message(log_buf)
285
286 !call neko_scratch_registry%relinquish_field(temp_indices)
287 call rhs%free()
288 call d_f_out%free()
289
290 end subroutine pde_filter_apply
291
293 subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, &
294 pctype)
295 class(pc_t), allocatable, target, intent(inout) :: pc
296 class(ksp_t), target, intent(inout) :: ksp
297 type(coef_t), target, intent(in) :: coef
298 type(dofmap_t), target, intent(in) :: dof
299 type(gs_t), target, intent(inout) :: gs
300 type(bc_list_t), target, intent(inout) :: bclst
301 character(len=*) :: pctype
302
303 call precon_factory(pc, pctype)
304
305 select type (pcp => pc)
306 type is (jacobi_t)
307 call pcp%init(coef, dof, gs)
308 type is (sx_jacobi_t)
309 call pcp%init(coef, dof, gs)
310 type is (device_jacobi_t)
311 call pcp%init(coef, dof, gs)
312 type is (hsmg_t)
313 if (len_trim(pctype) .gt. 4) then
314 if (index(pctype, '+') .eq. 5) then
315 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst, &
316 trim(pctype(6:)))
317 else
318 call neko_error('Unknown coarse grid solver')
319 end if
320 else
321 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
322 end if
323 end select
324
325 call ksp%set_pc(pc)
326
327 end subroutine filter_precon_factory
328
329end 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)
Multiplication by constant c .
subroutine, public device_subcol3(a_d, b_d, c_d, n)
Returns .
subroutine, public device_cfill(a_d, c, n)
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 registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines a field.
Definition field.f90:34
Filter to be applied to a scalar field.
Definition filter.f90:38
Gather-scatter.
Krylov preconditioner.
Definition pc_hsmg.f90:61
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:65
integer, parameter, public log_size
Definition log.f90:42
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:56
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:58
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_attributes(this, coef)
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:109
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
type(scratch_registry_t) function init(dof, size, expansion_size)
Constructor, optionally taking initial registry and expansion size as argument.
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:47
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Defines a jacobi preconditioner.
Base abstract class for filter.
Definition filter.f90:48
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:68
A Neumann boundary condition for scalar fields. Sets the flux of the field to the chosen value.
Definition neumann.f90:51
A PDE based filter mapping $\rho \mapsto \tilde{\rho}$, see Lazarov & O. Sigmund 2010,...
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.