38 use json_module,
only: json_file
46 use precon,
only:
pc_t, precon_factory, precon_destroy
77 class(
ax_t),
allocatable :: ax
81 class(
ksp_t),
allocatable :: ksp_filt
83 class(
pc_t),
allocatable :: pc_filt
91 real(kind=
rp) :: abstol_filt
95 character(len=:),
allocatable :: ksp_solver
97 character(len=:),
allocatable :: precon_type_filt
98 integer :: ksp_n, n, i
106 procedure, pass(this) :: init_from_attributes => &
119 type(json_file),
intent(inout) :: json
120 type(
coef_t),
intent(in) :: coef
123 call json_get(json,
"filter.radius", this%r)
133 this%precon_type_filt,
'jacobi')
135 call this%init_base(json, coef)
143 type(
coef_t),
intent(in) :: coef
146 n = this%coef%dof%size()
149 call this%bclst_filt%init()
152 call ax_helm_factory(this%Ax, full_formulation = .false.)
155 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
156 this%ksp_max_iter, this%abstol_filt)
160 this%coef, this%coef%dof, &
162 this%bclst_filt, this%precon_type_filt)
170 if (
allocated(this%Ax))
then
174 if (
allocated(this%ksp_filt))
then
175 call krylov_solver_destroy(this%ksp_filt)
176 deallocate(this%ksp_filt)
179 if (
allocated(this%pc_filt))
then
180 call precon_destroy(this%pc_filt)
181 deallocate(this%pc_filt)
184 call this%bclst_filt%free()
186 call this%free_base()
195 type(
field_t),
intent(in) :: F_in
196 type(
field_t),
intent(inout) :: F_out
200 character(len=LOG_SIZE) :: log_buf
203 n = this%coef%dof%size()
209 call rhs%init(this%coef%dof)
210 call d_f_out%init(this%coef%dof)
232 this%coef%h1(i,1,1,1) = this%r**2
234 this%coef%h2(i,1,1,1) = 1.0_rp
237 this%coef%ifh2 = .true.
243 call this%Ax%compute(rhs%x, d_f_out%x, this%coef, this%coef%msh, &
252 rhs%x(i,1,1,1) = f_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
258 call this%coef%gs_h%op(rhs, gs_op_add)
261 call this%bclst_filt%apply_scalar(rhs%x, n)
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)
274 call this%pc_filt%update()
279 write(log_buf,
'(A,A,A)')
'Iterations: ',&
280 'Start residual: ',
'Final residual:'
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
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
303 call precon_factory(pc, pctype)
305 select type (pcp => pc)
307 call pcp%init(coef, dof, gs)
309 call pcp%init(coef, dof, gs)
311 call pcp%init(coef, dof, gs)
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, &
321 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
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.
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.
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.
Filter to be applied to a scalar field.
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
integer, parameter, public ksp_max_iter
Maximum number of iters.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
integer, parameter neko_bcknd_device
Defines a Neumann boundary condition.
integer, parameter, public rp
Global precision used in computations.
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.
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
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.
Base type for a matrix-vector product providing .
A list of allocatable `bc_t`. Follows the standard interface of lists.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Defines a jacobi preconditioner.
Base abstract class for filter.
Defines a jacobi preconditioner.
Type for storing initial and final residuals in a Krylov solver.
Base abstract type for a canonical Krylov method, solving .
A Neumann boundary condition for scalar fields. Sets the flux of the field to the chosen value.
A PDE based filter mapping $\rho \mapsto \tilde{\rho}$, see Lazarov & O. Sigmund 2010,...
Abstract type to compute pressure residual.
Defines a canonical Krylov preconditioner.
Defines a jacobi preconditioner for SX-Aurora.