ddx 0.6.0
Libary for domain-decomposition methods for polarizable continuum models
|
Core routines and parameters specific to gradients. More...
Functions/Subroutines | |
subroutine | contract_grad_l (params, constants, isph, sigma, xi, basloc, dbsloc, vplm, vcos, vsin, fx) |
Compute the gradients of the ddCOSMO matrix. More... | |
subroutine | contract_gradi_lik (params, constants, isph, sigma, xi, basloc, dbsloc, vplm, vcos, vsin, fx) |
Contribution to the gradients of the ddCOSMO matrix. More... | |
subroutine | contract_gradi_lji (params, constants, isph, sigma, xi, basloc, dbsloc, vplm, vcos, vsin, fx) |
Contribution to the gradients of the ddCOSMO matrix. More... | |
subroutine | contract_grad_u (params, constants, isph, xi, phi, fx) |
Gradient of the characteristic function U. More... | |
subroutine | contract_grad_b (params, constants, isph, Xe, Xadj_e, force) |
Subroutine to compute contraction of B matrix. More... | |
subroutine | contract_grad_c (params, constants, workspace, Xr, Xe, Xadj_r_sgrid, Xadj_e_sgrid, Xadj_r, Xadj_e, force, diff_re, ddx_error) |
Subroutine to compute contraction of C matrix. More... | |
subroutine | contract_grad_f (params, constants, workspace, sol_adj, sol_sgrid, gradpsi, normal_hessian_cav, force, state, ddx_error) |
Subroutine to compute contraction of F0. More... | |
subroutine | contract_gradi_bik (params, constants, isph, Xe, Xadj_e, force) |
Subroutine to compute K^A counterpart for the HSP equation. Similar to contract_gradi_Lik. More... | |
subroutine | contract_gradi_bji (params, constants, isph, Xe, Xadj_e, force) |
Subroutine to compute K^A+K^C counterpart for the HSP equation. Similar to contract_gradi_Lji. More... | |
subroutine | contract_grad_c_worker1 (params, constants, workspace, Xadj_r_sgrid, Xadj_e_sgrid, diff_re, force, ddx_error) |
Subroutine to calculate the third derivative term in C1_C2 matrix, namely the derivative of PU_i. More... | |
subroutine | contract_grad_c_worker2 (params, constants, workspace, Xr, Xe, Xadj_r_sgrid, Xadj_e_sgrid, Xadj_r, Xadj_e, force, diff_re, ddx_error) |
Subroutine to compute the derivative of U_i(x_in) and \bf{k}^j_l0(x_in)Y^j_l0m0(x_in) More... | |
subroutine | contract_grad_f_worker1 (params, constants, workspace, sol_adj, sol_sgrid, gradpsi, force, ddx_error) |
Subroutine to compute the derivative of U_i(x_in) and \bf{k}^j_l0(x_in)Y^j_l0m0(x_in) in F0. More... | |
subroutine | contract_grad_f_worker2 (params, constants, gradpsi, normal_hessian_cav, force, state, ddx_error) |
subroutine | gradr_sph (params, constants, isph, vplm, vcos, vsin, basloc, dbsloc, g, ygrid, fx) |
Sphere contribution to the ddPCM matrix gradient using N^2 code. More... | |
subroutine | gradr_fmm (params, constants, workspace, g, ygrid, fx) |
Compute the ddPCM matrix gradient using FMMS (2 matvecs) More... | |
subroutine | gradr (params, constants, workspace, g, ygrid, fx) |
Gradient of the ddPCM matrix. More... | |
subroutine | gradr_dense (params, constants, workspace, g, ygrid, fx) |
Gradient of the ddPCM matrix using N^2 code. More... | |
subroutine | zeta_grad (params, constants, state, e_cav, forces) |
Force term: interaction of the external electric field with the zeta intermediate. This routine is called by the gradients of ddCOSMO ddPCM and ddLPB. More... | |
Core routines and parameters specific to gradients.
ddX software
subroutine ddx_gradients::contract_grad_l | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | sigma, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | xi, | ||
real(dp), dimension(constants % nbasis), intent(inout) | basloc, | ||
real(dp), dimension(3, constants % nbasis), intent(inout) | dbsloc, | ||
real(dp), dimension(constants % nbasis), intent(inout) | vplm, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vcos, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vsin, | ||
real(dp), dimension(3), intent(inout) | fx | ||
) |
Compute the gradients of the ddCOSMO matrix.
Definition at line 18 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_gradi_lik | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | sigma, | ||
real(dp), dimension(params % ngrid), intent(in) | xi, | ||
real(dp), dimension(constants % nbasis), intent(inout) | basloc, | ||
real(dp), dimension(3, constants % nbasis), intent(inout) | dbsloc, | ||
real(dp), dimension(constants % nbasis), intent(inout) | vplm, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vcos, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vsin, | ||
real(dp), dimension(3), intent(inout) | fx | ||
) |
Contribution to the gradients of the ddCOSMO matrix.
Definition at line 35 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_gradi_lji | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | sigma, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | xi, | ||
real(dp), dimension(constants % nbasis), intent(inout) | basloc, | ||
real(dp), dimension(3, constants % nbasis), intent(inout) | dbsloc, | ||
real(dp), dimension(constants % nbasis), intent(inout) | vplm, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vcos, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vsin, | ||
real(dp), dimension(3), intent(inout) | fx | ||
) |
Contribution to the gradients of the ddCOSMO matrix.
Definition at line 103 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_u | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | xi, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | phi, | ||
real(dp), dimension(3), intent(inout) | fx | ||
) |
Gradient of the characteristic function U.
Definition at line 210 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_b | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xe, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_e, | ||
real(dp), dimension(3), intent(inout) | force | ||
) |
Subroutine to compute contraction of B matrix.
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | isph | : Index of sphere |
[in] | Xe | : Solution vector Xe |
[in] | Xadj_e | : Adjoint solution on evaluated on grid points Xadj_e_sgrid |
[out] | force | : Force of adjoint part |
Definition at line 253 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_c | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xr, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xe, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_r_sgrid, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_e_sgrid, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xadj_r, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xadj_e, | ||
real(dp), dimension(3, params % nsph), intent(inout) | force, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(out) | diff_re, | ||
type(ddx_error_type), intent(inout) | ddx_error | ||
) |
Subroutine to compute contraction of C matrix.
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | workspace | : Input workspace |
[in] | Xr | : Solution of the Laplace problem |
[in] | Xe | : Solution of the HSP problem |
[in] | Xadj_r_sgrid | : Solution of the Adjoint Laplace problem evaluated at the grid |
[in] | Xadj_e_sgrid | : Solution of the Adjoint HSP problem evaluated at the grid |
[in] | Xadj_r | : Adjoint solution of the Laplace problem |
[in] | Xadj_e | : Adjoint solution of the HSP problem |
[in,out] | force | : Force |
[out] | diff_re | : epsilon_1/epsilon_2 * l'/r_j[Xr]_jl'm'
|
[in,out] | ddx_error | ddX error |
Definition at line 284 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_f | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | sol_adj, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | sol_sgrid, | ||
real(dp), dimension(3, constants % ncav), intent(in) | gradpsi, | ||
real(dp), dimension(3, constants % ncav), intent(in) | normal_hessian_cav, | ||
real(dp), dimension(3, params % nsph), intent(inout) | force, | ||
type(ddx_state_type), intent(inout) | state, | ||
type(ddx_error_type), intent(inout) | ddx_error | ||
) |
Subroutine to compute contraction of F0.
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | workspace | : Input workspace |
[in] | sol_sgrid | : Solution of the Adjoint problem evaluated at the grid |
[in] | sol_adj | : Adjoint solution |
[in] | normal_hessian_cav | : Normal of the Hessian evaluated at cavity points |
[in] | icav_g | : Index of outside cavity point |
[out] | force | : Force |
[in,out] | ddx_error | ddX error |
Definition at line 325 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_gradi_bik | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xe, | ||
real(dp), dimension(params % ngrid), intent(in) | Xadj_e, | ||
real(dp), dimension(3), intent(inout) | force | ||
) |
Subroutine to compute K^A counterpart for the HSP equation. Similar to contract_gradi_Lik.
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | isph | : Index of sphere |
[in] | Xe | : Solution vector Xe |
[in] | Xadj_e | : Adjoint solution on evaluated on grid points Xadj_e_sgrid |
[out] | force | : Force |
Definition at line 364 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_gradi_bji | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xe, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_e, | ||
real(dp), dimension(3), intent(inout) | force | ||
) |
Subroutine to compute K^A+K^C counterpart for the HSP equation. Similar to contract_gradi_Lji.
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | workspace | : Input workspace |
[in] | isph | : Index of sphere |
[in] | Xe | : Solution vector Xe |
[in] | Xadj_e | : Adjoint solution on evaluated on grid points Xadj_e_sgrid |
[out] | force | : Force of adjoint part |
Definition at line 452 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_c_worker1 | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_r_sgrid, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_e_sgrid, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | diff_re, | ||
real(dp), dimension(3, params % nsph), intent(inout) | force, | ||
type(ddx_error_type), intent(inout) | ddx_error | ||
) |
Subroutine to calculate the third derivative term in C1_C2 matrix, namely the derivative of PU_i.
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | workspace | : Input workspace |
[in] | Xadj_r_sgrid | : Adjoint Laplace solution evaluated at grid point |
[in] | Xadj_e_sgrid | : Adjoint HSP solution evaluated at grid point |
[in] | diff_re | : l'/r_j[Xr]_jl'm' -(i'_l'(r_j)/i_l'(r_j))[Xe]_jl'm' |
[out] | force | : Force |
[in,out] | ddx_error | ddX error |
Definition at line 588 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_c_worker2 | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xr, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xe, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_r_sgrid, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | Xadj_e_sgrid, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xadj_r, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | Xadj_e, | ||
real(dp), dimension(3, params % nsph), intent(inout) | force, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(out) | diff_re, | ||
type(ddx_error_type), intent(inout) | ddx_error | ||
) |
Subroutine to compute the derivative of U_i(x_in) and \bf{k}^j_l0(x_in)Y^j_l0m0(x_in)
contract_grad_C_worker2 : Force Derivative of U_i^e(x_in), k_l0, and Y_l0m0
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | workspace | : Input workspace |
[in] | ksph | : Derivative with respect to x_k |
[in] | Xr | : Solution of the Laplace problem |
[in] | Xe | : Solution of the HSP problem |
[in] | Xadj_r_sgrid | : Solution of the Adjoint Laplace problem evaluated at the grid |
[in] | Xadj_e_sgrid | : Solution of the Adjoint HSP problem evaluated at the grid |
[in] | Xadj_r | : Adjoint solution of the Laplace problem |
[in] | Xadj_e | : Adjoint solution of the HSP problem |
[in,out] | force | : Force |
[out] | diff_re | : epsilon_1/epsilon_2 * l'/r_j[Xr]_jl'm'
|
[in,out] | ddx_error | ddX error |
Definition at line 869 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_f_worker1 | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | sol_adj, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | sol_sgrid, | ||
real(dp), dimension(3, constants % ncav), intent(in) | gradpsi, | ||
real(dp), dimension(3, params % nsph), intent(inout) | force, | ||
type(ddx_error_type), intent(inout) | ddx_error | ||
) |
Subroutine to compute the derivative of U_i(x_in) and \bf{k}^j_l0(x_in)Y^j_l0m0(x_in) in F0.
contract_grad_f_worker1 : Force Derivative of U_i^e(x_in), k_l0, and Y_l0m0 in F0
[in] | params | : Input parameter file |
[in] | constants | : Input constants file |
[in] | workspace | : Input workspace |
[in] | ksph | : Derivative with respect to x_k |
[in] | sol_sgrid | : Solution of the Adjoint problem evaluated at the grid |
[in] | sol_adj | : Adjoint solution |
[in] | gradpsi | : Gradient of Psi_0 |
[in,out] | force | : Force |
[in,out] | ddx_error | ddX error |
Definition at line 1242 of file ddx_gradients.f90.
subroutine ddx_gradients::contract_grad_f_worker2 | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
real(dp), dimension(3, constants % ncav), intent(in) | gradpsi, | ||
real(dp), dimension(3, constants % ncav), intent(in) | normal_hessian_cav, | ||
real(dp), dimension(3, params % nsph), intent(inout) | force, | ||
type(ddx_state_type), intent(inout) | state, | ||
type(ddx_error_type), intent(inout) | ddx_error | ||
) |
Definition at line 1581 of file ddx_gradients.f90.
subroutine ddx_gradients::gradr_sph | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
integer, intent(in) | isph, | ||
real(dp), dimension(constants % nbasis), intent(inout) | vplm, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vcos, | ||
real(dp), dimension(params % lmax+1), intent(inout) | vsin, | ||
real(dp), dimension(constants % nbasis), intent(inout) | basloc, | ||
real(dp), dimension(3, constants % nbasis), intent(inout) | dbsloc, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | g, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | ygrid, | ||
real(dp), dimension(3), intent(inout) | fx | ||
) |
Sphere contribution to the ddPCM matrix gradient using N^2 code.
Definition at line 1636 of file ddx_gradients.f90.
subroutine ddx_gradients::gradr_fmm | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | g, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | ygrid, | ||
real(dp), dimension(3, params % nsph), intent(out) | fx | ||
) |
Compute the ddPCM matrix gradient using FMMS (2 matvecs)
Definition at line 1941 of file ddx_gradients.f90.
subroutine ddx_gradients::gradr | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | g, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | ygrid, | ||
real(dp), dimension(3, params % nsph), intent(out) | fx | ||
) |
Gradient of the ddPCM matrix.
Definition at line 2174 of file ddx_gradients.f90.
subroutine ddx_gradients::gradr_dense | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_workspace_type), intent(inout) | workspace, | ||
real(dp), dimension(constants % nbasis, params % nsph), intent(in) | g, | ||
real(dp), dimension(params % ngrid, params % nsph), intent(in) | ygrid, | ||
real(dp), dimension(3, params % nsph), intent(out) | fx | ||
) |
Gradient of the ddPCM matrix using N^2 code.
Definition at line 2194 of file ddx_gradients.f90.
subroutine ddx_gradients::zeta_grad | ( | type(ddx_params_type), intent(in) | params, |
type(ddx_constants_type), intent(in) | constants, | ||
type(ddx_state_type), intent(inout) | state, | ||
real(dp), dimension(3, constants % ncav), intent(in) | e_cav, | ||
real(dp), dimension(3, params % nsph), intent(inout) | forces | ||
) |
Force term: interaction of the external electric field with the zeta intermediate. This routine is called by the gradients of ddCOSMO ddPCM and ddLPB.
Definition at line 2219 of file ddx_gradients.f90.