DL_MG
|
auxiliary subroutines More...
Functions/Subroutines | |
subroutine | dl_mg_init_nonlin (eq_type, lambda, temp, n, c, q, rho, steric_weight, der_pot, expcap, linearised, use_fas, neutralisation_method, neutralisation_ion_ratios, ierror) |
Initialises the non-linear parameters. It can be called multiple times, e.g. if beta needs to be changed in continuation calculation. More... | |
subroutine | compute_rho_sum (rho) |
subroutine | compute_steric_weight_sum (has_stw, steric_weight) |
computes the accesible volume More... | |
subroutine | test_charge_neutrality (rho, pot, steric_weight) |
called by the PBE solver More... | |
subroutine | remove_zero_mode_source (mg, eq_type, remove_mode, rho, steric_weight, der_pot) |
removes zero modes from source term. Needed when using full PBC More... | |
subroutine | remove_zero_mode_potential (mg, eq_type) |
removes zero mode from the potential. Needed when using full PBC at the moment this is called only inside vcycle for mgp it migh need to be extend for arbitrary arrays with optional args N.B. the potential shift is applied also to the halo points to save a halo exchange. This is so because this sub is called immediately after mg_relax which does a full halo exchange. If used in other location the halo shift might need to be adjusted More... | |
real(wp) function | compute_norm (mg, component, scale) |
subroutine | return_betamu_gamma_ion_concentrations (n, betamu, ion_concentration, neutralisation_ratios, steric_weight_average) |
needed in dl_mg_solver_pbe if the arguments betamu or gammas are present set them to the values used in code added the ion concentrartions used internally (it might differ from the input value for counter-ion neutralisation, when using PBC) More... | |
subroutine | check_error_code (ierror) |
checks if the first error code bubbled up to the top level or if it was replaced More... | |
subroutine | mg_convergence_test (mg, level, test_tag, iter, converged, f_norm, res_norm, sol) |
convergence tests, for details see dl_mg_convergence_params More... | |
subroutine | set_boundary_values (mg, v, v_type, val, aval, av_type) |
subroutine | set_mass_term (mg, der_pot, ierr) |
computes mass term at non-zero potential used in linearized PBEs (Newton and deferred correction) More... | |
auxiliary subroutines
subroutine dl_mg_utils::check_error_code | ( | integer, intent(inout), optional | ierror | ) |
checks if the first error code bubbled up to the top level or if it was replaced
real(wp) function dl_mg_utils::compute_norm | ( | type(mg_t), intent(in), target | mg, |
character(len=1), intent(in) | component, | ||
logical, intent(in), optional | scale | ||
) |
subroutine dl_mg_utils::compute_rho_sum | ( | real(wp), dimension(isx:,isy:,isz:), intent(in) | rho | ) |
subroutine dl_mg_utils::compute_steric_weight_sum | ( | logical, intent(out) | has_stw, |
real(wp), dimension(isx:, isy:, isz:), intent(in), optional | steric_weight | ||
) |
computes the accesible volume
\[ ~ \int_V steric\_weight(r) dv / (V/(nx*ny*nz)) \]
Attention: this sum is not multipliplied with dv = dx * dy * dz because it's used as a normalisation factor for quatites which are also not multiplied with dv
subroutine dl_mg_utils::dl_mg_init_nonlin | ( | integer, intent(out) | eq_type, |
real(wp), intent(in) | lambda, | ||
real(wp), intent(in) | temp, | ||
integer, intent(in) | n, | ||
real(wp), dimension(n), intent(in) | c, | ||
real(wp), dimension(n), intent(in) | q, | ||
real(wp), dimension(:,:,:), intent(in) | rho, | ||
real(wp), dimension(:,:,:), intent(in), optional | steric_weight, | ||
real(wp), dimension(:,:,:), intent(in), optional | der_pot, | ||
real(wp), intent(in), optional | expcap, | ||
logical, intent(in), optional | linearised, | ||
logical, intent(in), optional | use_fas, | ||
integer, intent(in), optional | neutralisation_method, | ||
real(wp), dimension(n), intent(in), optional | neutralisation_ion_ratios, | ||
integer, intent(out), optional | ierror | ||
) |
Initialises the non-linear parameters. It can be called multiple times, e.g. if beta needs to be changed in continuation calculation.
This subroutine needes to be called also in the case of linearised PBE.
Note: from physical point of view it is sensible to assume that the solvent is electrically neutral, i.e. \( \sum_i c_i q_i = 0 \), however the neutrality condition is not checked.
for arguments description see dl_mg::dl_mg_solver_pbe
subroutine dl_mg_utils::mg_convergence_test | ( | type (mg_t), intent(inout), target | mg, |
integer, intent(in) | level, | ||
integer, intent(in) | test_tag, | ||
integer, intent(in) | iter, | ||
integer, intent(out) | converged, | ||
real(wp), intent(in), optional | f_norm, | ||
real(wp), intent(in), optional | res_norm, | ||
real(wp), dimension(isx:,isy:,isz:), intent(in), optional, target | sol | ||
) |
convergence tests, for details see dl_mg_convergence_params
subroutine dl_mg_utils::remove_zero_mode_potential | ( | type(mg_t), intent(inout), target | mg, |
integer, intent(in) | eq_type | ||
) |
removes zero mode from the potential. Needed when using full PBC at the moment this is called only inside vcycle for mgp it migh need to be extend for arbitrary arrays with optional args N.B. the potential shift is applied also to the halo points to save a halo exchange. This is so because this sub is called immediately after mg_relax which does a full halo exchange. If used in other location the halo shift might need to be adjusted
subroutine dl_mg_utils::remove_zero_mode_source | ( | type(mg_t), intent(inout), target | mg, |
integer, intent(in) | eq_type, | ||
integer, intent(in) | remove_mode, | ||
real(wp), dimension(isx:, isy:, isz:), intent(inout), optional, target | rho, | ||
real(wp), dimension(isx:, isy:, isz:), intent(in), optional, target | steric_weight, | ||
real(wp), dimension(isx:, isy:, isz:), intent(in), optional | der_pot | ||
) |
removes zero modes from source term. Needed when using full PBC
[in] | remove_mode | 1 (first call) or 2 for call in inner loops (defco or newton) |
subroutine dl_mg_utils::return_betamu_gamma_ion_concentrations | ( | integer, intent(in) | n, |
real(wp), dimension(n), intent(out), optional | betamu, | ||
real(wp), dimension(n), intent(out), optional | ion_concentration, | ||
real(wp), dimension(:), intent(out), optional | neutralisation_ratios, | ||
real(wp), intent(out), optional | steric_weight_average | ||
) |
needed in dl_mg_solver_pbe if the arguments betamu or gammas are present set them to the values used in code added the ion concentrartions used internally (it might differ from the input value for counter-ion neutralisation, when using PBC)
subroutine dl_mg_utils::set_boundary_values | ( | type(mg_t), intent(inout) | mg, |
real(wp), dimension(:,:,:), intent(inout) | v, | ||
integer, intent(in) | v_type, | ||
real(wp), intent(in), optional | val, | ||
real(wp), dimension(:,:,:), intent(in), optional | aval, | ||
integer, intent(in), optional | av_type | ||
) |
subroutine dl_mg_utils::set_mass_term | ( | type(mg_t), dimension(:), intent(inout) | mg, |
real(wp), dimension(isx:, isy:, isz:), intent(in) | der_pot, | ||
integer, intent(out), optional | ierr | ||
) |
computes mass term at non-zero potential used in linearized PBEs (Newton and deferred correction)
subroutine dl_mg_utils::test_charge_neutrality | ( | real(wp), dimension(:,:,:), intent(in) | rho, |
real(wp), dimension(:,:,:), intent(in) | pot, | ||
real(wp), dimension(:,:,:), intent(in), optional | steric_weight | ||
) |
called by the PBE solver