DL_MG
|
Higher order derivatives defect correction driver and its subroutines. More...
Functions/Subroutines | |
subroutine, public | dl_mg_defco_defect_corr_solver (eq_type, eps_full, eps_half, alpha, rho, fd_order, pot, der_pot, steric_weight, use_pot_in, use_damping, res, ierror) |
This subroutine computes the solution for the Poisson or Poisson-Boltzmann Equation with higher order approximation for the finite difference derivatives. It uses multigrid to obtain the solution for the second order finite differences and defect correction method to improve it to the specified order. Typically this increases solution accuracy but don't forget that higher order means does not imply higher precision in all cases. More... | |
subroutine | compute_current_defect (fd_order, eq_type, alpha, rho, eps, v_i, defect, defect_norm, defect_stpt, defect_norm_stpt, defect_bt, defect_norm_bt, steric_weight) |
subroutine | compute_boltzmann_defect_term (bt, pot, use_linear, steric_weight) |
Compute the Boltzmann term to add to the defect for the Poisson equation. More... | |
subroutine | compute_mu_defco (pot, steric_weight) |
computes the correction to chemical pontential in defco loop relevant when using PBE with PBC More... | |
subroutine | shift_pot_and_mu_defco (pot, steric_weight) |
subroutine | set_mg_grids_eps_and_steric (eq_type, eps_full, eps_half, steric_weight) |
allocate mg data structure for the mulitgrid solver and populate the data arrays for the fields that do not change during the computation: eps and steric weight More... | |
subroutine | free_mg |
Higher order derivatives defect correction driver and its subroutines.
James C. Womack, Lucian Anton 2016-17
subroutine dl_mg_defco::compute_boltzmann_defect_term | ( | real(kind=wp), dimension(isx:, isy:, isz:), intent(out) | bt, |
real(kind=wp), dimension(isx:, isy:, isz:), intent(in) | pot, | ||
logical, intent(in) | use_linear, | ||
real(kind=wp), dimension(isx:, isy:, isz:), intent(in), optional | steric_weight | ||
) |
Compute the Boltzmann term to add to the defect for the Poisson equation.
For the full non-linear PB equation the additional term in the defect is
\[ -\lambda \sum_{i} c_{i} q_{i} \beta W \exp( -\beta q_{i} \phi^{(i)})\]
and for the linearized PB equation, the additional term is
\[ \lambda {\lambda_{D}^{2}} W \phi^{(i)} \]
where \( c_{i}\) and \( q_{i} \) are the concentration and charge for ion type \(i\), \(\beta = 1/kT\), \( W\) is the steric weight \(\exp(-\beta V)]\), and \(\lambda_{D}\) is the Debye length.
If optional argument steric_weight is present, then the defect is evaluated using the steric potential, if not, then the steric potential is assumed to be 0.0 everywhere.
Adapted from implicit_solvent_boltzmann_defect from ONETEP's is_solvation_boltzmann module for inclusion in DL_MG. See original documentation in source.
J. C. Womack, 2016
subroutine dl_mg_defco::compute_current_defect | ( | integer, intent(in) | fd_order, |
integer, intent(in) | eq_type, | ||
real(kind=wp), intent(in) | alpha, | ||
real(kind=wp), dimension(1:,1:,1:), intent(in) | rho, | ||
real(kind=wp), dimension(1:,1:,1:), intent(in) | eps, | ||
real(kind=wp), dimension(1:,1:,1:), intent(in) | v_i, | ||
real(kind=wp), dimension(1:,1:,1:), intent(out) | defect, | ||
real(kind=wp), intent(out) | defect_norm, | ||
real(kind=wp), dimension(1:,1:,1:), intent(out) | defect_stpt, | ||
real(kind=wp), intent(out) | defect_norm_stpt, | ||
real(kind=wp), dimension(1:,1:,1:), intent(out) | defect_bt, | ||
real(kind=wp), intent(out) | defect_norm_bt, | ||
real(kind=wp), dimension(1:,1:,1:), intent(in), optional | steric_weight | ||
) |
subroutine dl_mg_defco::compute_mu_defco | ( | real(wp), dimension(isx:, isy:, isz:), intent(inout) | pot, |
real(wp), dimension(isx:, isy:, isz:), intent(in), optional | steric_weight | ||
) |
computes the correction to chemical pontential in defco loop relevant when using PBE with PBC
subroutine, public dl_mg_defco::dl_mg_defco_defect_corr_solver | ( | integer, intent(in) | eq_type, |
real(wp), dimension(:,:,:), intent(in) | eps_full, | ||
real(wp), dimension(:,:,:,:), intent(in) | eps_half, | ||
real(wp), intent(in) | alpha, | ||
real(wp), dimension(:,:,:), intent(in), target | rho, | ||
integer, intent(in) | fd_order, | ||
real(wp), dimension(:,:,:), intent(inout) | pot, | ||
real(wp), dimension(:,:,:), intent(in), optional | der_pot, | ||
real(wp), dimension(:,:,:), intent(in), optional | steric_weight, | ||
logical, intent(in), optional | use_pot_in, | ||
logical, intent(in), optional | use_damping, | ||
real(wp), dimension(:,:,:), intent(out), optional, target | res, | ||
integer, intent(out), optional | ierror | ||
) |
This subroutine computes the solution for the Poisson or Poisson-Boltzmann Equation with higher order approximation for the finite difference derivatives. It uses multigrid to obtain the solution for the second order finite differences and defect correction method to improve it to the specified order. Typically this increases solution accuracy but don't forget that higher order means does not imply higher precision in all cases.
[in] | eq_type | selector for the algorithm to be used |
[in] | eps_full | relative permitivity on the grid points ignored if the finite difference order is 2 |
[in] | eps_half | relative permitivity, the values correspond to the points located halfway between the grid points in each direction |
[in] | alpha | multiplicative constant for \(\rho\) defined by the units used |
[in] | rho | r.h.s. term (charge density) |
[in] | fd_order | finite difference order used in defect correction |
[in,out] | pot | the potential we seek |
[in] | der_pot | potential function at which the derivative of the Boltzmann term is computed (if not present defaults to 0) |
[in] | steric_weight | steric weight in Boltzmann term, i.e. \( e^{-\beta V_{steric}(r)} \) |
[in] | use_pot_in | start iterating from this potential, if not present starts from 0 inside the domain |
[in] | use_damping | if present and true use error damping in defect correction loop |
[out] | res | residual |
[out] | ierror | error flag |
subroutine dl_mg_defco::free_mg |
subroutine dl_mg_defco::set_mg_grids_eps_and_steric | ( | integer, intent(in) | eq_type, |
real(wp), dimension(isx:,isy:,isz:), intent(in) | eps_full, | ||
real(wp), dimension(isx:,isy:,isz:,:), intent(in), target | eps_half, | ||
real(wp), dimension(isx:,isy:,isz:), intent(in), optional | steric_weight | ||
) |
allocate mg data structure for the mulitgrid solver and populate the data arrays for the fields that do not change during the computation: eps and steric weight
subroutine dl_mg_defco::shift_pot_and_mu_defco | ( | real(wp), dimension(isx:, isy:, isz:), intent(inout) | pot, |
real(wp), dimension(isx:, isy:, isz:), intent(in), optional | steric_weight | ||
) |