38  use json_module, 
only: json_file
 
   45  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), 
intent(in) :: coef
 
  121    call json_get(json, 
"filter.radius", this%r)
 
  131         this%precon_type_filt, 
'jacobi')
 
  133    call this%init_base(json, coef)
 
 
  141    type(
coef_t), 
intent(in) :: coef
 
  144    n = this%coef%dof%size()
 
  147    call this%bclst_filt%init()
 
  150    call ax_helm_factory(this%Ax, full_formulation = .false.)
 
  153    call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
 
  154         this%ksp_max_iter, this%abstol_filt)
 
  158         this%coef, this%coef%dof, &
 
  160         this%bclst_filt, this%precon_type_filt)
 
 
  168    if (
allocated(this%Ax)) 
then 
  172    if (
allocated(this%ksp_filt)) 
then 
  173       call this%ksp_filt%free()
 
  174       deallocate(this%ksp_filt)
 
  177    if (
allocated(this%pc_filt)) 
then 
  178       call precon_destroy(this%pc_filt)
 
  179       deallocate(this%pc_filt)
 
  182    call this%bclst_filt%free()
 
  184    call this%free_base()
 
 
  193    type(
field_t), 
intent(in) :: F_in
 
  194    type(
field_t), 
intent(inout) :: F_out
 
  198    character(len=LOG_SIZE) :: log_buf
 
  201    n = this%coef%dof%size()
 
  207    call rhs%init(this%coef%dof)
 
  208    call d_f_out%init(this%coef%dof)
 
  230          this%coef%h1(i,1,1,1) = this%r**2
 
  232          this%coef%h2(i,1,1,1) = 1.0_rp
 
  235    this%coef%ifh2 = .true.
 
  241    call this%Ax%compute(rhs%x, d_f_out%x, this%coef, this%coef%msh, &
 
  250          rhs%x(i,1,1,1) = f_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
 
  256    call this%coef%gs_h%op(rhs, gs_op_add)
 
  259    call this%bclst_filt%apply_scalar(rhs%x, n)
 
  263    this%ksp_results(1) = &
 
  264         this%ksp_filt%solve(this%Ax, d_f_out, rhs%x, n, this%coef, &
 
  265         this%bclst_filt, this%coef%gs_h)
 
  272    call this%pc_filt%update()
 
  277    write(log_buf, 
'(A,A,A)') 
'Iterations:   ',&
 
  278         'Start residual:     ', 
'Final residual:' 
  280    write(log_buf, 
'(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
 
  281         this%ksp_results%res_start, this%ksp_results%res_final
 
 
  293    class(
pc_t), 
allocatable, 
target, 
intent(inout) :: pc
 
  294    class(
ksp_t), 
target, 
intent(inout) :: ksp
 
  295    type(
coef_t), 
target, 
intent(in) :: coef
 
  296    type(
dofmap_t), 
target, 
intent(in) :: dof
 
  297    type(
gs_t), 
target, 
intent(inout) :: gs
 
  298    type(
bc_list_t), 
target, 
intent(inout) :: bclst
 
  299    character(len=*) :: pctype
 
  301    call precon_factory(pc, pctype)
 
  303    select type (pcp => pc)
 
  305       call pcp%init(coef, dof, gs)
 
  307       call pcp%init(coef, dof, gs)
 
  309       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 .
 
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_components(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), 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.