55 use json_module,
only : json_file
67 use,
intrinsic :: iso_c_binding
74 integer :: smoother_itrs = 10
77 type(
gs_t),
pointer :: gs_h
98 class(
ax_t),
allocatable :: ax
103 procedure, pass(this) :: init_from_components => &
114 class(
phmg_t),
intent(inout),
target :: this
115 type(
coef_t),
intent(in),
target :: coef
116 type(
bc_list_t),
intent(inout),
target :: bclst
117 type(json_file),
intent(inout) :: phmg_params
118 integer :: crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree
119 integer :: smoother_itrs
120 character(len=:),
allocatable :: cheby_acc
121 integer,
allocatable :: pcrs_sched(:)
136 crs_tamg_cheby_degree, 4)
138 if (phmg_params%valid_path(
'pcoarsening_schedule'))
then
139 call json_get(phmg_params,
'pcoarsening_schedule', pcrs_sched)
141 allocate(pcrs_sched(2))
147 call this%init_from_components(coef, bclst, smoother_itrs, &
148 cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree,&
154 cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree, &
156 class(
phmg_t),
intent(inout),
target :: this
157 type(
coef_t),
intent(in),
target :: coef
158 type(
bc_list_t),
intent(inout),
target :: bclst
159 integer,
intent(in) :: smoother_itrs
160 character(len=:),
allocatable :: cheby_acc
161 integer,
intent(in) :: crs_tamg_lvls, crs_tamg_itrs
162 integer,
intent(in) :: crs_tamg_cheby_degree
163 integer,
intent(in),
allocatable :: pcrs_sched(:)
164 integer :: lx_crs, lx_mid
165 integer,
allocatable :: lx_lvls(:)
166 integer :: n, i, j, st
167 class(
bc_t),
pointer :: bc_j
168 logical :: use_jacobi, use_cheby
174 this%nlvls =
size(pcrs_sched) + 1
175 allocate(lx_lvls(0:this%nlvls - 1))
176 lx_lvls(1:) = pcrs_sched + 1
178 allocate(this%phmg_hrchy%lvl(0:this%nlvls - 1))
180 this%phmg_hrchy%lvl(0)%lvl = 0
181 this%phmg_hrchy%lvl(0)%smoother_itrs = smoother_itrs
182 this%phmg_hrchy%lvl(0)%Xh => coef%Xh
183 this%phmg_hrchy%lvl(0)%coef => coef
184 this%phmg_hrchy%lvl(0)%dm_Xh => coef%dof
185 this%phmg_hrchy%lvl(0)%gs_h => coef%gs_h
187 do i = 1, this%nlvls - 1
188 allocate(this%phmg_hrchy%lvl(i)%Xh)
189 allocate(this%phmg_hrchy%lvl(i)%dm_Xh)
190 allocate(this%phmg_hrchy%lvl(i)%gs_h)
191 allocate(this%phmg_hrchy%lvl(i)%coef)
193 this%phmg_hrchy%lvl(i)%lvl = i
194 this%phmg_hrchy%lvl(i)%smoother_itrs = smoother_itrs
195 call this%phmg_hrchy%lvl(i)%Xh%init(
gll, lx_lvls(i), lx_lvls(i), &
197 call this%phmg_hrchy%lvl(i)%dm_Xh%init(coef%msh, &
198 this%phmg_hrchy%lvl(i)%Xh)
199 call this%phmg_hrchy%lvl(i)%gs_h%init(this%phmg_hrchy%lvl(i)%dm_Xh)
200 call this%phmg_hrchy%lvl(i)%coef%init(this%phmg_hrchy%lvl(i)%gs_h)
203 do i = 0, this%nlvls - 1
204 call this%phmg_hrchy%lvl(i)%r%init(this%phmg_hrchy%lvl(i)%dm_Xh)
205 call this%phmg_hrchy%lvl(i)%w%init(this%phmg_hrchy%lvl(i)%dm_Xh)
206 call this%phmg_hrchy%lvl(i)%z%init(this%phmg_hrchy%lvl(i)%dm_Xh)
208 this%phmg_hrchy%lvl(i)%coef%ifh2 = coef%ifh2
209 call copy(this%phmg_hrchy%lvl(i)%coef%h1, coef%h1, &
210 this%phmg_hrchy%lvl(i)%dm_Xh%size())
212 call this%phmg_hrchy%lvl(i)%bc%init_base(this%phmg_hrchy%lvl(i)%coef)
213 if (bclst%size() .gt. 0 )
then
214 do j = 1, bclst%size()
216 call this%phmg_hrchy%lvl(i)%bc%mark_facets(bc_j%marked_facet)
219 call this%phmg_hrchy%lvl(i)%bc%finalize()
220 call this%phmg_hrchy%lvl(i)%bc%set_g(0.0_rp)
221 call this%phmg_hrchy%lvl(i)%bclst%init()
222 call this%phmg_hrchy%lvl(i)%bclst%append(this%phmg_hrchy%lvl(i)%bc)
225 if (trim(cheby_acc) .eq.
"schwarz")
then
226 call this%phmg_hrchy%lvl(i)%schwarz%init( &
227 this%phmg_hrchy%lvl(i)%Xh, &
228 this%phmg_hrchy%lvl(i)%dm_Xh, &
229 this%phmg_hrchy%lvl(i)%gs_h, &
230 this%phmg_hrchy%lvl(i)%bclst, &
235 call this%phmg_hrchy%lvl(i)%device_jacobi%init(&
236 this%phmg_hrchy%lvl(i)%coef, &
237 this%phmg_hrchy%lvl(i)%dm_Xh, &
238 this%phmg_hrchy%lvl(i)%gs_h)
240 call this%phmg_hrchy%lvl(i)%jacobi%init(&
241 this%phmg_hrchy%lvl(i)%coef, &
242 this%phmg_hrchy%lvl(i)%dm_Xh, &
243 this%phmg_hrchy%lvl(i)%gs_h)
247 if (trim(cheby_acc) .eq.
"jacobi")
then
248 call this%phmg_hrchy%lvl(i)%cheby_device%init( &
249 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs, &
250 this%phmg_hrchy%lvl(i)%device_jacobi)
253 call this%phmg_hrchy%lvl(i)%cheby_device%init( &
254 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs)
256 if (trim(cheby_acc) .eq.
"schwarz")
then
257 this%phmg_hrchy%lvl(i)%cheby_device%schwarz => &
258 this%phmg_hrchy%lvl(i)%schwarz
263 if (trim(cheby_acc) .eq.
"jacobi")
then
264 call this%phmg_hrchy%lvl(i)%cheby%init( &
265 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs, &
266 this%phmg_hrchy%lvl(i)%jacobi)
269 call this%phmg_hrchy%lvl(i)%cheby%init( &
270 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs)
272 if (trim(cheby_acc) .eq.
"schwarz")
then
273 this%phmg_hrchy%lvl(i)%cheby%schwarz => &
274 this%phmg_hrchy%lvl(i)%schwarz
285 call ax_helm_factory(this%ax, full_formulation = .false.)
288 allocate(this%intrp(this%nlvls - 1))
289 do i = 1, this%nlvls -1
290 call this%intrp(i)%init(this%phmg_hrchy%lvl(i-1)%Xh, &
291 this%phmg_hrchy%lvl(i)%Xh)
294 call this%amg_solver%init(this%ax, this%phmg_hrchy%lvl(this%nlvls -1)%Xh, &
295 this%phmg_hrchy%lvl(this%nlvls -1)%coef, this%msh, &
296 this%phmg_hrchy%lvl(this%nlvls-1)%gs_h, crs_tamg_lvls, &
297 this%phmg_hrchy%lvl(this%nlvls -1)%bclst, &
298 crs_tamg_itrs, crs_tamg_cheby_degree)
303 class(
phmg_t),
intent(inout) :: this
307 class(
phmg_t),
intent(inout) :: this
308 integer,
intent(in) :: n
309 real(kind=
rp),
dimension(n),
intent(inout) :: z
310 real(kind=
rp),
dimension(n),
intent(inout) :: r
311 type(c_ptr) :: z_d, r_d
315 associate( mglvl => this%phmg_hrchy%lvl)
330 call copy(mglvl(0)%r%x, r, n)
332 mglvl(0)%z%x = 0.0_rp
333 mglvl(0)%w%x = 0.0_rp
337 call copy(z, mglvl(0)%z%x, n)
345 class(
phmg_t),
intent(inout) :: this
350 class(
phmg_t),
intent(inout) :: this
351 type(ksp_monitor_t) :: ksp_results
352 character(len=2) :: lvl_name
355 associate(mg => this%phmg_hrchy%lvl, intrp => this%intrp, &
356 msh => this%msh, ax => this%Ax)
357 do lvl = 0, this%nlvls-2
358 write(lvl_name,
'(I0)') lvl
359 call profiler_start_region(
"PHMG_level_" // trim(lvl_name))
360 associate(z => mg(lvl)%z, r => mg(lvl)%r, w => mg(lvl)%w)
364 if (neko_bcknd_device .eq. 1)
then
365 mg(lvl)%cheby_device%zero_initial_guess = .true.
366 ksp_results = mg(lvl)%cheby_device%solve(ax, z, &
367 r%x, mg(lvl)%dm_Xh%size(), &
368 mg(lvl)%coef, mg(lvl)%bclst, &
369 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
371 mg(lvl)%cheby%zero_initial_guess = .true.
372 ksp_results = mg(lvl)%cheby%solve(ax, z, &
373 r%x, mg(lvl)%dm_Xh%size(), &
374 mg(lvl)%coef, mg(lvl)%bclst, &
375 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
381 call ax%compute(w%x, z%x, mg(lvl)%coef, msh, mg(lvl)%Xh)
382 call mg(lvl)%gs_h%op(w%x, mg(lvl)%dm_Xh%size(), gs_op_add, glb_cmd_event)
383 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
384 call mg(lvl)%bclst%apply_scalar(w%x, mg(lvl)%dm_Xh%size())
386 if (neko_bcknd_device .eq. 1)
then
387 call device_sub3(w%x_d, r%x_d, w%x_d, mg(lvl)%dm_Xh%size())
395 if (neko_bcknd_device .eq. 1)
then
396 call device_col2(w%x_d, mg(lvl)%coef%mult_d, mg(lvl)%dm_Xh%size())
398 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
401 call intrp(lvl+1)%map(mg(lvl+1)%r%x, w%x, msh%nelv, mg(lvl+1)%Xh)
403 call mg(lvl+1)%gs_h%op(mg(lvl+1)%r%x, mg(lvl+1)%dm_Xh%size(), &
404 gs_op_add, glb_cmd_event)
405 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
407 call mg(lvl+1)%bclst%apply_scalar( &
409 mg(lvl+1)%dm_Xh%size())
411 if (neko_bcknd_device .eq. 1)
then
412 call device_rzero(mg(lvl+1)%z%x_d, mg(lvl+1)%dm_Xh%size())
414 mg(lvl+1)%z%x = 0.0_rp
417 call profiler_end_region(
"PHMG_level_" // trim(lvl_name))
420 call profiler_start_region(
'PHMG_coarse-solve' )
424 call this%amg_solver%solve(mg(this%nlvls-1)%z%x, &
425 mg(this%nlvls-1)%r%x, &
426 mg(this%nlvls-1)%dm_Xh%size())
427 call profiler_end_region(
'PHMG_coarse-solve' )
429 do lvl = (this%nlvls-2), 0, -1
430 write(lvl_name,
'(I0)') lvl
431 call profiler_start_region(
"PHMG_level_" // trim(lvl_name))
432 associate(z => mg(lvl)%z, r => mg(lvl)%r, w => mg(lvl)%w)
436 call intrp(lvl+1)%map(w%x, mg(lvl+1)%z%x, msh%nelv, mg(lvl)%Xh)
438 call mg(lvl)%gs_h%op(w%x, mg(lvl)%dm_Xh%size(), gs_op_add, glb_cmd_event)
439 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
441 if (neko_bcknd_device .eq. 1)
then
442 call device_col2(w%x_d, mg(lvl)%coef%mult_d, mg(lvl)%dm_Xh%size())
444 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
450 if (neko_bcknd_device .eq. 1)
then
451 call device_add2(z%x_d, w%x_d, mg(lvl)%dm_Xh%size())
459 if (neko_bcknd_device .eq. 1)
then
460 ksp_results = mg(lvl)%cheby_device%solve(ax, z, &
461 r%x, mg(lvl)%dm_Xh%size(), &
462 mg(lvl)%coef, mg(lvl)%bclst, &
463 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
465 ksp_results = mg(lvl)%cheby%solve(ax, z, &
466 r%x, mg(lvl)%dm_Xh%size(), &
467 mg(lvl)%coef, mg(lvl)%bclst, &
468 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
471 call profiler_end_region(
"PHMG_level_" // trim(lvl_name))
488 class(ax_t),
intent(inout) :: Ax
489 type(mesh_t),
intent(inout) :: msh
490 type(field_t),
intent(inout) :: z, r, w
491 integer,
intent(in) :: n, lvl
492 integer :: i, iblk, ni, niblk
494 ni = mg%smoother_itrs
495 if (neko_bcknd_device .eq. 1)
then
497 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
498 call mg%gs_h%op(w%x, n, gs_op_add, glb_cmd_event)
499 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
500 call mg%bclst%apply_scalar(w%x, n)
501 call device_sub3(w%x_d, r%x_d, w%x_d, n)
503 call mg%device_jacobi%solve(w%x, w%x, n)
505 call device_add2s2(z%x_d, w%x_d, 0.6_rp, n)
509 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
510 call mg%gs_h%op(w%x, n, gs_op_add)
511 call mg%bclst%apply_scalar(w%x, n)
512 call sub3(w%x, r%x, w%x, n)
514 call mg%jacobi%solve(w%x, w%x, n)
516 call add2s2(z%x, w%x, 0.6_rp, n)
525 class(ax_t),
intent(inout) :: Ax
526 type(mesh_t),
intent(inout) :: msh
527 type(field_t) :: z, r, w
529 character(len=LOG_SIZE) :: log_buf
530 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
531 call mg%gs_h%op(w%x, mg%dm_Xh%size(), gs_op_add)
532 call mg%bclst%apply_scalar(w%x, mg%dm_Xh%size())
533 call device_sub3(w%x_d, r%x_d, w%x_d, mg%dm_Xh%size())
534 val = device_glsc2(w%x_d, w%x_d, mg%dm_Xh%size())
536 write(log_buf,
'(A15,I4,F12.6)')
'PRESMOO - PRE', lvl, val
537 else if (typ .eq. 2)
then
538 write(log_buf,
'(A15,I4,F12.6)')
'PRESMOO -POST', lvl, val
539 else if (typ .eq. 3)
then
540 write(log_buf,
'(A15,I4,F12.6)')
'POSTSMOO- PRE', lvl, val
541 else if (typ .eq. 4)
then
542 write(log_buf,
'(A15,I4,F12.6)')
'POSTSMOO-POST', lvl, val
543 else if (typ .eq. 5)
then
544 write(log_buf,
'(A15,I4,F12.6)')
'TAMG - PRE', lvl, val
545 else if (typ .eq. 6)
then
546 write(log_buf,
'(A15,I4,F12.6)')
'TAMG -POST', lvl, val
548 write(log_buf,
'(A15,I4,F12.6)')
'RESID', lvl, val
550 call neko_log%message(log_buf)
554 integer,
intent(in) :: nlvls
555 integer,
intent(in) :: smoo_type
558 character(len=LOG_SIZE) :: log_buf, smoo_name
560 call neko_log%section(
'PHMG')
562 if (smoo_type .eq. 1)
then
563 write(smoo_name,
'(A16)')
'CHEBY-acc JACOBI'
564 else if (smoo_type .eq. 2)
then
565 write(smoo_name,
'(A17)')
'CHEBY-acc SCHWARZ'
567 write(smoo_name,
'(A5)')
'CHEBY'
570 write(log_buf,
'(A28,I2,A8)') &
571 'Creating PHMG hierarchy with', &
573 call neko_log%message(log_buf)
577 write(log_buf,
'(A8,I2,A8,I2)') &
578 '-- level', i,
'-- lx:',
phmg%lvl(i)%Xh%lx
579 call neko_log%message(log_buf)
581 if (i .eq. clvl)
then
582 write(log_buf,
'(A19,A20)') &
584 call neko_log%message(log_buf)
586 write(log_buf,
'(A22,A20)') &
589 call neko_log%message(log_buf)
591 write(log_buf,
'(A28,I2)') &
593 phmg%lvl(i)%smoother_itrs
594 call neko_log%message(log_buf)
598 call neko_log%end_section()
__device__ T solve(const T u, const T y, const T guess, const T nu, const T kappa, const T B)
Return the device pointer for an associated Fortran array.
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.
Defines a boundary condition.
Chebyshev preconditioner.
Chebyshev preconditioner.
Jacobi preconditioner accelerator backend.
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_invcol2(a_d, b_d, n, strm)
Vector division .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
Device abstraction, common interface for various accelerators.
subroutine, public device_stream_wait_event(stream, event, flags)
Synchronize a device stream with an event.
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Defines a dirichlet boundary condition.
Defines a mapping of the degrees of freedom.
Routines to interpolate between different spaces.
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
subroutine, public sub3(a, b, c, n)
Vector subtraction .
subroutine, public add2(a, b, n)
Vector addition .
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Hybrid ph-multigrid preconditioner.
subroutine phmg_jacobi_smoother(z, r, w, mg, msh, ax, n, lvl)
Wraps jacobi solve as a residual update relaxation method.
subroutine phmg_init_from_components(this, coef, bclst, smoother_itrs, cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree, pcrs_sched)
subroutine phmg_solve(this, z, r, n)
subroutine phmg_init(this, coef, bclst, phmg_params)
subroutine phmg_update(this)
subroutine phmg_free(this)
subroutine phmg_mg_cycle(this)
subroutine phmg_resid_monitor(z, r, w, mg, msh, ax, lvl, typ)
subroutine print_phmg_info(nlvls, smoo_type, phmg)
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.
Overlapping schwarz solves.
Defines a function space.
integer, parameter, public gll
Implements multigrid using the TreeAMG hierarchy structure. USE:
Base type for a matrix-vector product providing .
Base type for a boundary condition.
A list of allocatable `bc_t`. Follows the standard interface of lists.
Defines a Chebyshev preconditioner.
Defines a Chebyshev preconditioner.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Interpolation between two space::space_t.
Defines a jacobi preconditioner.
Type for storing initial and final residuals in a Krylov solver.
Base abstract type for a canonical Krylov method, solving .
Defines a canonical Krylov preconditioner.
The function space for the SEM solution fields.
Type for the TreeAMG solver.