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