38 use json_module,
only : json_file
44 use precon,
only :
pc_t, precon_factory, precon_destroy
75 class(
ax_t),
allocatable :: ax
79 class(
ksp_t),
allocatable :: ksp_filt
81 class(
pc_t),
allocatable :: pc_filt
89 real(kind=
rp) :: abstol_filt
93 character(len=:),
allocatable :: ksp_solver
95 character(len=:),
allocatable :: precon_type_filt
96 integer :: ksp_n, n, i
104 procedure, pass(this) :: init_from_components => &
117 type(json_file),
intent(inout) :: json
118 type(
coef_t),
target,
intent(in) :: coef
119 real(kind=
rp) :: r, tol
121 character(len=:),
allocatable :: ksp_solver, precon_type
130 call this%init_from_components(coef, r, tol, max_iter, ksp_solver, &
137 ksp_solver, precon_type)
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
145 call this%init_base(coef)
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
154 n = this%coef%dof%size()
157 call this%bclst_filt%init()
160 call ax_helm_factory(this%Ax, full_formulation = .false.)
163 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
164 this%ksp_max_iter, this%abstol_filt)
168 this%coef, this%coef%dof, this%coef%gs_h, this%bclst_filt, &
169 this%precon_type_filt)
177 if (
allocated(this%Ax))
then
181 if (
allocated(this%ksp_filt))
then
182 call this%ksp_filt%free()
183 deallocate(this%ksp_filt)
186 if (
allocated(this%pc_filt))
then
187 call precon_destroy(this%pc_filt)
188 deallocate(this%pc_filt)
191 if (
allocated(this%ksp_solver))
then
192 deallocate(this%ksp_solver)
195 if (
allocated(this%precon_type_filt))
then
196 deallocate(this%precon_type_filt)
199 call this%bclst_filt%free()
201 call this%free_base()
211 type(
field_t),
intent(in) :: F_in
212 type(
field_t),
intent(inout) :: F_out
214 type(
field_t),
pointer :: RHS, d_F_out
215 character(len=LOG_SIZE) :: log_buf
216 integer :: temp_indices(2)
218 n = this%coef%dof%size()
234 call device_cfill(this%coef%h1_d, (this%r / (2.0_rp * sqrt(3.0_rp)))**2, n)
238 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
240 this%coef%h2 = 1.0_rp
242 this%coef%ifh2 = .true.
248 call this%Ax%compute(rhs%x, d_f_out%x, this%coef, this%coef%msh, &
257 rhs%x(i,1,1,1) = f_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
263 call this%coef%gs_h%op(rhs, gs_op_add)
266 call this%bclst_filt%apply_scalar(rhs%x, n)
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)
279 call this%pc_filt%update()
284 write(log_buf,
'(A,A,A)')
'Iterations: ',
'Start residual: ', &
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
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
306 call precon_factory(pc, pctype)
308 select type (pcp => pc)
310 call pcp%init(coef, dof, gs)
312 call pcp%init(coef, dof, gs)
314 call pcp%init(coef, dof, gs)
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, 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.
subroutine, public field_copy(a, b, n)
Copy a vector .
subroutine, public field_add3(a, b, c, n)
Vector addition .
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_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.
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 solution fields.
type(registry_t), target, public neko_registry
Global field registry.
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.
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. Sets the flux of the field to the chosen values.
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .
Abstract type to compute pressure residual.
Defines a canonical Krylov preconditioner.
Defines a jacobi preconditioner for SX-Aurora.