Neko 1.99.1
A portable framework for high-order spectral element flow simulations
|
Implements multigrid using the TreeAMG hierarchy structure. USE: More...
Data Types | |
type | tamg_solver_t |
Type for the TreeAMG solver. More... | |
type | tamg_wrk_t |
Functions/Subroutines | |
subroutine | tamg_mg_init (this, ax, xh, coef, msh, gs_h, nlvls, blst, max_iter, cheby_degree) |
Initialization of the TreeAMG multigrid solver. | |
subroutine | tamg_mg_solve (this, z, r, n) |
Solver function for the TreeAMG solver object. | |
recursive subroutine | tamg_mg_cycle (x, b, n, lvl, amg, mgstuff, zero_initial_guess) |
Recrsive multigrid cycle for the TreeAMG solver object. | |
recursive subroutine | tamg_mg_cycle_d (x, b, x_d, b_d, n, lvl, amg, mgstuff, zero_initial_guess) |
Recrsive multigrid cycle for the TreeAMG solver object on device. | |
subroutine | calc_resid (r, x, b, amg, lvl, n) |
Wrapper function to calculate residyal. | |
subroutine | print_preagg_info (lvl, nagg) |
subroutine | print_resid_info (r, x, b, r_d, x_d, b_d, amg, lvl, n) |
subroutine | fill_lvl_map (amg) |
type(tamg_hierarchy_t) :: amg type(tamg_solver_t) :: amg_solver
call amginit(ax, Xh, coef, msh, gs_h, 3) call amg_solverinit(amg, niter)
call amg_solversolve(xx, f, n)
|
private |
r | The residual to be returned |
x | The current solution |
b | The right-hand side |
amg | The TreeAMG object |
lvl | Current level of the cycle |
n | Number of dofs |
Definition at line 362 of file tree_amg_multigrid.f90.
|
private |
|
private |
Definition at line 386 of file tree_amg_multigrid.f90.
|
private |
x | The solution to be returned |
b | The right-hand side |
n | Number of dofs |
lvl | Current level of the cycle |
amg | The TreeAMG object |
mgstuff | The Solver object. TODO: rename this |
-------—<! SMOOTH <! -------—<!
Is coarsest grid.
-------—<! Residual <! -------—<!
-------—<! Restrict <! -------—<!
----------------—<! Call Coarse solve <! ----------------—<!
-------—<! Project <! -------—<!
-------—<! Correct <! -------—<!
-------—<! SMOOTH <! -------—<!
Definition at line 239 of file tree_amg_multigrid.f90.
|
private |
x | The solution to be returned |
b | The right-hand side |
n | Number of dofs |
lvl | Current level of the cycle |
amg | The TreeAMG object |
mgstuff | The Solver object. TODO: rename this |
-------—<! SMOOTH <! -------—<!
Is coarsest grid.
-------—<! Residual <! -------—<!
-------—<! Restrict <! -------—<!
----------------—<! Call Coarse solve <! ----------------—<!
-------—<! Project <! -------—<!
-------—<! Correct <! -------—<!
-------—<! SMOOTH <! -------—<!
Definition at line 298 of file tree_amg_multigrid.f90.
subroutine tree_amg_multigrid::tamg_mg_init | ( | class(tamg_solver_t), intent(inout), target | this, |
class(ax_t), intent(in), target | ax, | ||
type(space_t), intent(in), target | xh, | ||
type(coef_t), intent(in), target | coef, | ||
type(mesh_t), intent(in), target | msh, | ||
type(gs_t), intent(in), target | gs_h, | ||
integer, intent(in) | nlvls, | ||
type(bc_list_t), intent(in), target | blst, | ||
integer, intent(in) | max_iter, | ||
integer, intent(in) | cheby_degree | ||
) |
ax | Finest level matvec operator |
Xh | Finest level field |
coef | Finest level coeff thing |
msh | Finest level mesh information |
gs_h | Finest level gather scatter operator |
nlvls | Number of levels for the TreeAMG hierarchy |
blst | Finest level BC list |
max_iter | Number of AMG iterations |
Create level 1 (neko elements are level 0)
Create the remaining levels
Create the end point
Allocate work space on each level
Definition at line 102 of file tree_amg_multigrid.f90.
|
private |
z | The solution to be returned |
r | The right-hand side |
n | Number of dofs |
Definition at line 194 of file tree_amg_multigrid.f90.