| 1 |
module charge_charge |
| 2 |
|
| 3 |
use force_globals |
| 4 |
use definitions |
| 5 |
use atype_module |
| 6 |
use vector_class |
| 7 |
use simulation |
| 8 |
use status |
| 9 |
#ifdef IS_MPI |
| 10 |
use mpiSimulation |
| 11 |
#endif |
| 12 |
implicit none |
| 13 |
|
| 14 |
PRIVATE |
| 15 |
|
| 16 |
real(kind=dp), parameter :: pre = 332.06508_DP |
| 17 |
logical, save :: haveChargeMap = .false. |
| 18 |
|
| 19 |
public::do_charge_pair |
| 20 |
|
| 21 |
type :: ChargeList |
| 22 |
real(kind=DP) :: charge = 0.0_DP |
| 23 |
end type ChargeList |
| 24 |
|
| 25 |
type(ChargeList), dimension(:), allocatable :: ChargeMap |
| 26 |
|
| 27 |
contains |
| 28 |
|
| 29 |
subroutine createChargeMap(status) |
| 30 |
integer :: nAtypes |
| 31 |
integer :: status |
| 32 |
integer :: i |
| 33 |
real (kind=DP) :: thisCharge |
| 34 |
logical :: thisProperty |
| 35 |
|
| 36 |
status = 0 |
| 37 |
|
| 38 |
nAtypes = getSize(atypes) |
| 39 |
|
| 40 |
if (nAtypes == 0) then |
| 41 |
status = -1 |
| 42 |
return |
| 43 |
end if |
| 44 |
|
| 45 |
if (.not. allocated(ChargeMap)) then |
| 46 |
allocate(ChargeMap(nAtypes)) |
| 47 |
endif |
| 48 |
|
| 49 |
do i = 1, nAtypes |
| 50 |
|
| 51 |
call getElementProperty(atypes, i, "is_Charge", thisProperty) |
| 52 |
|
| 53 |
if (thisProperty) then |
| 54 |
call getElementProperty(atypes, i, "charge", thisCharge) |
| 55 |
ChargeMap(i)%charge = thisCharge |
| 56 |
endif |
| 57 |
|
| 58 |
end do |
| 59 |
|
| 60 |
haveChargeMap = .true. |
| 61 |
|
| 62 |
end subroutine createChargeMap |
| 63 |
|
| 64 |
subroutine do_charge_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
| 65 |
pot, f, do_pot) |
| 66 |
|
| 67 |
logical :: do_pot |
| 68 |
|
| 69 |
integer atom1, atom2, me1, me2, id1, id2 |
| 70 |
integer :: localError |
| 71 |
real(kind=dp) :: rij, r2, q1, q2, sw, vpair |
| 72 |
real(kind=dp) :: drdx, drdy, drdz, dudr, fx, fy, fz |
| 73 |
real(kind=dp) :: vterm |
| 74 |
|
| 75 |
real( kind = dp ) :: pot |
| 76 |
real( kind = dp ), dimension(3) :: d, fpair |
| 77 |
real( kind = dp ), dimension(3,nLocal) :: f |
| 78 |
|
| 79 |
if (.not.haveChargeMap) then |
| 80 |
localError = 0 |
| 81 |
call createChargeMap(localError) |
| 82 |
if ( localError .ne. 0 ) then |
| 83 |
call handleError("charge-charge", "ChargeMap creation failed!") |
| 84 |
return |
| 85 |
end if |
| 86 |
endif |
| 87 |
|
| 88 |
#ifdef IS_MPI |
| 89 |
me1 = atid_Row(atom1) |
| 90 |
me2 = atid_Col(atom2) |
| 91 |
#else |
| 92 |
me1 = atid(atom1) |
| 93 |
me2 = atid(atom2) |
| 94 |
#endif |
| 95 |
|
| 96 |
q1 = ChargeMap(me1)%charge |
| 97 |
q2 = ChargeMap(me2)%charge |
| 98 |
|
| 99 |
vterm = pre * q1 * q2 / rij |
| 100 |
|
| 101 |
dudr = -sw * vterm / rij |
| 102 |
|
| 103 |
drdx = d(1) / rij |
| 104 |
drdy = d(2) / rij |
| 105 |
drdz = d(3) / rij |
| 106 |
|
| 107 |
fx = dudr * drdx |
| 108 |
fy = dudr * drdy |
| 109 |
fz = dudr * drdz |
| 110 |
|
| 111 |
vpair = vpair + vterm |
| 112 |
|
| 113 |
#ifdef IS_MPI |
| 114 |
if (do_pot) then |
| 115 |
pot_Row(atom1) = pot_Row(atom1) + sw*vterm*0.5 |
| 116 |
pot_Col(atom2) = pot_Col(atom2) + sw*vterm*0.5 |
| 117 |
endif |
| 118 |
|
| 119 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
| 120 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
| 121 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
| 122 |
|
| 123 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
| 124 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
| 125 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
| 126 |
|
| 127 |
#else |
| 128 |
|
| 129 |
if (do_pot) pot = pot + sw*vterm |
| 130 |
|
| 131 |
f(1,atom1) = f(1,atom1) + fx |
| 132 |
f(2,atom1) = f(2,atom1) + fy |
| 133 |
f(3,atom1) = f(3,atom1) + fz |
| 134 |
|
| 135 |
f(1,atom2) = f(1,atom2) - fx |
| 136 |
f(2,atom2) = f(2,atom2) - fy |
| 137 |
f(3,atom2) = f(3,atom2) - fz |
| 138 |
|
| 139 |
#endif |
| 140 |
|
| 141 |
#ifdef IS_MPI |
| 142 |
id1 = AtomRowToGlobal(atom1) |
| 143 |
id2 = AtomColToGlobal(atom2) |
| 144 |
#else |
| 145 |
id1 = atom1 |
| 146 |
id2 = atom2 |
| 147 |
#endif |
| 148 |
|
| 149 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
| 150 |
|
| 151 |
fpair(1) = fpair(1) + fx |
| 152 |
fpair(2) = fpair(2) + fy |
| 153 |
fpair(3) = fpair(3) + fz |
| 154 |
|
| 155 |
endif |
| 156 |
return |
| 157 |
end subroutine do_charge_pair |
| 158 |
|
| 159 |
end module charge_charge |