Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
phmg.f90
Go to the documentation of this file.
1! Copyright (c) 2024-2025, 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!
34module phmg
35 use num_types, only : rp
36 use precon, only : pc_t
37 use gather_scatter, only : gs_t, gs_op_add
38 use space, only : space_t, gll
39 use dofmap, only : dofmap_t
40 use field, only : field_t
41 use coefs, only : coef_t
42 use mesh, only : mesh_t
43 use bc, only : bc_t
44 use bc_list, only : bc_list_t
45 use dirichlet, only : dirichlet_t
46 use utils, only : neko_error
47 use cheby, only : cheby_t
49 use jacobi, only : jacobi_t
51 use schwarz, only : schwarz_t
52 use ax_product, only : ax_t, ax_helm_factory
55 use json_module, only : json_file
57 use math, only : copy, col2, add2, sub3, add2s2
63 use krylov, only : ksp_t, ksp_monitor_t, ksp_max_iter, &
64 krylov_solver_factory
66 use logger, only : neko_log, log_size
67 use, intrinsic :: iso_c_binding
68 implicit none
69 private
70
71
72 type, private :: phmg_lvl_t
73 integer :: lvl = -1
74 integer :: smoother_itrs = 10
75 type(space_t), pointer :: xh
76 type(dofmap_t), pointer :: dm_xh
77 type(gs_t), pointer :: gs_h
79 type(cheby_t) :: cheby
83 type(coef_t), pointer :: coef
84 type(bc_list_t) :: bclst
85 type(dirichlet_t) :: bc
86 type(field_t) :: r, w, z
87 end type phmg_lvl_t
88
89 type, public :: phmg_hrchy_t
90 type(phmg_lvl_t), allocatable :: lvl(:)
91 end type phmg_hrchy_t
92
93
94 type, public, extends(pc_t) :: phmg_t
95 type(tamg_solver_t) :: amg_solver
96 integer :: nlvls
97 type(phmg_hrchy_t) :: phmg_hrchy
98 class(ax_t), allocatable :: ax
99 type(interpolator_t), allocatable :: intrp(:)
100 type(mesh_t), pointer :: msh
101 contains
102 procedure, pass(this) :: init => phmg_init
103 procedure, pass(this) :: init_from_components => &
105 procedure, pass(this) :: free => phmg_free
106 procedure, pass(this) :: solve => phmg_solve
107 procedure, pass(this) :: update => phmg_update
108 procedure, private, pass(this) :: mg_cycle => phmg_mg_cycle
109 end type phmg_t
110
111contains
112
113 subroutine phmg_init(this, coef, bclst, phmg_params)
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(:)
122
123 call json_get_or_default(phmg_params, 'smoother_iterations', &
124 smoother_itrs, 3)
125
126 call json_get_or_default(phmg_params, 'smoother_cheby_acc', &
127 cheby_acc, "jacobi")
128
129 call json_get_or_default(phmg_params, 'coarse_grid.levels', &
130 crs_tamg_lvls, 3)
131
132 call json_get_or_default(phmg_params, 'coarse_grid.iterations', &
133 crs_tamg_itrs, 1)
134
135 call json_get_or_default(phmg_params, 'coarse_grid.cheby_degree', &
136 crs_tamg_cheby_degree, 4)
137
138 if (phmg_params%valid_path('pcoarsening_schedule')) then
139 call json_get(phmg_params, 'pcoarsening_schedule', pcrs_sched)
140 else
141 allocate(pcrs_sched(2))
142 pcrs_sched(1) = 3
143 pcrs_sched(2) = 1
144 end if
145
146
147 call this%init_from_components(coef, bclst, smoother_itrs, &
148 cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree,&
149 pcrs_sched)
150
151 end subroutine phmg_init
152
153 subroutine phmg_init_from_components(this, coef, bclst, smoother_itrs, &
154 cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree, &
155 pcrs_sched)
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
169 use_jacobi = .true.
170 use_cheby = .true.
171
172 this%msh => coef%msh
173
174 this%nlvls = size(pcrs_sched) + 1
175 allocate(lx_lvls(0:this%nlvls - 1))
176 lx_lvls(1:) = pcrs_sched + 1
177
178 allocate(this%phmg_hrchy%lvl(0:this%nlvls - 1))
179
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
186
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)
192
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), &
196 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)
201 end do
202
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)
207
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())
211
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()
215 bc_j => bclst%get(j)
216 call this%phmg_hrchy%lvl(i)%bc%mark_facets(bc_j%marked_facet)
217 end do
218 end if
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)
223
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, &
231 coef%msh)
232 end if
233
234 if (neko_bcknd_device .eq. 1) then
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)
239 else
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)
244 end if
245
246 if (neko_bcknd_device .eq. 1) then
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)
251 st = 1
252 else
253 call this%phmg_hrchy%lvl(i)%cheby_device%init( &
254 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs)
255 st = 0
256 if (trim(cheby_acc) .eq. "schwarz") then
257 this%phmg_hrchy%lvl(i)%cheby_device%schwarz => &
258 this%phmg_hrchy%lvl(i)%schwarz
259 st = 2
260 end if
261 end if
262 else
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)
267 st = 1
268 else
269 call this%phmg_hrchy%lvl(i)%cheby%init( &
270 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs)
271 st = 0
272 if (trim(cheby_acc) .eq. "schwarz") then
273 this%phmg_hrchy%lvl(i)%cheby%schwarz => &
274 this%phmg_hrchy%lvl(i)%schwarz
275 st = 2
276 end if
277 end if
278 end if
279
280 end do
281
282 call print_phmg_info(this%nlvls, st, this%phmg_hrchy)
283
284 ! Create backend specific Ax operator
285 call ax_helm_factory(this%ax, full_formulation = .false.)
286
287 ! Interpolator Fine + mg levels
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)
292 end do
293
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)
299
300 end subroutine phmg_init_from_components
301
302 subroutine phmg_free(this)
303 class(phmg_t), intent(inout) :: this
304 end subroutine phmg_free
305
306 subroutine phmg_solve(this, z, r, n)
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
312 type(ksp_monitor_t) :: ksp_results
313
314 call profiler_start_region('PHMG_solve', 8)
315 associate( mglvl => this%phmg_hrchy%lvl)
316 if (neko_bcknd_device .eq. 1) then
317 z_d = device_get_ptr(z)
318 r_d = device_get_ptr(r)
319 !We should not work with the input
320 call device_copy(mglvl(0)%r%x_d, r_d, n)
321
322 call device_rzero(mglvl(0)%z%x_d, n)
323 call device_rzero(mglvl(0)%w%x_d, n)
324
325 call this%mg_cycle()
326
327 call device_copy(z_d, mglvl(0)%z%x_d, n)
328 else
329 !We should not work with the input
330 call copy(mglvl(0)%r%x, r, n)
331
332 mglvl(0)%z%x = 0.0_rp
333 mglvl(0)%w%x = 0.0_rp
334
335 call this%mg_cycle()
336
337 call copy(z, mglvl(0)%z%x, n)
338 end if
339 end associate
340 call profiler_end_region('PHMG_solve', 8)
341
342 end subroutine phmg_solve
343
344 subroutine phmg_update(this)
345 class(phmg_t), intent(inout) :: this
346 end subroutine phmg_update
347
348
349 subroutine phmg_mg_cycle(this)
350 class(phmg_t), intent(inout) :: this
351 type(ksp_monitor_t) :: ksp_results
352 character(len=2) :: lvl_name
353 integer :: lvl
354
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)
361 !------------!
362 ! SMOOTH !
363 !------------!
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)
370 else
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)
376 end if
377
378 !------------!
379 ! Residual !
380 !------------!
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())
385
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())
388 else
389 w%x = r%x - w%x
390 end if
391
392 !------------!
393 ! Restrict !
394 !------------!
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())
397 else
398 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
399 end if
400
401 call intrp(lvl+1)%map(mg(lvl+1)%r%x, w%x, msh%nelv, mg(lvl+1)%Xh)
402
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)
406
407 call mg(lvl+1)%bclst%apply_scalar( &
408 mg(lvl+1)%r%x, &
409 mg(lvl+1)%dm_Xh%size())
410
411 if (neko_bcknd_device .eq. 1) then
412 call device_rzero(mg(lvl+1)%z%x_d, mg(lvl+1)%dm_Xh%size())
413 else
414 mg(lvl+1)%z%x = 0.0_rp
415 end if
416 end associate
417 call profiler_end_region( "PHMG_level_" // trim(lvl_name))
418 end do
419
420 call profiler_start_region( 'PHMG_coarse-solve' )
421 !------------!
422 ! SOLVE !
423 !------------!
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' )
428
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)
433 !------------!
434 ! Project !
435 !------------!
436 call intrp(lvl+1)%map(w%x, mg(lvl+1)%z%x, msh%nelv, mg(lvl)%Xh)
437
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)
440
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())
443 else
444 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
445 end if
446
447 !------------!
448 ! Correct !
449 !------------!
450 if (neko_bcknd_device .eq. 1) then
451 call device_add2(z%x_d, w%x_d, mg(lvl)%dm_Xh%size())
452 else
453 z%x = z%x + w%x
454 end if
455
456 !------------!
457 ! SMOOTH !
458 !------------!
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)
464 else
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)
469 end if
470 end associate
471 call profiler_end_region( "PHMG_level_" // trim(lvl_name))
472 end do
473 end associate
474
475 end subroutine phmg_mg_cycle
476
486 subroutine phmg_jacobi_smoother(z, r, w, mg, msh, Ax, n, lvl)
487 type(phmg_lvl_t) :: mg
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
493
494 ni = mg%smoother_itrs
495 if (neko_bcknd_device .eq. 1) then
496 do i = 1, ni
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)
502
503 call mg%device_jacobi%solve(w%x, w%x, n)
504
505 call device_add2s2(z%x_d, w%x_d, 0.6_rp, n)
506 end do
507 else
508 do i = 1, ni
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)
513
514 call mg%jacobi%solve(w%x, w%x, n)
515
516 call add2s2(z%x, w%x, 0.6_rp, n)
517 end do
518 end if
519 end subroutine phmg_jacobi_smoother
520
521
522 subroutine phmg_resid_monitor(z, r, w, mg, msh, Ax, lvl, typ)
523 integer :: lvl, typ
524 type(phmg_lvl_t) :: mg
525 class(ax_t), intent(inout) :: Ax
526 type(mesh_t), intent(inout) :: msh
527 type(field_t) :: z, r, w
528 real(kind=rp) :: val
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())
535 if (typ .eq. 1) then
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
547 else
548 write(log_buf, '(A15,I4,F12.6)') 'RESID', lvl, val
549 end if
550 call neko_log%message(log_buf)
551 end subroutine phmg_resid_monitor
552
553 subroutine print_phmg_info(nlvls, smoo_type, phmg)
554 integer, intent(in) :: nlvls
555 integer, intent(in) :: smoo_type
556 type(phmg_hrchy_t) :: phmg
557 integer :: i, clvl
558 character(len=LOG_SIZE) :: log_buf, smoo_name
559
560 call neko_log%section('PHMG')
561
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'
566 else
567 write(smoo_name, '(A5)') 'CHEBY'
568 end if
569
570 write(log_buf, '(A28,I2,A8)') &
571 'Creating PHMG hierarchy with', &
572 nlvls, 'levels.'
573 call neko_log%message(log_buf)
574
575 clvl = nlvls - 1
576 do i = 0, nlvls-1
577 write(log_buf, '(A8,I2,A8,I2)') &
578 '-- level', i, '-- lx:', phmg%lvl(i)%Xh%lx
579 call neko_log%message(log_buf)
580
581 if (i .eq. clvl) then
582 write(log_buf, '(A19,A20)') &
583 'Solve:', 'tAMG'
584 call neko_log%message(log_buf)
585 else
586 write(log_buf, '(A22,A20)') &
587 'Smoother:', &
588 trim(smoo_name)
589 call neko_log%message(log_buf)
590
591 write(log_buf, '(A28,I2)') &
592 'Smoother Iters:', &
593 phmg%lvl(i)%smoother_itrs
594 call neko_log%message(log_buf)
595 end if
596 end do
597
598 call neko_log%end_section()
599
600 end subroutine print_phmg_info
601
602end module phmg
__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.
Definition device.F90:107
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
Defines a boundary condition.
Definition bc.f90:34
Chebyshev preconditioner.
Chebyshev preconditioner.
Definition cheby.f90:34
Coefficients.
Definition coef.f90:34
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.
Definition device.F90:34
subroutine, public device_stream_wait_event(stream, event, flags)
Synchronize a device stream with an event.
Definition device.F90:1362
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:51
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:62
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Gather-scatter.
Routines to interpolate between different spaces.
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:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:781
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:812
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Hybrid ph-multigrid preconditioner.
Definition phmg.f90:34
subroutine phmg_jacobi_smoother(z, r, w, mg, msh, ax, n, lvl)
Wraps jacobi solve as a residual update relaxation method.
Definition phmg.f90:487
subroutine phmg_init_from_components(this, coef, bclst, smoother_itrs, cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree, pcrs_sched)
Definition phmg.f90:156
subroutine phmg_solve(this, z, r, n)
Definition phmg.f90:307
subroutine phmg_init(this, coef, bclst, phmg_params)
Definition phmg.f90:114
subroutine phmg_update(this)
Definition phmg.f90:345
subroutine phmg_free(this)
Definition phmg.f90:303
subroutine phmg_mg_cycle(this)
Definition phmg.f90:350
subroutine phmg_resid_monitor(z, r, w, mg, msh, ax, lvl, typ)
Definition phmg.f90:523
subroutine print_phmg_info(nlvls, smoo_type, phmg)
Definition phmg.f90:554
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:107
Overlapping schwarz solves.
Definition schwarz.f90:61
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
Implements multigrid using the TreeAMG hierarchy structure. USE:
Utilities.
Definition utils.f90:35
Base type for a matrix-vector product providing .
Definition ax.f90:43
Base type for a boundary condition.
Definition bc.f90:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Defines a Chebyshev preconditioner.
Definition cheby.f90:52
Defines a Chebyshev preconditioner.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:61
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
Interpolation between two space::space_t.
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:73
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
The function space for the SEM solution fields.
Definition space.f90:63
Type for the TreeAMG solver.