DL_MG
|
Store of convergence and iteration control parameters: More...
Data Types | |
type | conv_params_t |
Functions/Subroutines | |
subroutine | dl_mg_set_conv_params (eps_half, eq_type, fd_order, tol_res_rel, tol_res_abs, tol_pot_rel, tol_pot_abs, tol_newton_res_rel, tol_newton_res_abs, tol_cg_res_rel, tol_cg_res_abs, tol_vcyc_res_rel, tol_vcyc_res_abs, tol_level_1_res_rel, tol_level_1_res_abs, tol_mu_rel, tol_mu_abs, mg_max_conv_rate, max_iters_defco, max_iters_newton, max_iters_cg, max_iters_vcyc, max_iters_level_1, v_iterations, omega_sor, use_cg) |
subroutine | write_convergence_params (eps_half, eq_type, fd_order, tol_res_rel, tol_res_abs, tol_pot_rel, tol_pot_abs, tol_newton_res_rel, tol_newton_res_abs, tol_cg_res_rel, tol_cg_res_abs, tol_vcyc_res_rel, tol_vcyc_res_abs, tol_level_1_res_rel, tol_level_1_res_abs, tol_mu_rel, tol_mu_abs, mg_max_conv_rate, max_iters_defco, max_iters_newton, max_iters_cg, max_iters_vcyc, max_iters_level_1, v_iterations, omega_sor, use_cg) |
Variables | |
real(wp) | tol_res_rel_ = 1.0e-2_wp |
real(wp) | tol_res_abs_ = 1.0e-2_wp |
real(wp) | tol_pot_rel_ = 1.0e-6_wp |
real(wp) | tol_pot_abs_ = 1.0e-6_wp |
real(wp) | tol_newton_res_rel_ = 1.0e-8_wp |
real(wp) | tol_newton_res_abs_ = 1.0e-8_wp |
real(wp) | tol_cg_res_rel_ = 1.0e-8_wp |
real(wp) | tol_cg_res_abs_ = 1.0e-8_wp |
real(wp) | tol_vcyc_res_rel_ = 1.0e-8_wp |
real(wp) | tol_vcyc_res_abs_ = 1.0e-8_wp |
real(wp) | tol_level_1_res_rel_ |
initialisation is done in the subroutine from below More... | |
real(wp) | tol_level_1_res_abs_ = 1.0e-12_wp |
real(wp) | tol_mu_rel_ = 1.0e-3_wp |
real(wp) | tol_mu_abs_ = 1.0e-3_wp |
integer | max_iters_defco_ = 30 |
maximum number of iterations for the higher order correction More... | |
integer | max_iters_newton_ = 30 |
.... Newton iteration More... | |
integer | max_iters_cg_ = 100 |
integer | max_iters_vcyc_ = 30 |
maximum number of V cycles (recommended value not more that 50) if conjugate gradient is on this is reset to 1, unless set explicitly via input argument or envinroment variable More... | |
integer | max_iters_level_1_ = 100 |
max number of iteration at the coarsest level More... | |
integer, dimension(2) | v_iterations_ = (/2, 1/) |
Number of smoother iteration to be applied at each level. First value is for the descending branch (coarsening), the second is for the ascending branch (prolongation). More... | |
real(wp) | omega_sor_ = 1.15 |
SOR relaxation parameter, the default is set tp 1.15 as recommended in Trottenberg's book. More... | |
real(wp) | mg_max_conv_rate_ = 0.85_wp |
if consecutive residuals ratio is larger than the value below it is assumed that round-off limit was reached. i.e. the computation cannot progress becase of roundoff errors, or perhaps the problem is too hard for this algorithm More... | |
logical, save | use_cg_ = .false. |
integer, parameter | tconv_level_1 =1 |
convergence tags used in subroutine mg_convergence_test: MG coarse level, MG V-cycle, Cojugate Gradient, Newton More... | |
integer, parameter | tconv_vcycle =2 |
integer, parameter | tconv_cg =3 |
integer, parameter | tconv_newton =4 |
integer, parameter | nconv_tags =4 |
type(conv_params_t), dimension(nconv_tags), save | conv_params |
used in mg_convergence_test More... | |
Store of convergence and iteration control parameters:
The governing convergence criteria is as follow:
A. if defect correction order is > 2 the following applies:
defect norm < MAX(tol_rel_res * (initial defect norm), tol_res_abs)
AND
||pot^(n) -pot^(n-1)|| < MAX(tol_res_pot * ||pot^(n-1)||, tol_abs_pot)
These two conditions try to ensure that the iteration process does not stops too early because of a norm reduction stall or small defect. The tol_abs_{res,pot} values ensure stopping in the case of zero solution or account for roundoff errors.
B. defect correction == 2
In this case the standard the multigrid convergence algorithm is used
residual norm < MAX(tol_rel_res * || f||, tol_res_abs)
where tol_res_rel, and tol_res_abs are transferrend to the selected equation (Poisson or or PB flavor) and their tolerance parameters are ignored if not present in the argument list. Also, tol_rel_pot and tol_abs_pot are ignored in this cases.
When the PBE equation is solved with periodic boundary conditions the chemical potential convergence is tested with the following condition:
abs(mu^(i) - mu^(i-1)) < MAX(tol_mu_rel * abs(mu^(i), tol_mu_abs)
Convergece is reached if the above condition and the Newton converge test are satisfied
subroutine dl_mg_convergence_params::dl_mg_set_conv_params | ( | real(wp), dimension(:,:,:,:), intent(in) | eps_half, |
integer, intent(in) | eq_type, | ||
integer, intent(in) | fd_order, | ||
real(wp), intent(in), optional | tol_res_rel, | ||
real(wp), intent(in), optional | tol_res_abs, | ||
real(wp), intent(in), optional | tol_pot_rel, | ||
real(wp), intent(in), optional | tol_pot_abs, | ||
real(wp), intent(in), optional | tol_newton_res_rel, | ||
real(wp), intent(in), optional | tol_newton_res_abs, | ||
real(wp), intent(in), optional | tol_cg_res_rel, | ||
real(wp), intent(in), optional | tol_cg_res_abs, | ||
real(wp), intent(in), optional | tol_vcyc_res_rel, | ||
real(wp), intent(in), optional | tol_vcyc_res_abs, | ||
real(wp), intent(in), optional | tol_level_1_res_rel, | ||
real(wp), intent(in), optional | tol_level_1_res_abs, | ||
real(wp), intent(in), optional | tol_mu_rel, | ||
real(wp), intent(in), optional | tol_mu_abs, | ||
real(wp), intent(in), optional | mg_max_conv_rate, | ||
integer, intent(in), optional | max_iters_defco, | ||
integer, intent(in), optional | max_iters_newton, | ||
integer, intent(in), optional | max_iters_cg, | ||
integer, intent(in), optional | max_iters_vcyc, | ||
integer, intent(in), optional | max_iters_level_1, | ||
integer, dimension(2), intent(in), optional | v_iterations, | ||
real(wp), intent(in), optional | omega_sor, | ||
logical, intent(in), optional | use_cg | ||
) |
[in] | eps_half | to check for internal |
set more convergence parameters if provided
subroutine dl_mg_convergence_params::write_convergence_params | ( | real(wp), dimension(:,:,:,:), intent(in) | eps_half, |
integer, intent(in) | eq_type, | ||
integer, intent(in) | fd_order, | ||
real(wp), intent(in) | tol_res_rel, | ||
real(wp), intent(in) | tol_res_abs, | ||
real(wp), intent(in) | tol_pot_rel, | ||
real(wp), intent(in) | tol_pot_abs, | ||
real(wp), intent(in) | tol_newton_res_rel, | ||
real(wp), intent(in) | tol_newton_res_abs, | ||
real(wp), intent(in) | tol_cg_res_rel, | ||
real(wp), intent(in) | tol_cg_res_abs, | ||
real(wp), intent(in) | tol_vcyc_res_rel, | ||
real(wp), intent(in) | tol_vcyc_res_abs, | ||
real(wp), intent(in) | tol_level_1_res_rel, | ||
real(wp), intent(in) | tol_level_1_res_abs, | ||
real(wp), intent(in) | tol_mu_rel, | ||
real(wp), intent(in) | tol_mu_abs, | ||
real(wp), intent(in) | mg_max_conv_rate, | ||
integer, intent(in) | max_iters_defco, | ||
integer, intent(in) | max_iters_newton, | ||
integer, intent(in) | max_iters_cg, | ||
integer, intent(in) | max_iters_vcyc, | ||
integer, intent(in) | max_iters_level_1, | ||
integer, dimension(2), intent(in) | v_iterations, | ||
real(wp), intent(in) | omega_sor, | ||
logical, intent(in) | use_cg | ||
) |
type(conv_params_t), dimension(nconv_tags), save dl_mg_convergence_params::conv_params |
used in mg_convergence_test
integer dl_mg_convergence_params::max_iters_cg_ = 100 |
integer dl_mg_convergence_params::max_iters_defco_ = 30 |
maximum number of iterations for the higher order correction
integer dl_mg_convergence_params::max_iters_level_1_ = 100 |
max number of iteration at the coarsest level
integer dl_mg_convergence_params::max_iters_newton_ = 30 |
.... Newton iteration
integer dl_mg_convergence_params::max_iters_vcyc_ = 30 |
maximum number of V cycles (recommended value not more that 50) if conjugate gradient is on this is reset to 1, unless set explicitly via input argument or envinroment variable
real(wp) dl_mg_convergence_params::mg_max_conv_rate_ = 0.85_wp |
if consecutive residuals ratio is larger than the value below it is assumed that round-off limit was reached. i.e. the computation cannot progress becase of roundoff errors, or perhaps the problem is too hard for this algorithm
integer, parameter dl_mg_convergence_params::nconv_tags =4 |
real(wp) dl_mg_convergence_params::omega_sor_ = 1.15 |
SOR relaxation parameter, the default is set tp 1.15 as recommended in Trottenberg's book.
integer, parameter dl_mg_convergence_params::tconv_cg =3 |
integer, parameter dl_mg_convergence_params::tconv_level_1 =1 |
convergence tags used in subroutine mg_convergence_test: MG coarse level, MG V-cycle, Cojugate Gradient, Newton
integer, parameter dl_mg_convergence_params::tconv_newton =4 |
integer, parameter dl_mg_convergence_params::tconv_vcycle =2 |
real(wp) dl_mg_convergence_params::tol_cg_res_abs_ = 1.0e-8_wp |
real(wp) dl_mg_convergence_params::tol_cg_res_rel_ = 1.0e-8_wp |
real(wp) dl_mg_convergence_params::tol_level_1_res_abs_ = 1.0e-12_wp |
real(wp) dl_mg_convergence_params::tol_level_1_res_rel_ |
initialisation is done in the subroutine from below
real(wp) dl_mg_convergence_params::tol_mu_abs_ = 1.0e-3_wp |
real(wp) dl_mg_convergence_params::tol_mu_rel_ = 1.0e-3_wp |
real(wp) dl_mg_convergence_params::tol_newton_res_abs_ = 1.0e-8_wp |
real(wp) dl_mg_convergence_params::tol_newton_res_rel_ = 1.0e-8_wp |
real(wp) dl_mg_convergence_params::tol_pot_abs_ = 1.0e-6_wp |
real(wp) dl_mg_convergence_params::tol_pot_rel_ = 1.0e-6_wp |
real(wp) dl_mg_convergence_params::tol_res_abs_ = 1.0e-2_wp |
real(wp) dl_mg_convergence_params::tol_res_rel_ = 1.0e-2_wp |
real(wp) dl_mg_convergence_params::tol_vcyc_res_abs_ = 1.0e-8_wp |
real(wp) dl_mg_convergence_params::tol_vcyc_res_rel_ = 1.0e-8_wp |
logical, save dl_mg_convergence_params::use_cg_ = .false. |
integer, dimension(2) dl_mg_convergence_params::v_iterations_ = (/2, 1/) |
Number of smoother iteration to be applied at each level. First value is for the descending branch (coarsening), the second is for the ascending branch (prolongation).