Computes the electric field of a trial set of induced point dipoles at polarizable sites. This is intended to be used as matrix-vector routine in the solution of the linear system.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(ommp_electrostatics_type), | intent(inout) | :: | eel |
Data structure for electrostatic part of the system |
||
character, | intent(in) | :: | in_kind | |||
logical, | intent(in) | :: | do_V |
Flag to control which properties have to be computed. |
||
logical, | intent(in) | :: | do_E |
Flag to control which properties have to be computed. |
||
logical, | intent(in) | :: | do_Egrd |
Flag to control which properties have to be computed. |
||
logical, | intent(in) | :: | do_EHes |
Flag to control which properties have to be computed. |
subroutine elec_prop_D2D(eel, in_kind, do_V, do_E, do_Egrd, do_EHes)
!! Computes the electric field of a trial set of induced point dipoles
!! at polarizable sites. This is intended to be used as matrix-vector
!! routine in the solution of the linear system.
implicit none
type(ommp_electrostatics_type), intent(inout) :: eel
!! Data structure for electrostatic part of the system
logical, intent(in) :: do_V, do_E, do_Egrd, do_EHes
!! Flag to control which properties have to be computed.
character, intent(in) :: in_kind
integer(ip) :: i, j, jpol, ipol, ij, idx, ikernel, knd
logical :: to_scale, to_do
real(rp) :: kernel(5), dr(3), tmpV, tmpE(3), tmpEgr(6), tmpHE(10), scalf
knd = 1 ! Default
if(in_kind == 'P') then
knd = _amoeba_P_
elseif(in_kind == 'D') then
knd = _amoeba_D_
elseif(eel%amoeba) then
call fatal_error("Unrecognized interaction '"//in_kind//"' in elec&
&_prop_D2D.")
else
knd = 1
end if
if(do_EHes) then
ikernel = 5
elseif(do_Egrd) then
ikernel = 4
elseif(do_E) then
ikernel = 3
elseif(do_V) then
ikernel = 2
else
return
end if
if(eel%use_fmm) then
call prepare_fmm_ipd(eel, knd)
!$omp parallel do default(shared) schedule(dynamic) &
!$omp private(i,j,ij,ipol,jpol,idx,dr,kernel,to_do,to_scale,scalf,tmpV,tmpE,tmpEgr,tmpHE)
do ipol=1, eel%pol_atoms
i = eel%polar_mm(ipol)
call cart_propfar_at_ipart(eel%fmm_ipd(knd), i, &
do_V, eel%V_D2D(ipol,knd), &
do_E, eel%E_D2D(:,ipol,knd), &
do_Egrd, eel%Egrd_D2D(:,ipol,knd), &
do_EHes, eel%EHes_D2D(:,ipol,knd))
! Near field is computed internally because dumped kernel is required
do ij=eel%fmm_near_field_list%ri(i), eel%fmm_near_field_list%ri(i+1)-1
j = eel%fmm_near_field_list%ci(ij)
jpol = eel%polar_mm(j)
! If the atom is not polarizable, skip
if(jpol < 1) cycle
!loop on target
to_do = .true.
to_scale = .false.
scalf = 1.0
! Check if the element should be scaled
do idx=eel%list_P_P%ri(ipol), eel%list_P_P%ri(ipol+1)-1
if(eel%list_P_P%ci(idx) == jpol) then
to_scale = .true.
exit
end if
end do
!If it should set the correct variables
if(to_scale) then
to_do = eel%todo_P_P(idx)
scalf = eel%scalef_P_P(idx)
end if
if(to_do) then
call damped_coulomb_kernel(eel, j, i,&
ikernel, kernel, dr)
if(do_V) tmpV = 0.0_rp
if(do_E) tmpE = 0.0_rp
if(do_Egrd) tmpEgr = 0.0_rp
if(do_EHes) tmpHE = 0.0_rp
call mu_elec_prop(eel%ipd(:,jpol,knd), dr, kernel, &
do_V, tmpV, &
do_E, tmpE, &
do_Egrd, tmpEgr, &
do_EHes, tmpHE)
if(to_scale) then
if(do_V) eel%V_D2D(ipol,knd) = eel%V_D2D(ipol,knd) + tmpV * scalf
if(do_E) eel%E_D2D(:, ipol,knd) = eel%E_D2D(:, ipol,knd) + tmpE * scalf
if(do_Egrd) eel%Egrd_D2D(:, ipol,knd) = eel%Egrd_D2D(:, ipol,knd) + tmpEgr * scalf
if(do_EHes) eel%EHes_D2D(:, ipol,knd) = eel%EHes_D2D(:, ipol,knd) + tmpHE * scalf
else
if(do_V) eel%V_D2D(ipol,knd) = eel%V_D2D(ipol,knd) + tmpV
if(do_E) eel%E_D2D(:, ipol,knd) = eel%E_D2D(:, ipol,knd) + tmpE
if(do_Egrd) eel%Egrd_D2D(:, ipol,knd) = eel%Egrd_D2D(:, ipol,knd) + tmpEgr
if(do_EHes) eel%EHes_D2D(:, ipol,knd) = eel%EHes_D2D(:, ipol,knd) + tmpHE
end if
end if
end do
end do
if(allocated(eel%list_P_P_fmm_far)) then
! Now remove screened interactions from far-field
do ipol=1, eel%pol_atoms
i = eel%polar_mm(ipol)
do idx=eel%list_P_P_fmm_far%ri(ipol), eel%list_P_P_fmm_far%ri(ipol+1)-1
jpol = eel%list_P_P_fmm_far%ci(idx)
j = eel%polar_mm(jpol)
scalf = 1.0 - eel%scalef_P_P_fmm_far(idx)
call damped_coulomb_kernel(eel, j, i,&
2, kernel(1:3), dr)
if(do_V) tmpV = 0.0_rp
if(do_E) tmpE = 0.0_rp
if(do_Egrd) tmpEgr = 0.0_rp
if(do_EHes) tmpHE = 0.0_rp
call mu_elec_prop(eel%ipd(:,jpol,knd), dr, kernel, &
do_V, tmpV, &
do_E, tmpE, &
do_Egrd, tmpEgr, &
do_EHes, tmpHE)
if(do_V) eel%V_D2D(ipol,knd) = eel%V_D2D(ipol,knd) + tmpV * scalf
if(do_E) eel%E_D2D(:, ipol,knd) = eel%E_D2D(:, ipol,knd) + tmpE * scalf
if(do_Egrd) eel%Egrd_D2D(:, ipol,knd) = eel%Egrd_D2D(:, ipol,knd) + tmpEgr * scalf
if(do_EHes) eel%EHes_D2D(:, ipol,knd) = eel%EHes_D2D(:, ipol,knd) + tmpHE * scalf
end do
end do
end if
else
!$omp parallel do default(shared) schedule(dynamic) &
!$omp private(i,j,idx,dr,kernel,to_do,to_scale,scalf,tmpV,tmpE,tmpEgr,tmpHE)
do j=1, eel%pol_atoms
do i=1, eel%pol_atoms
if(j == i) cycle
!loop on target
to_do = .true.
to_scale = .false.
scalf = 1.0
! Check if the element should be scaled
do idx=eel%list_P_P%ri(i), eel%list_P_P%ri(i+1)-1
if(eel%list_P_P%ci(idx) == j) then
to_scale = .true.
exit
end if
end do
!If it should set the correct variables
if(to_scale) then
to_do = eel%todo_P_P(idx)
scalf = eel%scalef_P_P(idx)
end if
if(to_do) then
call damped_coulomb_kernel(eel, eel%polar_mm(i), &
eel%polar_mm(j),&
ikernel, kernel, dr)
if(do_V) tmpV = 0.0_rp
if(do_E) tmpE = 0.0_rp
if(do_Egrd) tmpEgr = 0.0_rp
if(do_EHes) tmpHE = 0.0_rp
call mu_elec_prop(eel%ipd(:,i,knd), dr, kernel, &
do_V, tmpV, &
do_E, tmpE, &
do_Egrd, tmpEgr, &
do_EHes, tmpHE)
if(to_scale) then
if(do_V) eel%V_D2D(j,knd) = eel%V_D2D(j,knd) + tmpV * scalf
if(do_E) eel%E_D2D(:, j,knd) = eel%E_D2D(:, j,knd) + tmpE * scalf
if(do_Egrd) eel%Egrd_D2D(:, j,knd) = eel%Egrd_D2D(:, j,knd) + tmpEgr * scalf
if(do_EHes) eel%EHes_D2D(:, j,knd) = eel%EHes_D2D(:, j,knd) + tmpHE * scalf
else
if(do_V) eel%V_D2D(j,knd) = eel%V_D2D(j,knd) + tmpV
if(do_E) eel%E_D2D(:, j,knd) = eel%E_D2D(:, j,knd) + tmpE
if(do_Egrd) eel%Egrd_D2D(:, j,knd) = eel%Egrd_D2D(:, j,knd) + tmpEgr
if(do_EHes) eel%EHes_D2D(:, j,knd) = eel%EHes_D2D(:, j,knd) + tmpHE
end if
end if
end do
end do
end if
end subroutine elec_prop_D2D