Neko 1.99.1
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
64 use krylov, only : ksp_t, ksp_monitor_t, ksp_max_iter, &
65 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 end type phmg_t
109
110contains
111
112 subroutine phmg_init(this, coef, bclst, phmg_params)
113 class(phmg_t), intent(inout), target :: this
114 type(coef_t), intent(in), target :: coef
115 type(bc_list_t), intent(inout), target :: bclst
116 type(json_file), intent(inout) :: phmg_params
117 integer :: crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree
118 integer :: smoother_itrs
119 character(len=:), allocatable :: cheby_acc
120 integer, allocatable :: pcrs_sched(:)
121
122 call json_get_or_default(phmg_params, 'smoother_iterations', &
123 smoother_itrs, 10)
124
125 call json_get_or_default(phmg_params, 'smoother_cheby_acc', &
126 cheby_acc, "jacobi")
127
128 call json_get_or_default(phmg_params, 'coarse_grid.levels', &
129 crs_tamg_lvls, 3)
130
131 call json_get_or_default(phmg_params, 'coarse_grid.iterations', &
132 crs_tamg_itrs, 1)
133
134 call json_get_or_default(phmg_params, 'coarse_grid.cheby_degree', &
135 crs_tamg_cheby_degree, 5)
136
137 if (phmg_params%valid_path('pcoarsening_schedule')) then
138 call json_get(phmg_params, 'pcoarsening_schedule', pcrs_sched)
139 else
140 allocate(pcrs_sched(2))
141 pcrs_sched(1) = 3
142 pcrs_sched(2) = 1
143 end if
144
145
146 call this%init_from_components(coef, bclst, smoother_itrs, &
147 cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree,&
148 pcrs_sched)
149
150 end subroutine phmg_init
151
152 subroutine phmg_init_from_components(this, coef, bclst, smoother_itrs, &
153 cheby_acc, crs_tamg_lvls, crs_tamg_itrs, crs_tamg_cheby_degree, &
154 pcrs_sched)
155 class(phmg_t), intent(inout), target :: this
156 type(coef_t), intent(in), target :: coef
157 type(bc_list_t), intent(inout), target :: bclst
158 integer, intent(in) :: smoother_itrs
159 character(len=:), allocatable :: cheby_acc
160 integer, intent(in) :: crs_tamg_lvls, crs_tamg_itrs
161 integer, intent(in) :: crs_tamg_cheby_degree
162 integer, intent(in), allocatable :: pcrs_sched(:)
163 integer :: lx_crs, lx_mid
164 integer, allocatable :: lx_lvls(:)
165 integer :: n, i, j, st
166 class(bc_t), pointer :: bc_j
167 logical :: use_jacobi, use_cheby
168 use_jacobi = .true.
169 use_cheby = .true.
170
171 this%msh => coef%msh
172
173 this%nlvls = size(pcrs_sched) + 1
174 allocate(lx_lvls(0:this%nlvls - 1))
175 lx_lvls(1:) = pcrs_sched + 1
176
177 allocate(this%phmg_hrchy%lvl(0:this%nlvls - 1))
178
179 this%phmg_hrchy%lvl(0)%lvl = 0
180 this%phmg_hrchy%lvl(0)%smoother_itrs = smoother_itrs
181 this%phmg_hrchy%lvl(0)%Xh => coef%Xh
182 this%phmg_hrchy%lvl(0)%coef => coef
183 this%phmg_hrchy%lvl(0)%dm_Xh => coef%dof
184 this%phmg_hrchy%lvl(0)%gs_h => coef%gs_h
185
186 do i = 1, this%nlvls - 1
187 allocate(this%phmg_hrchy%lvl(i)%Xh)
188 allocate(this%phmg_hrchy%lvl(i)%dm_Xh)
189 allocate(this%phmg_hrchy%lvl(i)%gs_h)
190 allocate(this%phmg_hrchy%lvl(i)%coef)
191
192 this%phmg_hrchy%lvl(i)%lvl = i
193 this%phmg_hrchy%lvl(i)%smoother_itrs = smoother_itrs
194 call this%phmg_hrchy%lvl(i)%Xh%init(gll, lx_lvls(i), lx_lvls(i), &
195 lx_lvls(i))
196 call this%phmg_hrchy%lvl(i)%dm_Xh%init(coef%msh, &
197 this%phmg_hrchy%lvl(i)%Xh)
198 call this%phmg_hrchy%lvl(i)%gs_h%init(this%phmg_hrchy%lvl(i)%dm_Xh)
199 call this%phmg_hrchy%lvl(i)%coef%init(this%phmg_hrchy%lvl(i)%gs_h)
200 end do
201
202 do i = 0, this%nlvls - 1
203 call this%phmg_hrchy%lvl(i)%r%init(this%phmg_hrchy%lvl(i)%dm_Xh)
204 call this%phmg_hrchy%lvl(i)%w%init(this%phmg_hrchy%lvl(i)%dm_Xh)
205 call this%phmg_hrchy%lvl(i)%z%init(this%phmg_hrchy%lvl(i)%dm_Xh)
206
207 this%phmg_hrchy%lvl(i)%coef%ifh2 = coef%ifh2
208 call copy(this%phmg_hrchy%lvl(i)%coef%h1, coef%h1, &
209 this%phmg_hrchy%lvl(i)%dm_Xh%size())
210
211 call this%phmg_hrchy%lvl(i)%bc%init_base(this%phmg_hrchy%lvl(i)%coef)
212 if (bclst%size() .gt. 0 ) then
213 do j = 1, bclst%size()
214 bc_j => bclst%get(j)
215 call this%phmg_hrchy%lvl(i)%bc%mark_facets(bc_j%marked_facet)
216 end do
217 end if
218 call this%phmg_hrchy%lvl(i)%bc%finalize()
219 call this%phmg_hrchy%lvl(i)%bc%set_g(0.0_rp)
220 call this%phmg_hrchy%lvl(i)%bclst%init()
221 call this%phmg_hrchy%lvl(i)%bclst%append(this%phmg_hrchy%lvl(i)%bc)
222
224 if (trim(cheby_acc) .eq. "schwarz") then
225 call this%phmg_hrchy%lvl(i)%schwarz%init( &
226 this%phmg_hrchy%lvl(i)%Xh, &
227 this%phmg_hrchy%lvl(i)%dm_Xh, &
228 this%phmg_hrchy%lvl(i)%gs_h, &
229 this%phmg_hrchy%lvl(i)%bclst, &
230 coef%msh)
231 end if
232
233 if (neko_bcknd_device .eq. 1) then
234 call this%phmg_hrchy%lvl(i)%device_jacobi%init(&
235 this%phmg_hrchy%lvl(i)%coef, &
236 this%phmg_hrchy%lvl(i)%dm_Xh, &
237 this%phmg_hrchy%lvl(i)%gs_h)
238 else
239 call this%phmg_hrchy%lvl(i)%jacobi%init(&
240 this%phmg_hrchy%lvl(i)%coef, &
241 this%phmg_hrchy%lvl(i)%dm_Xh, &
242 this%phmg_hrchy%lvl(i)%gs_h)
243 end if
244
245 if (neko_bcknd_device .eq. 1) then
246 if (trim(cheby_acc) .eq. "jacobi") then
247 call this%phmg_hrchy%lvl(i)%cheby_device%init( &
248 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs, &
249 this%phmg_hrchy%lvl(i)%device_jacobi)
250 st = 1
251 else
252 call this%phmg_hrchy%lvl(i)%cheby_device%init( &
253 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs)
254 st = 0
255 if (trim(cheby_acc) .eq. "schwarz") then
256 this%phmg_hrchy%lvl(i)%cheby_device%schwarz => &
257 this%phmg_hrchy%lvl(i)%schwarz
258 st = 2
259 end if
260 end if
261 else
262 if (trim(cheby_acc) .eq. "jacobi") then
263 call this%phmg_hrchy%lvl(i)%cheby%init( &
264 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs, &
265 this%phmg_hrchy%lvl(i)%jacobi)
266 st = 1
267 else
268 call this%phmg_hrchy%lvl(i)%cheby%init( &
269 this%phmg_hrchy%lvl(i)%dm_Xh%size(), smoother_itrs)
270 st = 0
271 if (trim(cheby_acc) .eq. "schwarz") then
272 this%phmg_hrchy%lvl(i)%cheby%schwarz => &
273 this%phmg_hrchy%lvl(i)%schwarz
274 st = 2
275 end if
276 end if
277 end if
278
279 end do
280
281 call print_phmg_info(this%nlvls, st, this%phmg_hrchy)
282
283 ! Create backend specific Ax operator
284 call ax_helm_factory(this%ax, full_formulation = .false.)
285
286 ! Interpolator Fine + mg levels
287 allocate(this%intrp(this%nlvls - 1))
288 do i = 1, this%nlvls -1
289 call this%intrp(i)%init(this%phmg_hrchy%lvl(i-1)%Xh, &
290 this%phmg_hrchy%lvl(i)%Xh)
291 end do
292
293 call this%amg_solver%init(this%ax, this%phmg_hrchy%lvl(this%nlvls -1)%Xh, &
294 this%phmg_hrchy%lvl(this%nlvls -1)%coef, this%msh, &
295 this%phmg_hrchy%lvl(this%nlvls-1)%gs_h, crs_tamg_lvls, &
296 this%phmg_hrchy%lvl(this%nlvls -1)%bclst, &
297 crs_tamg_itrs, crs_tamg_cheby_degree)
298
299 end subroutine phmg_init_from_components
300
301 subroutine phmg_free(this)
302 class(phmg_t), intent(inout) :: this
303 end subroutine phmg_free
304
305 subroutine phmg_solve(this, z, r, n)
306 class(phmg_t), intent(inout) :: this
307 integer, intent(in) :: n
308 real(kind=rp), dimension(n), intent(inout) :: z
309 real(kind=rp), dimension(n), intent(inout) :: r
310 type(c_ptr) :: z_d, r_d
311 type(ksp_monitor_t) :: ksp_results
312
313
314 associate( mglvl => this%phmg_hrchy%lvl)
315 if (neko_bcknd_device .eq. 1) then
316 z_d = device_get_ptr(z)
317 r_d = device_get_ptr(r)
318 !We should not work with the input
319 call device_copy(mglvl(0)%r%x_d, r_d, n)
320 call device_rzero(mglvl(0)%z%x_d, n)
321 call device_rzero(mglvl(0)%w%x_d, n)
322 call phmg_mg_cycle(mglvl(0)%z, mglvl(0)%r, mglvl(0)%w, 0, &
323 this%nlvls -1, mglvl, this%intrp, this%msh, this%Ax, &
324 this%amg_solver)
325
326 call device_copy(z_d, mglvl(0)%z%x_d, n)
327 else
328 !We should not work with the input
329 call copy(mglvl(0)%r%x, r, n)
330
331 mglvl(0)%z%x = 0.0_rp
332 mglvl(0)%w%x = 0.0_rp
333
334 call phmg_mg_cycle(mglvl(0)%z, mglvl(0)%r, mglvl(0)%w, 0, &
335 this%nlvls -1, mglvl, this%intrp, this%msh, this%Ax, &
336 this%amg_solver)
337
338 call copy(z, mglvl(0)%z%x, n)
339 end if
340 end associate
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 recursive subroutine phmg_mg_cycle(z, r, w, lvl, clvl, &
350 mg, intrp, msh, Ax, amg_solver)
351 type(ksp_monitor_t) :: ksp_results
352 integer :: lvl, clvl
353 type(phmg_lvl_t) :: mg(0:clvl)
354 type(interpolator_t) :: intrp(1:clvl)
355 type(tamg_solver_t), intent(inout) :: amg_solver
356 class(ax_t), intent(inout) :: ax
357 type(mesh_t), intent(inout) :: msh
358 type(field_t) :: z, r, w
359 integer :: i
360 logical :: use_jacobi
361 real(kind=rp) :: val
362
363 use_jacobi = .false.
364 call profiler_start_region('PHMG_cycle', 8)
368 call profiler_start_region('PHMG_PreSmooth', 9)
369 if (use_jacobi) then
370 call phmg_jacobi_smoother(z, r, w, mg(lvl), msh, ax, &
371 mg(lvl)%dm_Xh%size(), lvl)
372 else
373 if (neko_bcknd_device .eq. 1) then
374 mg(lvl)%cheby_device%zero_initial_guess = .true.
375 ksp_results = mg(lvl)%cheby_device%solve(ax, z, &
376 r%x, mg(lvl)%dm_Xh%size(), &
377 mg(lvl)%coef, mg(lvl)%bclst, &
378 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
379 else
380 mg(lvl)%cheby%zero_initial_guess = .true.
381 ksp_results = mg(lvl)%cheby%solve(ax, z, &
382 r%x, mg(lvl)%dm_Xh%size(), &
383 mg(lvl)%coef, mg(lvl)%bclst, &
384 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
385 end if
386 end if
387 call profiler_end_region('PHMG_PreSmooth', 9)
388
392 call ax%compute(w%x, z%x, mg(lvl)%coef, msh, mg(lvl)%Xh)
393 call mg(lvl)%gs_h%op(w%x, mg(lvl)%dm_Xh%size(), gs_op_add, glb_cmd_event)
394 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
395 call mg(lvl)%bclst%apply_scalar(w%x, mg(lvl)%dm_Xh%size())
396
397 if (neko_bcknd_device .eq. 1) then
398 call device_sub3(w%x_d, r%x_d, w%x_d, mg(lvl)%dm_Xh%size())
399 else
400 w%x = r%x - w%x
401 end if
402
406 if (neko_bcknd_device .eq. 1) then
407 call device_col2(w%x_d, mg(lvl)%coef%mult_d, mg(lvl)%dm_Xh%size())
408 else
409 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
410 end if
411
412 call profiler_start_region('PHMG_map_to_coarse', 9)
413 call intrp(lvl+1)%map(mg(lvl+1)%r%x, w%x, msh%nelv, mg(lvl+1)%Xh)
414 call profiler_end_region('PHMG_map_to_coarse', 9)
415
416 call mg(lvl+1)%gs_h%op(mg(lvl+1)%r%x, mg(lvl+1)%dm_Xh%size(), &
417 gs_op_add, glb_cmd_event)
418 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
419
420 call mg(lvl+1)%bclst%apply_scalar( &
421 mg(lvl+1)%r%x, &
422 mg(lvl+1)%dm_Xh%size())
426 if (neko_bcknd_device .eq. 1) then
427 call device_rzero(mg(lvl+1)%z%x_d, mg(lvl+1)%dm_Xh%size())
428 else
429 mg(lvl+1)%z%x = 0.0_rp
430 end if
431 if (lvl+1 .eq. clvl) then
432 call profiler_start_region('PHMG_tAMG_coarse_grid', 9)
433 call amg_solver%solve(mg(lvl+1)%z%x, &
434 mg(lvl+1)%r%x, &
435 mg(lvl+1)%dm_Xh%size())
436 call profiler_end_region('PHMG_tAMG_coarse_grid', 9)
437
438 else
439 call phmg_mg_cycle(mg(lvl+1)%z, mg(lvl+1)%r, mg(lvl+1)%w, lvl+1, &
440 clvl, mg, intrp, msh, ax, amg_solver)
441 end if
442
446 call profiler_start_region('PHMG_map_to_fine', 9)
447 call intrp(lvl+1)%map(w%x, mg(lvl+1)%z%x, msh%nelv, mg(lvl)%Xh)
448 call profiler_end_region('PHMG_map_to_fine', 9)
449
450 call mg(lvl)%gs_h%op(w%x, mg(lvl)%dm_Xh%size(), gs_op_add, glb_cmd_event)
451 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
452
453 if (neko_bcknd_device .eq. 1) then
454 call device_col2(w%x_d, mg(lvl)%coef%mult_d, mg(lvl)%dm_Xh%size())
455 else
456 call col2(w%x, mg(lvl)%coef%mult, mg(lvl)%dm_Xh%size())
457 end if
458
462 if (neko_bcknd_device .eq. 1) then
463 call device_add2(z%x_d, w%x_d, mg(lvl)%dm_Xh%size())
464 else
465 z%x = z%x + w%x
466 end if
467
471 call profiler_start_region('PHMG_PostSmooth', 9)
472 if (use_jacobi) then
473 call phmg_jacobi_smoother(z, r, w, mg(lvl), msh, ax, &
474 mg(lvl)%dm_Xh%size(), lvl)
475 else
476 if (neko_bcknd_device .eq. 1) then
477 ksp_results = mg(lvl)%cheby_device%solve(ax, z, &
478 r%x, mg(lvl)%dm_Xh%size(), &
479 mg(lvl)%coef, mg(lvl)%bclst, &
480 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
481 else
482 ksp_results = mg(lvl)%cheby%solve(ax, z, &
483 r%x, mg(lvl)%dm_Xh%size(), &
484 mg(lvl)%coef, mg(lvl)%bclst, &
485 mg(lvl)%gs_h, niter = mg(lvl)%smoother_itrs)
486 end if
487 end if
488 call profiler_end_region('PHMG_PostSmooth', 9)
489
490 call profiler_end_region('PHMG_cycle', 8)
491 end subroutine phmg_mg_cycle
492
502 subroutine phmg_jacobi_smoother(z, r, w, mg, msh, Ax, n, lvl)
503 type(phmg_lvl_t) :: mg
504 class(ax_t), intent(inout) :: Ax
505 type(mesh_t), intent(inout) :: msh
506 type(field_t), intent(inout) :: z, r, w
507 integer, intent(in) :: n, lvl
508 integer :: i, iblk, ni, niblk
509
510 ni = mg%smoother_itrs
511 if (neko_bcknd_device .eq. 1) then
512 do i = 1, ni
513 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
514 call mg%gs_h%op(w%x, n, gs_op_add, glb_cmd_event)
515 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
516 call mg%bclst%apply_scalar(w%x, n)
517 call device_sub3(w%x_d, r%x_d, w%x_d, n)
518
519 call mg%device_jacobi%solve(w%x, w%x, n)
520
521 call device_add2s2(z%x_d, w%x_d, 0.6_rp, n)
522 end do
523 else
524 do i = 1, ni
525 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
526 call mg%gs_h%op(w%x, n, gs_op_add)
527 call mg%bclst%apply_scalar(w%x, n)
528 call sub3(w%x, r%x, w%x, n)
529
530 call mg%jacobi%solve(w%x, w%x, n)
531
532 call add2s2(z%x, w%x, 0.6_rp, n)
533 end do
534 end if
535 end subroutine phmg_jacobi_smoother
536
537
538 subroutine phmg_resid_monitor(z, r, w, mg, msh, Ax, lvl, typ)
539 integer :: lvl, typ
540 type(phmg_lvl_t) :: mg
541 class(ax_t), intent(inout) :: Ax
542 type(mesh_t), intent(inout) :: msh
543 type(field_t) :: z, r, w
544 real(kind=rp) :: val
545 character(len=LOG_SIZE) :: log_buf
546 call ax%compute(w%x, z%x, mg%coef, msh, mg%Xh)
547 call mg%gs_h%op(w%x, mg%dm_Xh%size(), gs_op_add)
548 call mg%bclst%apply_scalar(w%x, mg%dm_Xh%size())
549 call device_sub3(w%x_d, r%x_d, w%x_d, mg%dm_Xh%size())
550 val = device_glsc2(w%x_d, w%x_d, mg%dm_Xh%size())
551 if (typ .eq. 1) then
552 write(log_buf, '(A15,I4,F12.6)') 'PRESMOO - PRE', lvl, val
553 else if (typ .eq. 2) then
554 write(log_buf, '(A15,I4,F12.6)') 'PRESMOO -POST', lvl, val
555 else if (typ .eq. 3) then
556 write(log_buf, '(A15,I4,F12.6)') 'POSTSMOO- PRE', lvl, val
557 else if (typ .eq. 4) then
558 write(log_buf, '(A15,I4,F12.6)') 'POSTSMOO-POST', lvl, val
559 else if (typ .eq. 5) then
560 write(log_buf, '(A15,I4,F12.6)') 'TAMG - PRE', lvl, val
561 else if (typ .eq. 6) then
562 write(log_buf, '(A15,I4,F12.6)') 'TAMG -POST', lvl, val
563 else
564 write(log_buf, '(A15,I4,F12.6)') 'RESID', lvl, val
565 end if
566 call neko_log%message(log_buf)
567 end subroutine phmg_resid_monitor
568
569 subroutine print_phmg_info(nlvls, smoo_type, phmg)
570 integer, intent(in) :: nlvls
571 integer, intent(in) :: smoo_type
572 type(phmg_hrchy_t) :: phmg
573 integer :: i, clvl
574 character(len=LOG_SIZE) :: log_buf, smoo_name
575
576 call neko_log%section('PHMG')
577
578 if (smoo_type .eq. 1) then
579 write(smoo_name, '(A16)') 'CHEBY-acc JACOBI'
580 else if (smoo_type .eq. 2) then
581 write(smoo_name, '(A17)') 'CHEBY-acc SCHWARZ'
582 else
583 write(smoo_name, '(A5)') 'CHEBY'
584 end if
585
586 write(log_buf, '(A28,I2,A8)') &
587 'Creating PHMG hierarchy with', &
588 nlvls, 'levels.'
589 call neko_log%message(log_buf)
590
591 clvl = nlvls - 1
592 do i = 0, nlvls-1
593 write(log_buf, '(A8,I2,A8,I2)') &
594 '-- level', i, '-- lx:', phmg%lvl(i)%Xh%lx
595 call neko_log%message(log_buf)
596
597 if (i .eq. clvl) then
598 write(log_buf, '(A19,A20)') &
599 'Solve:', 'tAMG'
600 call neko_log%message(log_buf)
601 else
602 write(log_buf, '(A22,A20)') &
603 'Smoother:', &
604 trim(smoo_name)
605 call neko_log%message(log_buf)
606
607 write(log_buf, '(A28,I2)') &
608 'Smoother Iters:', &
609 phmg%lvl(i)%smoother_itrs
610 call neko_log%message(log_buf)
611 end if
612 end do
613
614 call neko_log%end_section()
615
616 end subroutine print_phmg_info
617
618end 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:101
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:1208
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:70
integer, parameter, public log_size
Definition log.f90:46
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:503
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:155
subroutine phmg_solve(this, z, r, n)
Definition phmg.f90:306
subroutine phmg_init(this, coef, bclst, phmg_params)
Definition phmg.f90:113
subroutine phmg_update(this)
Definition phmg.f90:345
subroutine phmg_free(this)
Definition phmg.f90:302
subroutine phmg_resid_monitor(z, r, w, mg, msh, ax, lvl, typ)
Definition phmg.f90:539
subroutine print_phmg_info(nlvls, smoo_type, phmg)
Definition phmg.f90:570
recursive subroutine phmg_mg_cycle(z, r, w, lvl, clvl, mg, intrp, msh, ax, amg_solver)
Definition phmg.f90:351
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:61
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:55
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:48
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.