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 |
real(kind=dp), save :: ecr = 0.0_DP |
16 |
real(kind=dp), save :: rt = 0.0_DP |
17 |
real(kind=dp), save :: pre = 0.0_DP |
18 |
logical, save :: haveCutoffs = .false. |
19 |
logical, save :: haveChargeMap = .false. |
20 |
|
21 |
public::setCutoffsCharge |
22 |
public::do_charge_pair |
23 |
|
24 |
type :: ChargeList |
25 |
real(kind=DP) :: charge = 0.0_DP |
26 |
end type ChargeList |
27 |
|
28 |
type(ChargeList), dimension(:), allocatable :: ChargeMap |
29 |
|
30 |
contains |
31 |
|
32 |
subroutine setCutoffsCharge(this_ecr, this_rt) |
33 |
real(kind=dp), intent(in) :: this_ecr, this_rt |
34 |
ecr = this_ecr |
35 |
rt = this_rt |
36 |
|
37 |
! pre converts from fundamental charge to kcal/mol |
38 |
pre = 332.06508_DP |
39 |
|
40 |
haveCutoffs = .true. |
41 |
|
42 |
return |
43 |
end subroutine setCutoffsCharge |
44 |
|
45 |
subroutine createChargeMap(status) |
46 |
integer :: nAtypes |
47 |
integer :: status |
48 |
integer :: i |
49 |
real (kind=DP) :: thisCharge |
50 |
logical :: thisProperty |
51 |
|
52 |
status = 0 |
53 |
|
54 |
nAtypes = getSize(atypes) |
55 |
|
56 |
if (nAtypes == 0) then |
57 |
status = -1 |
58 |
return |
59 |
end if |
60 |
|
61 |
if (.not. allocated(ChargeMap)) then |
62 |
allocate(ChargeMap(nAtypes)) |
63 |
endif |
64 |
|
65 |
do i = 1, nAtypes |
66 |
|
67 |
call getElementProperty(atypes, i, "is_Charge", thisProperty) |
68 |
|
69 |
if (thisProperty) then |
70 |
call getElementProperty(atypes, i, "charge", thisCharge) |
71 |
ChargeMap(i)%charge = thisCharge |
72 |
endif |
73 |
|
74 |
end do |
75 |
|
76 |
haveChargeMap = .true. |
77 |
|
78 |
end subroutine createChargeMap |
79 |
|
80 |
|
81 |
subroutine do_charge_pair(atom1, atom2, d, rij, r2, pot, f, & |
82 |
do_pot, do_stress) |
83 |
|
84 |
logical :: do_pot, do_stress |
85 |
|
86 |
integer atom1, atom2, me1, me2, id1, id2 |
87 |
integer :: localError |
88 |
real(kind=dp) :: rij, r2, q1, q2 |
89 |
real(kind=dp) :: drdx, drdy, drdz, dudr, fx, fy, fz |
90 |
real(kind=dp) :: taper, dtdr, vterm |
91 |
|
92 |
real( kind = dp ) :: pot |
93 |
real( kind = dp ), dimension(3) :: d |
94 |
real( kind = dp ), dimension(3,nLocal) :: f |
95 |
|
96 |
if (.not. haveCutoffs) then |
97 |
write(default_error,*) 'charge-charge does not have cutoffs set!' |
98 |
return |
99 |
endif |
100 |
|
101 |
if (.not.haveChargeMap) then |
102 |
localError = 0 |
103 |
call createChargeMap(localError) |
104 |
if ( localError .ne. 0 ) then |
105 |
call handleError("charge-charge", "ChargeMap creation failed!") |
106 |
return |
107 |
end if |
108 |
endif |
109 |
|
110 |
#ifdef IS_MPI |
111 |
me1 = atid_Row(atom1) |
112 |
me2 = atid_Col(atom2) |
113 |
#else |
114 |
me1 = atid(atom1) |
115 |
me2 = atid(atom2) |
116 |
#endif |
117 |
|
118 |
|
119 |
q1 = ChargeMap(me1)%charge |
120 |
q2 = ChargeMap(me2)%charge |
121 |
|
122 |
if (rij.le.ecr) then |
123 |
|
124 |
if (rij.lt.rt) then |
125 |
taper = 1.0d0 |
126 |
dtdr = 0.0d0 |
127 |
else |
128 |
taper = (ecr + 2.0d0*rij - 3.0d0*rt)*(ecr-rij)**2/ ((ecr-rt)**3) |
129 |
dtdr = 6.0d0*(rij*rij - rij*rt - rij*ecr +ecr*rt)/((ecr-rt)**3) |
130 |
endif |
131 |
|
132 |
vterm = pre * q1 * q2 / rij |
133 |
|
134 |
!!if(rij .le. 1.5) then |
135 |
!! write(*,*) atom1, atom2, q1, q2, rij, vterm |
136 |
!!endif |
137 |
|
138 |
dudr = vterm * ( dtdr - taper / rij ) |
139 |
|
140 |
drdx = d(1) / rij |
141 |
drdy = d(2) / rij |
142 |
drdz = d(3) / rij |
143 |
|
144 |
fx = dudr * drdx |
145 |
fy = dudr * drdy |
146 |
fz = dudr * drdz |
147 |
|
148 |
|
149 |
#ifdef IS_MPI |
150 |
if (do_pot) then |
151 |
pot_Row(atom1) = pot_Row(atom1) + vterm*taper*0.5 |
152 |
pot_Col(atom2) = pot_Col(atom2) + vterm*taper*0.5 |
153 |
endif |
154 |
|
155 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
156 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
157 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
158 |
|
159 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
160 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
161 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
162 |
|
163 |
#else |
164 |
|
165 |
if (do_pot) pot = pot + vterm*taper |
166 |
|
167 |
f(1,atom1) = f(1,atom1) + fx |
168 |
f(2,atom1) = f(2,atom1) + fy |
169 |
f(3,atom1) = f(3,atom1) + fz |
170 |
|
171 |
f(1,atom2) = f(1,atom2) - fx |
172 |
f(2,atom2) = f(2,atom2) - fy |
173 |
f(3,atom2) = f(3,atom2) - fz |
174 |
|
175 |
#endif |
176 |
|
177 |
if (do_stress) then |
178 |
|
179 |
#ifdef IS_MPI |
180 |
id1 = tagRow(atom1) |
181 |
id2 = tagColumn(atom2) |
182 |
#else |
183 |
id1 = atom1 |
184 |
id2 = atom2 |
185 |
#endif |
186 |
|
187 |
|
188 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
189 |
|
190 |
! because the d vector is the rj - ri vector, and |
191 |
! because fx, fy, fz are the force on atom i, we need a |
192 |
! negative sign here: |
193 |
|
194 |
tau_Temp(1) = tau_Temp(1) - d(1) * fx |
195 |
tau_Temp(2) = tau_Temp(2) - d(1) * fy |
196 |
tau_Temp(3) = tau_Temp(3) - d(1) * fz |
197 |
tau_Temp(4) = tau_Temp(4) - d(2) * fx |
198 |
tau_Temp(5) = tau_Temp(5) - d(2) * fy |
199 |
tau_Temp(6) = tau_Temp(6) - d(2) * fz |
200 |
tau_Temp(7) = tau_Temp(7) - d(3) * fx |
201 |
tau_Temp(8) = tau_Temp(8) - d(3) * fy |
202 |
tau_Temp(9) = tau_Temp(9) - d(3) * fz |
203 |
|
204 |
virial_Temp = virial_Temp + & |
205 |
(tau_Temp(1) + tau_Temp(5) + tau_Temp(9)) |
206 |
|
207 |
endif |
208 |
endif |
209 |
|
210 |
endif |
211 |
return |
212 |
end subroutine do_charge_pair |
213 |
|
214 |
end module charge_charge |