1 |
gezelter |
411 |
!! |
2 |
|
|
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
|
|
!! |
4 |
|
|
!! The University of Notre Dame grants you ("Licensee") a |
5 |
|
|
!! non-exclusive, royalty free, license to use, modify and |
6 |
|
|
!! redistribute this software in source and binary code form, provided |
7 |
|
|
!! that the following conditions are met: |
8 |
|
|
!! |
9 |
|
|
!! 1. Acknowledgement of the program authors must be made in any |
10 |
|
|
!! publication of scientific results based in part on use of the |
11 |
|
|
!! program. An acceptable form of acknowledgement is citation of |
12 |
|
|
!! the article in which the program was described (Matthew |
13 |
|
|
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
|
|
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
|
|
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
|
|
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
|
|
!! |
18 |
|
|
!! 2. Redistributions of source code must retain the above copyright |
19 |
|
|
!! notice, this list of conditions and the following disclaimer. |
20 |
|
|
!! |
21 |
|
|
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
|
|
!! notice, this list of conditions and the following disclaimer in the |
23 |
|
|
!! documentation and/or other materials provided with the |
24 |
|
|
!! distribution. |
25 |
|
|
!! |
26 |
|
|
!! This software is provided "AS IS," without a warranty of any |
27 |
|
|
!! kind. All express or implied conditions, representations and |
28 |
|
|
!! warranties, including any implied warranty of merchantability, |
29 |
|
|
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
|
|
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
|
|
!! be liable for any damages suffered by licensee as a result of |
32 |
|
|
!! using, modifying or distributing the software or its |
33 |
|
|
!! derivatives. In no event will the University of Notre Dame or its |
34 |
|
|
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
|
|
!! direct, indirect, special, consequential, incidental or punitive |
36 |
|
|
!! damages, however caused and regardless of the theory of liability, |
37 |
|
|
!! arising out of the use of or inability to use software, even if the |
38 |
|
|
!! University of Notre Dame has been advised of the possibility of |
39 |
|
|
!! such damages. |
40 |
|
|
!! |
41 |
|
|
|
42 |
|
|
module electrostatic_module |
43 |
gezelter |
507 |
|
44 |
gezelter |
411 |
use force_globals |
45 |
|
|
use definitions |
46 |
|
|
use atype_module |
47 |
|
|
use vector_class |
48 |
|
|
use simulation |
49 |
|
|
use status |
50 |
|
|
#ifdef IS_MPI |
51 |
|
|
use mpiSimulation |
52 |
|
|
#endif |
53 |
|
|
implicit none |
54 |
|
|
|
55 |
|
|
PRIVATE |
56 |
|
|
|
57 |
gezelter |
602 |
#define __FORTRAN90 |
58 |
|
|
#include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" |
59 |
|
|
|
60 |
gezelter |
434 |
!! these prefactors convert the multipole interactions into kcal / mol |
61 |
|
|
!! all were computed assuming distances are measured in angstroms |
62 |
|
|
!! Charge-Charge, assuming charges are measured in electrons |
63 |
gezelter |
411 |
real(kind=dp), parameter :: pre11 = 332.0637778_dp |
64 |
gezelter |
434 |
!! Charge-Dipole, assuming charges are measured in electrons, and |
65 |
|
|
!! dipoles are measured in debyes |
66 |
|
|
real(kind=dp), parameter :: pre12 = 69.13373_dp |
67 |
|
|
!! Dipole-Dipole, assuming dipoles are measured in debyes |
68 |
|
|
real(kind=dp), parameter :: pre22 = 14.39325_dp |
69 |
|
|
!! Charge-Quadrupole, assuming charges are measured in electrons, and |
70 |
|
|
!! quadrupoles are measured in 10^-26 esu cm^2 |
71 |
|
|
!! This unit is also known affectionately as an esu centi-barn. |
72 |
|
|
real(kind=dp), parameter :: pre14 = 69.13373_dp |
73 |
gezelter |
411 |
|
74 |
gezelter |
602 |
!! variables to handle different summation methods for long-range electrostatics: |
75 |
|
|
integer, save :: summationMethod = NONE |
76 |
chrisfen |
603 |
logical, save :: summationMethodChecked = .false. |
77 |
gezelter |
602 |
real(kind=DP), save :: defaultCutoff = 0.0_DP |
78 |
|
|
logical, save :: haveDefaultCutoff = .false. |
79 |
|
|
real(kind=DP), save :: dampingAlpha = 0.0_DP |
80 |
|
|
logical, save :: haveDampingAlpha = .false. |
81 |
|
|
real(kind=DP), save :: dielectric = 0.0_DP |
82 |
|
|
logical, save :: haveDielectric = .false. |
83 |
|
|
real(kind=DP), save :: constERFC = 0.0_DP |
84 |
|
|
real(kind=DP), save :: constEXP = 0.0_DP |
85 |
|
|
logical, save :: haveDWAconstants = .false. |
86 |
chrisfen |
607 |
real(kind=dp), save :: rcuti = 0.0_dp |
87 |
|
|
real(kind=dp), save :: rcuti2 = 0.0_dp |
88 |
|
|
real(kind=dp), save :: rcuti3 = 0.0_dp |
89 |
|
|
real(kind=dp), save :: rcuti4 = 0.0_dp |
90 |
chrisfen |
626 |
logical, save :: is_Undamped_Wolf = .false. |
91 |
|
|
logical, save :: is_Damped_Wolf = .false. |
92 |
gezelter |
602 |
|
93 |
chuckv |
632 |
#ifdef __IFC |
94 |
|
|
! error function for ifc version > 7. |
95 |
chuckv |
631 |
double precision, external :: derfc |
96 |
chuckv |
632 |
#endif |
97 |
|
|
|
98 |
gezelter |
602 |
public :: setElectrostaticSummationMethod |
99 |
|
|
public :: setElectrostaticCutoffRadius |
100 |
|
|
public :: setDampedWolfAlpha |
101 |
|
|
public :: setReactionFieldDielectric |
102 |
gezelter |
411 |
public :: newElectrostaticType |
103 |
|
|
public :: setCharge |
104 |
|
|
public :: setDipoleMoment |
105 |
|
|
public :: setSplitDipoleDistance |
106 |
|
|
public :: setQuadrupoleMoments |
107 |
|
|
public :: doElectrostaticPair |
108 |
|
|
public :: getCharge |
109 |
|
|
public :: getDipoleMoment |
110 |
chrisfen |
445 |
public :: pre22 |
111 |
chuckv |
492 |
public :: destroyElectrostaticTypes |
112 |
gezelter |
411 |
|
113 |
|
|
type :: Electrostatic |
114 |
|
|
integer :: c_ident |
115 |
|
|
logical :: is_Charge = .false. |
116 |
|
|
logical :: is_Dipole = .false. |
117 |
|
|
logical :: is_SplitDipole = .false. |
118 |
|
|
logical :: is_Quadrupole = .false. |
119 |
chrisfen |
532 |
logical :: is_Tap = .false. |
120 |
gezelter |
411 |
real(kind=DP) :: charge = 0.0_DP |
121 |
|
|
real(kind=DP) :: dipole_moment = 0.0_DP |
122 |
|
|
real(kind=DP) :: split_dipole_distance = 0.0_DP |
123 |
|
|
real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP |
124 |
|
|
end type Electrostatic |
125 |
|
|
|
126 |
|
|
type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap |
127 |
|
|
|
128 |
|
|
contains |
129 |
|
|
|
130 |
gezelter |
602 |
subroutine setElectrostaticSummationMethod(the_ESM) |
131 |
|
|
integer, intent(in) :: the_ESM |
132 |
|
|
|
133 |
|
|
if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then |
134 |
|
|
call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method") |
135 |
|
|
endif |
136 |
|
|
|
137 |
chrisfen |
610 |
summationMethod = the_ESM |
138 |
chrisfen |
626 |
|
139 |
|
|
if (summationMethod .eq. UNDAMPED_WOLF) is_Undamped_Wolf = .true. |
140 |
|
|
if (summationMethod .eq. DAMPED_WOLF) is_Damped_Wolf = .true. |
141 |
gezelter |
602 |
end subroutine setElectrostaticSummationMethod |
142 |
|
|
|
143 |
|
|
subroutine setElectrostaticCutoffRadius(thisRcut) |
144 |
|
|
real(kind=dp), intent(in) :: thisRcut |
145 |
|
|
defaultCutoff = thisRcut |
146 |
|
|
haveDefaultCutoff = .true. |
147 |
|
|
end subroutine setElectrostaticCutoffRadius |
148 |
|
|
|
149 |
|
|
subroutine setDampedWolfAlpha(thisAlpha) |
150 |
|
|
real(kind=dp), intent(in) :: thisAlpha |
151 |
|
|
dampingAlpha = thisAlpha |
152 |
|
|
haveDampingAlpha = .true. |
153 |
|
|
end subroutine setDampedWolfAlpha |
154 |
|
|
|
155 |
|
|
subroutine setReactionFieldDielectric(thisDielectric) |
156 |
|
|
real(kind=dp), intent(in) :: thisDielectric |
157 |
|
|
dielectric = thisDielectric |
158 |
|
|
haveDielectric = .true. |
159 |
|
|
end subroutine setReactionFieldDielectric |
160 |
|
|
|
161 |
gezelter |
411 |
subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, & |
162 |
chrisfen |
532 |
is_SplitDipole, is_Quadrupole, is_Tap, status) |
163 |
gezelter |
507 |
|
164 |
gezelter |
411 |
integer, intent(in) :: c_ident |
165 |
|
|
logical, intent(in) :: is_Charge |
166 |
|
|
logical, intent(in) :: is_Dipole |
167 |
|
|
logical, intent(in) :: is_SplitDipole |
168 |
|
|
logical, intent(in) :: is_Quadrupole |
169 |
chrisfen |
532 |
logical, intent(in) :: is_Tap |
170 |
gezelter |
411 |
integer, intent(out) :: status |
171 |
|
|
integer :: nAtypes, myATID, i, j |
172 |
|
|
|
173 |
|
|
status = 0 |
174 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
175 |
gezelter |
507 |
|
176 |
gezelter |
411 |
!! Be simple-minded and assume that we need an ElectrostaticMap that |
177 |
|
|
!! is the same size as the total number of atom types |
178 |
|
|
|
179 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
180 |
gezelter |
507 |
|
181 |
gezelter |
411 |
nAtypes = getSize(atypes) |
182 |
gezelter |
507 |
|
183 |
gezelter |
411 |
if (nAtypes == 0) then |
184 |
|
|
status = -1 |
185 |
|
|
return |
186 |
|
|
end if |
187 |
gezelter |
507 |
|
188 |
gezelter |
411 |
if (.not. allocated(ElectrostaticMap)) then |
189 |
|
|
allocate(ElectrostaticMap(nAtypes)) |
190 |
|
|
endif |
191 |
gezelter |
507 |
|
192 |
gezelter |
411 |
end if |
193 |
|
|
|
194 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
195 |
|
|
status = -1 |
196 |
|
|
return |
197 |
|
|
endif |
198 |
gezelter |
507 |
|
199 |
gezelter |
411 |
! set the values for ElectrostaticMap for this atom type: |
200 |
|
|
|
201 |
|
|
ElectrostaticMap(myATID)%c_ident = c_ident |
202 |
|
|
ElectrostaticMap(myATID)%is_Charge = is_Charge |
203 |
|
|
ElectrostaticMap(myATID)%is_Dipole = is_Dipole |
204 |
|
|
ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole |
205 |
|
|
ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole |
206 |
chrisfen |
532 |
ElectrostaticMap(myATID)%is_Tap = is_Tap |
207 |
gezelter |
507 |
|
208 |
gezelter |
411 |
end subroutine newElectrostaticType |
209 |
|
|
|
210 |
|
|
subroutine setCharge(c_ident, charge, status) |
211 |
|
|
integer, intent(in) :: c_ident |
212 |
|
|
real(kind=dp), intent(in) :: charge |
213 |
|
|
integer, intent(out) :: status |
214 |
|
|
integer :: myATID |
215 |
|
|
|
216 |
|
|
status = 0 |
217 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
218 |
|
|
|
219 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
220 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!") |
221 |
|
|
status = -1 |
222 |
|
|
return |
223 |
|
|
end if |
224 |
|
|
|
225 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
226 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!") |
227 |
|
|
status = -1 |
228 |
|
|
return |
229 |
|
|
endif |
230 |
|
|
|
231 |
|
|
if (.not.ElectrostaticMap(myATID)%is_Charge) then |
232 |
|
|
call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!") |
233 |
|
|
status = -1 |
234 |
|
|
return |
235 |
gezelter |
507 |
endif |
236 |
gezelter |
411 |
|
237 |
|
|
ElectrostaticMap(myATID)%charge = charge |
238 |
|
|
end subroutine setCharge |
239 |
|
|
|
240 |
|
|
subroutine setDipoleMoment(c_ident, dipole_moment, status) |
241 |
|
|
integer, intent(in) :: c_ident |
242 |
|
|
real(kind=dp), intent(in) :: dipole_moment |
243 |
|
|
integer, intent(out) :: status |
244 |
|
|
integer :: myATID |
245 |
|
|
|
246 |
|
|
status = 0 |
247 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
248 |
|
|
|
249 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
250 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!") |
251 |
|
|
status = -1 |
252 |
|
|
return |
253 |
|
|
end if |
254 |
|
|
|
255 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
256 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!") |
257 |
|
|
status = -1 |
258 |
|
|
return |
259 |
|
|
endif |
260 |
|
|
|
261 |
|
|
if (.not.ElectrostaticMap(myATID)%is_Dipole) then |
262 |
|
|
call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!") |
263 |
|
|
status = -1 |
264 |
|
|
return |
265 |
|
|
endif |
266 |
|
|
|
267 |
|
|
ElectrostaticMap(myATID)%dipole_moment = dipole_moment |
268 |
|
|
end subroutine setDipoleMoment |
269 |
|
|
|
270 |
|
|
subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status) |
271 |
|
|
integer, intent(in) :: c_ident |
272 |
|
|
real(kind=dp), intent(in) :: split_dipole_distance |
273 |
|
|
integer, intent(out) :: status |
274 |
|
|
integer :: myATID |
275 |
|
|
|
276 |
|
|
status = 0 |
277 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
278 |
|
|
|
279 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
280 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!") |
281 |
|
|
status = -1 |
282 |
|
|
return |
283 |
|
|
end if |
284 |
|
|
|
285 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
286 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!") |
287 |
|
|
status = -1 |
288 |
|
|
return |
289 |
|
|
endif |
290 |
|
|
|
291 |
|
|
if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then |
292 |
|
|
call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!") |
293 |
|
|
status = -1 |
294 |
|
|
return |
295 |
|
|
endif |
296 |
|
|
|
297 |
|
|
ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance |
298 |
|
|
end subroutine setSplitDipoleDistance |
299 |
|
|
|
300 |
|
|
subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status) |
301 |
|
|
integer, intent(in) :: c_ident |
302 |
|
|
real(kind=dp), intent(in), dimension(3) :: quadrupole_moments |
303 |
|
|
integer, intent(out) :: status |
304 |
|
|
integer :: myATID, i, j |
305 |
|
|
|
306 |
|
|
status = 0 |
307 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
308 |
|
|
|
309 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
310 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!") |
311 |
|
|
status = -1 |
312 |
|
|
return |
313 |
|
|
end if |
314 |
|
|
|
315 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
316 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!") |
317 |
|
|
status = -1 |
318 |
|
|
return |
319 |
|
|
endif |
320 |
|
|
|
321 |
|
|
if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then |
322 |
|
|
call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!") |
323 |
|
|
status = -1 |
324 |
|
|
return |
325 |
|
|
endif |
326 |
gezelter |
507 |
|
327 |
gezelter |
411 |
do i = 1, 3 |
328 |
gezelter |
507 |
ElectrostaticMap(myATID)%quadrupole_moments(i) = & |
329 |
|
|
quadrupole_moments(i) |
330 |
|
|
enddo |
331 |
gezelter |
411 |
|
332 |
|
|
end subroutine setQuadrupoleMoments |
333 |
|
|
|
334 |
gezelter |
507 |
|
335 |
gezelter |
411 |
function getCharge(atid) result (c) |
336 |
|
|
integer, intent(in) :: atid |
337 |
|
|
integer :: localError |
338 |
|
|
real(kind=dp) :: c |
339 |
gezelter |
507 |
|
340 |
gezelter |
411 |
if (.not.allocated(ElectrostaticMap)) then |
341 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!") |
342 |
|
|
return |
343 |
|
|
end if |
344 |
gezelter |
507 |
|
345 |
gezelter |
411 |
if (.not.ElectrostaticMap(atid)%is_Charge) then |
346 |
|
|
call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!") |
347 |
|
|
return |
348 |
|
|
endif |
349 |
gezelter |
507 |
|
350 |
gezelter |
411 |
c = ElectrostaticMap(atid)%charge |
351 |
|
|
end function getCharge |
352 |
|
|
|
353 |
|
|
function getDipoleMoment(atid) result (dm) |
354 |
|
|
integer, intent(in) :: atid |
355 |
|
|
integer :: localError |
356 |
|
|
real(kind=dp) :: dm |
357 |
gezelter |
507 |
|
358 |
gezelter |
411 |
if (.not.allocated(ElectrostaticMap)) then |
359 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!") |
360 |
|
|
return |
361 |
|
|
end if |
362 |
gezelter |
507 |
|
363 |
gezelter |
411 |
if (.not.ElectrostaticMap(atid)%is_Dipole) then |
364 |
|
|
call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!") |
365 |
|
|
return |
366 |
|
|
endif |
367 |
gezelter |
507 |
|
368 |
gezelter |
411 |
dm = ElectrostaticMap(atid)%dipole_moment |
369 |
|
|
end function getDipoleMoment |
370 |
|
|
|
371 |
gezelter |
602 |
subroutine checkSummationMethod() |
372 |
|
|
|
373 |
chrisfen |
607 |
if (.not.haveDefaultCutoff) then |
374 |
|
|
call handleError("checkSummationMethod", "no Default Cutoff set!") |
375 |
|
|
endif |
376 |
|
|
|
377 |
|
|
rcuti = 1.0d0 / defaultCutoff |
378 |
|
|
rcuti2 = rcuti*rcuti |
379 |
|
|
rcuti3 = rcuti2*rcuti |
380 |
|
|
rcuti4 = rcuti2*rcuti2 |
381 |
|
|
|
382 |
gezelter |
602 |
if (summationMethod .eq. DAMPED_WOLF) then |
383 |
|
|
if (.not.haveDWAconstants) then |
384 |
|
|
|
385 |
|
|
if (.not.haveDampingAlpha) then |
386 |
|
|
call handleError("checkSummationMethod", "no Damping Alpha set!") |
387 |
|
|
endif |
388 |
|
|
|
389 |
chrisfen |
603 |
if (.not.haveDefaultCutoff) then |
390 |
|
|
call handleError("checkSummationMethod", "no Default Cutoff set!") |
391 |
|
|
endif |
392 |
|
|
|
393 |
|
|
constEXP = exp(-dampingAlpha*dampingAlpha*defaultCutoff*defaultCutoff) |
394 |
chuckv |
631 |
constERFC = derfc(dampingAlpha*defaultCutoff) |
395 |
gezelter |
602 |
|
396 |
|
|
haveDWAconstants = .true. |
397 |
|
|
endif |
398 |
|
|
endif |
399 |
|
|
|
400 |
chrisfen |
603 |
if (summationMethod .eq. REACTION_FIELD) then |
401 |
|
|
if (.not.haveDielectric) then |
402 |
|
|
call handleError("checkSummationMethod", "no reaction field Dielectric set!") |
403 |
|
|
endif |
404 |
|
|
endif |
405 |
|
|
|
406 |
|
|
summationMethodChecked = .true. |
407 |
gezelter |
602 |
end subroutine checkSummationMethod |
408 |
|
|
|
409 |
|
|
|
410 |
|
|
|
411 |
gezelter |
411 |
subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, & |
412 |
chrisfen |
607 |
vpair, fpair, pot, eFrame, f, t, do_pot) |
413 |
gezelter |
507 |
|
414 |
gezelter |
411 |
logical, intent(in) :: do_pot |
415 |
gezelter |
507 |
|
416 |
gezelter |
411 |
integer, intent(in) :: atom1, atom2 |
417 |
|
|
integer :: localError |
418 |
|
|
|
419 |
chrisfen |
607 |
real(kind=dp), intent(in) :: rij, r2, sw |
420 |
gezelter |
411 |
real(kind=dp), intent(in), dimension(3) :: d |
421 |
|
|
real(kind=dp), intent(inout) :: vpair |
422 |
|
|
real(kind=dp), intent(inout), dimension(3) :: fpair |
423 |
|
|
|
424 |
chrisfen |
626 |
real( kind = dp ) :: pot |
425 |
gezelter |
411 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
426 |
|
|
real( kind = dp ), dimension(3,nLocal) :: f |
427 |
|
|
real( kind = dp ), dimension(3,nLocal) :: t |
428 |
gezelter |
507 |
|
429 |
gezelter |
439 |
real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i |
430 |
|
|
real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j |
431 |
|
|
real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i |
432 |
|
|
real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j |
433 |
gezelter |
411 |
|
434 |
|
|
logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole |
435 |
|
|
logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole |
436 |
chrisfen |
532 |
logical :: i_is_Tap, j_is_Tap |
437 |
gezelter |
411 |
integer :: me1, me2, id1, id2 |
438 |
|
|
real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j |
439 |
gezelter |
439 |
real (kind=dp) :: qxx_i, qyy_i, qzz_i |
440 |
|
|
real (kind=dp) :: qxx_j, qyy_j, qzz_j |
441 |
|
|
real (kind=dp) :: cx_i, cy_i, cz_i |
442 |
|
|
real (kind=dp) :: cx_j, cy_j, cz_j |
443 |
|
|
real (kind=dp) :: cx2, cy2, cz2 |
444 |
gezelter |
411 |
real (kind=dp) :: ct_i, ct_j, ct_ij, a1 |
445 |
gezelter |
421 |
real (kind=dp) :: riji, ri, ri2, ri3, ri4 |
446 |
chrisfen |
597 |
real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2 |
447 |
gezelter |
421 |
real (kind=dp) :: xhat, yhat, zhat |
448 |
gezelter |
411 |
real (kind=dp) :: dudx, dudy, dudz |
449 |
chrisfen |
626 |
real (kind=dp) :: scale, sc2, bigR |
450 |
gezelter |
411 |
|
451 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
452 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!") |
453 |
|
|
return |
454 |
|
|
end if |
455 |
|
|
|
456 |
gezelter |
602 |
if (.not.summationMethodChecked) then |
457 |
|
|
call checkSummationMethod() |
458 |
chrisfen |
611 |
|
459 |
gezelter |
602 |
endif |
460 |
|
|
|
461 |
|
|
|
462 |
gezelter |
411 |
#ifdef IS_MPI |
463 |
|
|
me1 = atid_Row(atom1) |
464 |
|
|
me2 = atid_Col(atom2) |
465 |
|
|
#else |
466 |
|
|
me1 = atid(atom1) |
467 |
|
|
me2 = atid(atom2) |
468 |
|
|
#endif |
469 |
|
|
|
470 |
|
|
!! some variables we'll need independent of electrostatic type: |
471 |
|
|
|
472 |
|
|
riji = 1.0d0 / rij |
473 |
|
|
|
474 |
gezelter |
421 |
xhat = d(1) * riji |
475 |
|
|
yhat = d(2) * riji |
476 |
|
|
zhat = d(3) * riji |
477 |
gezelter |
411 |
|
478 |
|
|
!! logicals |
479 |
|
|
i_is_Charge = ElectrostaticMap(me1)%is_Charge |
480 |
|
|
i_is_Dipole = ElectrostaticMap(me1)%is_Dipole |
481 |
|
|
i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole |
482 |
|
|
i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole |
483 |
chrisfen |
532 |
i_is_Tap = ElectrostaticMap(me1)%is_Tap |
484 |
gezelter |
411 |
|
485 |
|
|
j_is_Charge = ElectrostaticMap(me2)%is_Charge |
486 |
|
|
j_is_Dipole = ElectrostaticMap(me2)%is_Dipole |
487 |
|
|
j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole |
488 |
|
|
j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole |
489 |
chrisfen |
532 |
j_is_Tap = ElectrostaticMap(me2)%is_Tap |
490 |
gezelter |
411 |
|
491 |
|
|
if (i_is_Charge) then |
492 |
|
|
q_i = ElectrostaticMap(me1)%charge |
493 |
|
|
endif |
494 |
gezelter |
507 |
|
495 |
gezelter |
411 |
if (i_is_Dipole) then |
496 |
|
|
mu_i = ElectrostaticMap(me1)%dipole_moment |
497 |
|
|
#ifdef IS_MPI |
498 |
gezelter |
439 |
uz_i(1) = eFrame_Row(3,atom1) |
499 |
|
|
uz_i(2) = eFrame_Row(6,atom1) |
500 |
|
|
uz_i(3) = eFrame_Row(9,atom1) |
501 |
gezelter |
411 |
#else |
502 |
gezelter |
439 |
uz_i(1) = eFrame(3,atom1) |
503 |
|
|
uz_i(2) = eFrame(6,atom1) |
504 |
|
|
uz_i(3) = eFrame(9,atom1) |
505 |
gezelter |
411 |
#endif |
506 |
gezelter |
439 |
ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat |
507 |
gezelter |
411 |
|
508 |
|
|
if (i_is_SplitDipole) then |
509 |
|
|
d_i = ElectrostaticMap(me1)%split_dipole_distance |
510 |
|
|
endif |
511 |
gezelter |
507 |
|
512 |
gezelter |
411 |
endif |
513 |
|
|
|
514 |
gezelter |
439 |
if (i_is_Quadrupole) then |
515 |
|
|
qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1) |
516 |
|
|
qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2) |
517 |
|
|
qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3) |
518 |
|
|
#ifdef IS_MPI |
519 |
|
|
ux_i(1) = eFrame_Row(1,atom1) |
520 |
|
|
ux_i(2) = eFrame_Row(4,atom1) |
521 |
|
|
ux_i(3) = eFrame_Row(7,atom1) |
522 |
|
|
uy_i(1) = eFrame_Row(2,atom1) |
523 |
|
|
uy_i(2) = eFrame_Row(5,atom1) |
524 |
|
|
uy_i(3) = eFrame_Row(8,atom1) |
525 |
|
|
uz_i(1) = eFrame_Row(3,atom1) |
526 |
|
|
uz_i(2) = eFrame_Row(6,atom1) |
527 |
|
|
uz_i(3) = eFrame_Row(9,atom1) |
528 |
|
|
#else |
529 |
|
|
ux_i(1) = eFrame(1,atom1) |
530 |
|
|
ux_i(2) = eFrame(4,atom1) |
531 |
|
|
ux_i(3) = eFrame(7,atom1) |
532 |
|
|
uy_i(1) = eFrame(2,atom1) |
533 |
|
|
uy_i(2) = eFrame(5,atom1) |
534 |
|
|
uy_i(3) = eFrame(8,atom1) |
535 |
|
|
uz_i(1) = eFrame(3,atom1) |
536 |
|
|
uz_i(2) = eFrame(6,atom1) |
537 |
|
|
uz_i(3) = eFrame(9,atom1) |
538 |
|
|
#endif |
539 |
|
|
cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat |
540 |
|
|
cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat |
541 |
|
|
cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat |
542 |
|
|
endif |
543 |
|
|
|
544 |
gezelter |
411 |
if (j_is_Charge) then |
545 |
|
|
q_j = ElectrostaticMap(me2)%charge |
546 |
|
|
endif |
547 |
gezelter |
507 |
|
548 |
gezelter |
411 |
if (j_is_Dipole) then |
549 |
|
|
mu_j = ElectrostaticMap(me2)%dipole_moment |
550 |
|
|
#ifdef IS_MPI |
551 |
gezelter |
439 |
uz_j(1) = eFrame_Col(3,atom2) |
552 |
|
|
uz_j(2) = eFrame_Col(6,atom2) |
553 |
|
|
uz_j(3) = eFrame_Col(9,atom2) |
554 |
gezelter |
411 |
#else |
555 |
gezelter |
439 |
uz_j(1) = eFrame(3,atom2) |
556 |
|
|
uz_j(2) = eFrame(6,atom2) |
557 |
|
|
uz_j(3) = eFrame(9,atom2) |
558 |
gezelter |
411 |
#endif |
559 |
chrisfen |
465 |
ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat |
560 |
gezelter |
411 |
|
561 |
|
|
if (j_is_SplitDipole) then |
562 |
|
|
d_j = ElectrostaticMap(me2)%split_dipole_distance |
563 |
|
|
endif |
564 |
|
|
endif |
565 |
|
|
|
566 |
gezelter |
439 |
if (j_is_Quadrupole) then |
567 |
|
|
qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1) |
568 |
|
|
qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2) |
569 |
|
|
qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3) |
570 |
|
|
#ifdef IS_MPI |
571 |
|
|
ux_j(1) = eFrame_Col(1,atom2) |
572 |
|
|
ux_j(2) = eFrame_Col(4,atom2) |
573 |
|
|
ux_j(3) = eFrame_Col(7,atom2) |
574 |
|
|
uy_j(1) = eFrame_Col(2,atom2) |
575 |
|
|
uy_j(2) = eFrame_Col(5,atom2) |
576 |
|
|
uy_j(3) = eFrame_Col(8,atom2) |
577 |
|
|
uz_j(1) = eFrame_Col(3,atom2) |
578 |
|
|
uz_j(2) = eFrame_Col(6,atom2) |
579 |
|
|
uz_j(3) = eFrame_Col(9,atom2) |
580 |
|
|
#else |
581 |
|
|
ux_j(1) = eFrame(1,atom2) |
582 |
|
|
ux_j(2) = eFrame(4,atom2) |
583 |
|
|
ux_j(3) = eFrame(7,atom2) |
584 |
|
|
uy_j(1) = eFrame(2,atom2) |
585 |
|
|
uy_j(2) = eFrame(5,atom2) |
586 |
|
|
uy_j(3) = eFrame(8,atom2) |
587 |
|
|
uz_j(1) = eFrame(3,atom2) |
588 |
|
|
uz_j(2) = eFrame(6,atom2) |
589 |
|
|
uz_j(3) = eFrame(9,atom2) |
590 |
|
|
#endif |
591 |
|
|
cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat |
592 |
|
|
cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat |
593 |
|
|
cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat |
594 |
|
|
endif |
595 |
chrisfen |
554 |
|
596 |
gezelter |
411 |
epot = 0.0_dp |
597 |
|
|
dudx = 0.0_dp |
598 |
|
|
dudy = 0.0_dp |
599 |
|
|
dudz = 0.0_dp |
600 |
|
|
|
601 |
gezelter |
439 |
dudux_i = 0.0_dp |
602 |
|
|
duduy_i = 0.0_dp |
603 |
|
|
duduz_i = 0.0_dp |
604 |
gezelter |
411 |
|
605 |
gezelter |
439 |
dudux_j = 0.0_dp |
606 |
|
|
duduy_j = 0.0_dp |
607 |
|
|
duduz_j = 0.0_dp |
608 |
gezelter |
411 |
|
609 |
|
|
if (i_is_Charge) then |
610 |
|
|
|
611 |
|
|
if (j_is_Charge) then |
612 |
gezelter |
507 |
|
613 |
chrisfen |
611 |
if (summationMethod .eq. UNDAMPED_WOLF) then |
614 |
chrisfen |
626 |
|
615 |
chrisfen |
597 |
vterm = pre11 * q_i * q_j * (riji - rcuti) |
616 |
|
|
vpair = vpair + vterm |
617 |
chrisfen |
626 |
epot = epot + sw*vterm |
618 |
chrisfen |
597 |
|
619 |
|
|
dudr = - sw * pre11 * q_i * q_j * (riji*riji*riji - rcuti2*rcuti) |
620 |
|
|
|
621 |
|
|
dudx = dudx + dudr * d(1) |
622 |
|
|
dudy = dudy + dudr * d(2) |
623 |
|
|
dudz = dudz + dudr * d(3) |
624 |
gezelter |
411 |
|
625 |
chrisfen |
597 |
else |
626 |
chrisfen |
626 |
|
627 |
chrisfen |
597 |
vterm = pre11 * q_i * q_j * riji |
628 |
|
|
vpair = vpair + vterm |
629 |
chrisfen |
626 |
epot = epot + sw*vterm |
630 |
chrisfen |
597 |
|
631 |
|
|
dudr = - sw * vterm * riji |
632 |
|
|
|
633 |
|
|
dudx = dudx + dudr * xhat |
634 |
|
|
dudy = dudy + dudr * yhat |
635 |
|
|
dudz = dudz + dudr * zhat |
636 |
|
|
|
637 |
|
|
endif |
638 |
|
|
|
639 |
gezelter |
411 |
endif |
640 |
|
|
|
641 |
|
|
if (j_is_Dipole) then |
642 |
|
|
|
643 |
chrisfen |
626 |
pref = pre12 * q_i * mu_j |
644 |
gezelter |
411 |
|
645 |
chrisfen |
611 |
if (summationMethod .eq. UNDAMPED_WOLF) then |
646 |
chrisfen |
597 |
ri2 = riji * riji |
647 |
|
|
ri3 = ri2 * riji |
648 |
gezelter |
507 |
|
649 |
chrisfen |
626 |
pref = pre12 * q_i * mu_j |
650 |
chrisfen |
597 |
vterm = - pref * ct_j * (ri2 - rcuti2) |
651 |
chrisfen |
626 |
vpair = vpair + vterm |
652 |
|
|
epot = epot + sw*vterm |
653 |
chrisfen |
597 |
|
654 |
|
|
!! this has a + sign in the () because the rij vector is |
655 |
|
|
!! r_j - r_i and the charge-dipole potential takes the origin |
656 |
|
|
!! as the point dipole, which is atom j in this case. |
657 |
|
|
|
658 |
chrisfen |
626 |
dudx = dudx - sw*pref * ( ri3*( uz_j(1) - 3.0d0*ct_j*xhat) & |
659 |
chrisfen |
597 |
- rcuti3*( uz_j(1) - 3.0d0*ct_j*d(1)*rcuti ) ) |
660 |
chrisfen |
626 |
dudy = dudy - sw*pref * ( ri3*( uz_j(2) - 3.0d0*ct_j*yhat) & |
661 |
chrisfen |
597 |
- rcuti3*( uz_j(2) - 3.0d0*ct_j*d(2)*rcuti ) ) |
662 |
chrisfen |
626 |
dudz = dudz - sw*pref * ( ri3*( uz_j(3) - 3.0d0*ct_j*zhat) & |
663 |
chrisfen |
597 |
- rcuti3*( uz_j(3) - 3.0d0*ct_j*d(3)*rcuti ) ) |
664 |
|
|
|
665 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) - sw*pref*( ri2*xhat - d(1)*rcuti3 ) |
666 |
|
|
duduz_j(2) = duduz_j(2) - sw*pref*( ri2*yhat - d(2)*rcuti3 ) |
667 |
|
|
duduz_j(3) = duduz_j(3) - sw*pref*( ri2*zhat - d(3)*rcuti3 ) |
668 |
gezelter |
411 |
|
669 |
chrisfen |
597 |
else |
670 |
|
|
if (j_is_SplitDipole) then |
671 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_j * d_j) |
672 |
|
|
ri = 1.0_dp / BigR |
673 |
|
|
scale = rij * ri |
674 |
|
|
else |
675 |
|
|
ri = riji |
676 |
|
|
scale = 1.0_dp |
677 |
|
|
endif |
678 |
|
|
|
679 |
|
|
ri2 = ri * ri |
680 |
|
|
ri3 = ri2 * ri |
681 |
|
|
sc2 = scale * scale |
682 |
chrisfen |
626 |
|
683 |
|
|
pref = pre12 * q_i * mu_j |
684 |
chrisfen |
597 |
vterm = - pref * ct_j * ri2 * scale |
685 |
chrisfen |
626 |
vpair = vpair + vterm |
686 |
|
|
epot = epot + sw*vterm |
687 |
chrisfen |
597 |
|
688 |
|
|
!! this has a + sign in the () because the rij vector is |
689 |
|
|
!! r_j - r_i and the charge-dipole potential takes the origin |
690 |
|
|
!! as the point dipole, which is atom j in this case. |
691 |
|
|
|
692 |
chrisfen |
626 |
dudx = dudx - sw*pref * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2) |
693 |
|
|
dudy = dudy - sw*pref * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2) |
694 |
|
|
dudz = dudz - sw*pref * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2) |
695 |
chrisfen |
597 |
|
696 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) - sw*pref * ri2 * xhat * scale |
697 |
|
|
duduz_j(2) = duduz_j(2) - sw*pref * ri2 * yhat * scale |
698 |
|
|
duduz_j(3) = duduz_j(3) - sw*pref * ri2 * zhat * scale |
699 |
gezelter |
411 |
|
700 |
chrisfen |
597 |
endif |
701 |
gezelter |
411 |
endif |
702 |
gezelter |
421 |
|
703 |
gezelter |
439 |
if (j_is_Quadrupole) then |
704 |
|
|
ri2 = riji * riji |
705 |
|
|
ri3 = ri2 * riji |
706 |
gezelter |
440 |
ri4 = ri2 * ri2 |
707 |
gezelter |
439 |
cx2 = cx_j * cx_j |
708 |
|
|
cy2 = cy_j * cy_j |
709 |
|
|
cz2 = cz_j * cz_j |
710 |
|
|
|
711 |
chrisfen |
611 |
if (summationMethod .eq. UNDAMPED_WOLF) then |
712 |
chrisfen |
626 |
pref = pre14 * q_i / 3.0_dp |
713 |
chrisfen |
597 |
vterm1 = pref * ri3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + & |
714 |
|
|
qyy_j * (3.0_dp*cy2 - 1.0_dp) + & |
715 |
|
|
qzz_j * (3.0_dp*cz2 - 1.0_dp) ) |
716 |
|
|
vterm2 = pref * rcuti3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + & |
717 |
|
|
qyy_j * (3.0_dp*cy2 - 1.0_dp) + & |
718 |
|
|
qzz_j * (3.0_dp*cz2 - 1.0_dp) ) |
719 |
chrisfen |
626 |
vpair = vpair + ( vterm1 - vterm2 ) |
720 |
|
|
epot = epot + sw*( vterm1 - vterm2 ) |
721 |
chrisfen |
597 |
|
722 |
|
|
dudx = dudx - (5.0_dp * & |
723 |
chrisfen |
626 |
(vterm1*riji*xhat - vterm2*rcuti2*d(1))) + sw*pref * ( & |
724 |
chrisfen |
597 |
(ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(1)) - & |
725 |
|
|
qxx_j*2.0_dp*(xhat - rcuti*d(1))) + & |
726 |
|
|
(ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(1)) - & |
727 |
|
|
qyy_j*2.0_dp*(xhat - rcuti*d(1))) + & |
728 |
|
|
(ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(1)) - & |
729 |
|
|
qzz_j*2.0_dp*(xhat - rcuti*d(1))) ) |
730 |
|
|
dudy = dudy - (5.0_dp * & |
731 |
chrisfen |
626 |
(vterm1*riji*yhat - vterm2*rcuti2*d(2))) + sw*pref * ( & |
732 |
chrisfen |
597 |
(ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(2)) - & |
733 |
|
|
qxx_j*2.0_dp*(yhat - rcuti*d(2))) + & |
734 |
|
|
(ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(2)) - & |
735 |
|
|
qyy_j*2.0_dp*(yhat - rcuti*d(2))) + & |
736 |
|
|
(ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(2)) - & |
737 |
|
|
qzz_j*2.0_dp*(yhat - rcuti*d(2))) ) |
738 |
|
|
dudz = dudz - (5.0_dp * & |
739 |
chrisfen |
626 |
(vterm1*riji*zhat - vterm2*rcuti2*d(3))) + sw*pref * ( & |
740 |
chrisfen |
597 |
(ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(3)) - & |
741 |
|
|
qxx_j*2.0_dp*(zhat - rcuti*d(3))) + & |
742 |
|
|
(ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(3)) - & |
743 |
|
|
qyy_j*2.0_dp*(zhat - rcuti*d(3))) + & |
744 |
|
|
(ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(3)) - & |
745 |
|
|
qzz_j*2.0_dp*(zhat - rcuti*d(3))) ) |
746 |
|
|
|
747 |
chrisfen |
626 |
dudux_j(1) = dudux_j(1) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*xhat) -& |
748 |
chrisfen |
597 |
rcuti4*(qxx_j*6.0_dp*cx_j*d(1))) |
749 |
chrisfen |
626 |
dudux_j(2) = dudux_j(2) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*yhat) -& |
750 |
chrisfen |
597 |
rcuti4*(qxx_j*6.0_dp*cx_j*d(2))) |
751 |
chrisfen |
626 |
dudux_j(3) = dudux_j(3) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*zhat) -& |
752 |
chrisfen |
597 |
rcuti4*(qxx_j*6.0_dp*cx_j*d(3))) |
753 |
|
|
|
754 |
chrisfen |
626 |
duduy_j(1) = duduy_j(1) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*xhat) -& |
755 |
chrisfen |
597 |
rcuti4*(qyy_j*6.0_dp*cx_j*d(1))) |
756 |
chrisfen |
626 |
duduy_j(2) = duduy_j(2) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*yhat) -& |
757 |
chrisfen |
597 |
rcuti4*(qyy_j*6.0_dp*cx_j*d(2))) |
758 |
chrisfen |
626 |
duduy_j(3) = duduy_j(3) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*zhat) -& |
759 |
chrisfen |
597 |
rcuti4*(qyy_j*6.0_dp*cx_j*d(3))) |
760 |
|
|
|
761 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*xhat) -& |
762 |
chrisfen |
597 |
rcuti4*(qzz_j*6.0_dp*cx_j*d(1))) |
763 |
chrisfen |
626 |
duduz_j(2) = duduz_j(2) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*yhat) -& |
764 |
chrisfen |
597 |
rcuti4*(qzz_j*6.0_dp*cx_j*d(2))) |
765 |
chrisfen |
626 |
duduz_j(3) = duduz_j(3) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*zhat) -& |
766 |
chrisfen |
597 |
rcuti4*(qzz_j*6.0_dp*cx_j*d(3))) |
767 |
|
|
|
768 |
|
|
else |
769 |
chrisfen |
626 |
pref = pre14 * q_i / 3.0_dp |
770 |
chrisfen |
597 |
vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + & |
771 |
|
|
qyy_j * (3.0_dp*cy2 - 1.0_dp) + & |
772 |
|
|
qzz_j * (3.0_dp*cz2 - 1.0_dp)) |
773 |
chrisfen |
626 |
vpair = vpair + vterm |
774 |
|
|
epot = epot + sw*vterm |
775 |
chrisfen |
597 |
|
776 |
chrisfen |
626 |
dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + sw*pref * ri4 * ( & |
777 |
chrisfen |
597 |
qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + & |
778 |
|
|
qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + & |
779 |
|
|
qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) ) |
780 |
chrisfen |
626 |
dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + sw*pref * ri4 * ( & |
781 |
chrisfen |
597 |
qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + & |
782 |
|
|
qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + & |
783 |
|
|
qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) ) |
784 |
chrisfen |
626 |
dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + sw*pref * ri4 * ( & |
785 |
chrisfen |
597 |
qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + & |
786 |
|
|
qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + & |
787 |
|
|
qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) ) |
788 |
|
|
|
789 |
chrisfen |
626 |
dudux_j(1) = dudux_j(1) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*xhat) |
790 |
|
|
dudux_j(2) = dudux_j(2) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*yhat) |
791 |
|
|
dudux_j(3) = dudux_j(3) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*zhat) |
792 |
chrisfen |
597 |
|
793 |
chrisfen |
626 |
duduy_j(1) = duduy_j(1) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*xhat) |
794 |
|
|
duduy_j(2) = duduy_j(2) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*yhat) |
795 |
|
|
duduy_j(3) = duduy_j(3) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*zhat) |
796 |
chrisfen |
597 |
|
797 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*xhat) |
798 |
|
|
duduz_j(2) = duduz_j(2) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*yhat) |
799 |
|
|
duduz_j(3) = duduz_j(3) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*zhat) |
800 |
chrisfen |
597 |
|
801 |
|
|
endif |
802 |
gezelter |
439 |
endif |
803 |
gezelter |
411 |
endif |
804 |
gezelter |
507 |
|
805 |
gezelter |
411 |
if (i_is_Dipole) then |
806 |
gezelter |
507 |
|
807 |
gezelter |
411 |
if (j_is_Charge) then |
808 |
chrisfen |
626 |
|
809 |
|
|
pref = pre12 * q_j * mu_i |
810 |
|
|
|
811 |
chrisfen |
611 |
if (summationMethod .eq. UNDAMPED_WOLF) then |
812 |
chrisfen |
597 |
ri2 = riji * riji |
813 |
|
|
ri3 = ri2 * riji |
814 |
gezelter |
507 |
|
815 |
chrisfen |
626 |
pref = pre12 * q_j * mu_i |
816 |
chrisfen |
597 |
vterm = pref * ct_i * (ri2 - rcuti2) |
817 |
chrisfen |
626 |
vpair = vpair + vterm |
818 |
|
|
epot = epot + sw*vterm |
819 |
chrisfen |
597 |
|
820 |
|
|
!! this has a + sign in the () because the rij vector is |
821 |
|
|
!! r_j - r_i and the charge-dipole potential takes the origin |
822 |
|
|
!! as the point dipole, which is atom j in this case. |
823 |
|
|
|
824 |
chrisfen |
626 |
dudx = dudx + sw*pref * ( ri3*( uz_i(1) - 3.0d0*ct_i*xhat) & |
825 |
chrisfen |
597 |
- rcuti3*( uz_i(1) - 3.0d0*ct_i*d(1)*rcuti ) ) |
826 |
chrisfen |
626 |
dudy = dudy + sw*pref * ( ri3*( uz_i(2) - 3.0d0*ct_i*yhat) & |
827 |
chrisfen |
597 |
- rcuti3*( uz_i(2) - 3.0d0*ct_i*d(2)*rcuti ) ) |
828 |
chrisfen |
626 |
dudz = dudz + sw*pref * ( ri3*( uz_i(3) - 3.0d0*ct_i*zhat) & |
829 |
chrisfen |
597 |
- rcuti3*( uz_i(3) - 3.0d0*ct_i*d(3)*rcuti ) ) |
830 |
|
|
|
831 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) - sw*pref*( ri2*xhat - d(1)*rcuti3 ) |
832 |
|
|
duduz_i(2) = duduz_i(2) - sw*pref*( ri2*yhat - d(2)*rcuti3 ) |
833 |
|
|
duduz_i(3) = duduz_i(3) - sw*pref*( ri2*zhat - d(3)*rcuti3 ) |
834 |
gezelter |
411 |
|
835 |
chrisfen |
597 |
else |
836 |
|
|
if (i_is_SplitDipole) then |
837 |
gezelter |
421 |
BigR = sqrt(r2 + 0.25_dp * d_i * d_i) |
838 |
|
|
ri = 1.0_dp / BigR |
839 |
chrisfen |
597 |
scale = rij * ri |
840 |
|
|
else |
841 |
gezelter |
421 |
ri = riji |
842 |
|
|
scale = 1.0_dp |
843 |
|
|
endif |
844 |
chrisfen |
597 |
|
845 |
|
|
ri2 = ri * ri |
846 |
|
|
ri3 = ri2 * ri |
847 |
|
|
sc2 = scale * scale |
848 |
chrisfen |
626 |
|
849 |
|
|
pref = pre12 * q_j * mu_i |
850 |
chrisfen |
597 |
vterm = pref * ct_i * ri2 * scale |
851 |
chrisfen |
626 |
vpair = vpair + vterm |
852 |
|
|
epot = epot + sw*vterm |
853 |
chrisfen |
597 |
|
854 |
chrisfen |
626 |
dudx = dudx + sw*pref * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2) |
855 |
|
|
dudy = dudy + sw*pref * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2) |
856 |
|
|
dudz = dudz + sw*pref * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2) |
857 |
chrisfen |
597 |
|
858 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref * ri2 * xhat * scale |
859 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref * ri2 * yhat * scale |
860 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref * ri2 * zhat * scale |
861 |
gezelter |
421 |
endif |
862 |
chrisfen |
597 |
endif |
863 |
chrisfen |
626 |
|
864 |
chrisfen |
597 |
if (j_is_Dipole) then |
865 |
gezelter |
421 |
|
866 |
chrisfen |
611 |
if (summationMethod .eq. UNDAMPED_WOLF) then |
867 |
chrisfen |
597 |
ri2 = riji * riji |
868 |
|
|
ri3 = ri2 * riji |
869 |
|
|
ri4 = ri2 * ri2 |
870 |
gezelter |
507 |
|
871 |
chrisfen |
626 |
pref = pre22 * mu_i * mu_j |
872 |
chrisfen |
597 |
vterm = pref * (ri3 - rcuti3) * (ct_ij - 3.0d0 * ct_i * ct_j) |
873 |
chrisfen |
626 |
vpair = vpair + vterm |
874 |
|
|
epot = epot + sw*vterm |
875 |
chrisfen |
597 |
|
876 |
|
|
a1 = 5.0d0 * ct_i * ct_j - ct_ij |
877 |
|
|
|
878 |
chrisfen |
626 |
dudx = dudx + sw*pref*3.0d0*ri4 & |
879 |
|
|
* (a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) & |
880 |
|
|
- sw*pref*3.0d0*rcuti4 & |
881 |
|
|
* (a1*rcuti*d(1)-ct_i*uz_j(1)-ct_j*uz_i(1)) |
882 |
|
|
dudy = dudy + sw*pref*3.0d0*ri4 & |
883 |
|
|
* (a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) & |
884 |
|
|
- sw*pref*3.0d0*rcuti4 & |
885 |
|
|
* (a1*rcuti*d(2)-ct_i*uz_j(2)-ct_j*uz_i(2)) |
886 |
|
|
dudz = dudz + sw*pref*3.0d0*ri4 & |
887 |
|
|
* (a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) & |
888 |
|
|
- sw*pref*3.0d0*rcuti4 & |
889 |
|
|
* (a1*rcuti*d(3)-ct_i*uz_j(3)-ct_j*uz_i(3)) |
890 |
chrisfen |
597 |
|
891 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref*(ri3*(uz_j(1)-3.0d0*ct_j*xhat) & |
892 |
chrisfen |
597 |
- rcuti3*(uz_j(1) - 3.0d0*ct_j*d(1)*rcuti)) |
893 |
chrisfen |
626 |
duduz_i(2) = duduz_i(2) + sw*pref*(ri3*(uz_j(2)-3.0d0*ct_j*yhat) & |
894 |
chrisfen |
597 |
- rcuti3*(uz_j(2) - 3.0d0*ct_j*d(2)*rcuti)) |
895 |
chrisfen |
626 |
duduz_i(3) = duduz_i(3) + sw*pref*(ri3*(uz_j(3)-3.0d0*ct_j*zhat) & |
896 |
chrisfen |
597 |
- rcuti3*(uz_j(3) - 3.0d0*ct_j*d(3)*rcuti)) |
897 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) + sw*pref*(ri3*(uz_i(1)-3.0d0*ct_i*xhat) & |
898 |
chrisfen |
597 |
- rcuti3*(uz_i(1) - 3.0d0*ct_i*d(1)*rcuti)) |
899 |
chrisfen |
626 |
duduz_j(2) = duduz_j(2) + sw*pref*(ri3*(uz_i(2)-3.0d0*ct_i*yhat) & |
900 |
chrisfen |
597 |
- rcuti3*(uz_i(2) - 3.0d0*ct_i*d(2)*rcuti)) |
901 |
chrisfen |
626 |
duduz_j(3) = duduz_j(3) + sw*pref*(ri3*(uz_i(3)-3.0d0*ct_i*zhat) & |
902 |
chrisfen |
597 |
- rcuti3*(uz_i(3) - 3.0d0*ct_i*d(3)*rcuti)) |
903 |
chrisfen |
626 |
|
904 |
chrisfen |
597 |
else |
905 |
|
|
if (i_is_SplitDipole) then |
906 |
|
|
if (j_is_SplitDipole) then |
907 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j) |
908 |
|
|
else |
909 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_i * d_i) |
910 |
|
|
endif |
911 |
|
|
ri = 1.0_dp / BigR |
912 |
|
|
scale = rij * ri |
913 |
|
|
else |
914 |
|
|
if (j_is_SplitDipole) then |
915 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_j * d_j) |
916 |
|
|
ri = 1.0_dp / BigR |
917 |
|
|
scale = rij * ri |
918 |
|
|
else |
919 |
|
|
ri = riji |
920 |
|
|
scale = 1.0_dp |
921 |
|
|
endif |
922 |
|
|
endif |
923 |
|
|
|
924 |
|
|
ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3) |
925 |
|
|
|
926 |
|
|
ri2 = ri * ri |
927 |
|
|
ri3 = ri2 * ri |
928 |
|
|
ri4 = ri2 * ri2 |
929 |
|
|
sc2 = scale * scale |
930 |
|
|
|
931 |
chrisfen |
626 |
pref = pre22 * mu_i * mu_j |
932 |
chrisfen |
597 |
vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2) |
933 |
chrisfen |
626 |
vpair = vpair + vterm |
934 |
|
|
epot = epot + sw*vterm |
935 |
chrisfen |
597 |
|
936 |
|
|
a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij |
937 |
|
|
|
938 |
chrisfen |
626 |
dudx = dudx + sw*pref*3.0d0*ri4*scale & |
939 |
|
|
*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) |
940 |
|
|
dudy = dudy + sw*pref*3.0d0*ri4*scale & |
941 |
|
|
*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) |
942 |
|
|
dudz = dudz + sw*pref*3.0d0*ri4*scale & |
943 |
|
|
*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) |
944 |
chrisfen |
597 |
|
945 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref*ri3 & |
946 |
|
|
*(uz_j(1) - 3.0d0*ct_j*xhat*sc2) |
947 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref*ri3 & |
948 |
|
|
*(uz_j(2) - 3.0d0*ct_j*yhat*sc2) |
949 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref*ri3 & |
950 |
|
|
*(uz_j(3) - 3.0d0*ct_j*zhat*sc2) |
951 |
chrisfen |
597 |
|
952 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) + sw*pref*ri3 & |
953 |
|
|
*(uz_i(1) - 3.0d0*ct_i*xhat*sc2) |
954 |
|
|
duduz_j(2) = duduz_j(2) + sw*pref*ri3 & |
955 |
|
|
*(uz_i(2) - 3.0d0*ct_i*yhat*sc2) |
956 |
|
|
duduz_j(3) = duduz_j(3) + sw*pref*ri3 & |
957 |
|
|
*(uz_i(3) - 3.0d0*ct_i*zhat*sc2) |
958 |
chrisfen |
597 |
endif |
959 |
gezelter |
411 |
endif |
960 |
|
|
endif |
961 |
gezelter |
439 |
|
962 |
|
|
if (i_is_Quadrupole) then |
963 |
|
|
if (j_is_Charge) then |
964 |
gezelter |
507 |
|
965 |
gezelter |
439 |
ri2 = riji * riji |
966 |
|
|
ri3 = ri2 * riji |
967 |
gezelter |
440 |
ri4 = ri2 * ri2 |
968 |
gezelter |
439 |
cx2 = cx_i * cx_i |
969 |
|
|
cy2 = cy_i * cy_i |
970 |
|
|
cz2 = cz_i * cz_i |
971 |
gezelter |
507 |
|
972 |
chrisfen |
611 |
if (summationMethod .eq. UNDAMPED_WOLF) then |
973 |
chrisfen |
626 |
pref = pre14 * q_j / 3.0_dp |
974 |
chrisfen |
597 |
vterm1 = pref * ri3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + & |
975 |
|
|
qyy_i * (3.0_dp*cy2 - 1.0_dp) + & |
976 |
|
|
qzz_i * (3.0_dp*cz2 - 1.0_dp) ) |
977 |
|
|
vterm2 = pref * rcuti3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + & |
978 |
|
|
qyy_i * (3.0_dp*cy2 - 1.0_dp) + & |
979 |
|
|
qzz_i * (3.0_dp*cz2 - 1.0_dp) ) |
980 |
chrisfen |
626 |
vpair = vpair + ( vterm1 - vterm2 ) |
981 |
|
|
epot = epot + sw*( vterm1 - vterm2 ) |
982 |
chrisfen |
597 |
|
983 |
chrisfen |
626 |
dudx = dudx - sw*(5.0_dp*(vterm1*riji*xhat-vterm2*rcuti2*d(1))) +& |
984 |
|
|
sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(1)) - & |
985 |
chrisfen |
597 |
qxx_i*2.0_dp*(xhat - rcuti*d(1))) + & |
986 |
|
|
(ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(1)) - & |
987 |
|
|
qyy_i*2.0_dp*(xhat - rcuti*d(1))) + & |
988 |
|
|
(ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(1)) - & |
989 |
|
|
qzz_i*2.0_dp*(xhat - rcuti*d(1))) ) |
990 |
chrisfen |
626 |
dudy = dudy - sw*(5.0_dp*(vterm1*riji*yhat-vterm2*rcuti2*d(2))) +& |
991 |
|
|
sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(2)) - & |
992 |
chrisfen |
597 |
qxx_i*2.0_dp*(yhat - rcuti*d(2))) + & |
993 |
|
|
(ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(2)) - & |
994 |
|
|
qyy_i*2.0_dp*(yhat - rcuti*d(2))) + & |
995 |
|
|
(ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(2)) - & |
996 |
|
|
qzz_i*2.0_dp*(yhat - rcuti*d(2))) ) |
997 |
chrisfen |
626 |
dudz = dudz - sw*(5.0_dp*(vterm1*riji*zhat-vterm2*rcuti2*d(3))) +& |
998 |
|
|
sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(3)) - & |
999 |
chrisfen |
597 |
qxx_i*2.0_dp*(zhat - rcuti*d(3))) + & |
1000 |
|
|
(ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(3)) - & |
1001 |
|
|
qyy_i*2.0_dp*(zhat - rcuti*d(3))) + & |
1002 |
|
|
(ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(3)) - & |
1003 |
|
|
qzz_i*2.0_dp*(zhat - rcuti*d(3))) ) |
1004 |
|
|
|
1005 |
chrisfen |
626 |
dudux_i(1) = dudux_i(1) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*xhat) -& |
1006 |
chrisfen |
597 |
rcuti4*(qxx_i*6.0_dp*cx_i*d(1))) |
1007 |
chrisfen |
626 |
dudux_i(2) = dudux_i(2) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*yhat) -& |
1008 |
chrisfen |
597 |
rcuti4*(qxx_i*6.0_dp*cx_i*d(2))) |
1009 |
chrisfen |
626 |
dudux_i(3) = dudux_i(3) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*zhat) -& |
1010 |
chrisfen |
597 |
rcuti4*(qxx_i*6.0_dp*cx_i*d(3))) |
1011 |
|
|
|
1012 |
chrisfen |
626 |
duduy_i(1) = duduy_i(1) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*xhat) -& |
1013 |
chrisfen |
597 |
rcuti4*(qyy_i*6.0_dp*cx_i*d(1))) |
1014 |
chrisfen |
626 |
duduy_i(2) = duduy_i(2) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*yhat) -& |
1015 |
chrisfen |
597 |
rcuti4*(qyy_i*6.0_dp*cx_i*d(2))) |
1016 |
chrisfen |
626 |
duduy_i(3) = duduy_i(3) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*zhat) -& |
1017 |
chrisfen |
597 |
rcuti4*(qyy_i*6.0_dp*cx_i*d(3))) |
1018 |
|
|
|
1019 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*xhat) -& |
1020 |
chrisfen |
597 |
rcuti4*(qzz_i*6.0_dp*cx_i*d(1))) |
1021 |
chrisfen |
626 |
duduz_i(2) = duduz_i(2) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*yhat) -& |
1022 |
chrisfen |
597 |
rcuti4*(qzz_i*6.0_dp*cx_i*d(2))) |
1023 |
chrisfen |
626 |
duduz_i(3) = duduz_i(3) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*zhat) -& |
1024 |
chrisfen |
597 |
rcuti4*(qzz_i*6.0_dp*cx_i*d(3))) |
1025 |
gezelter |
507 |
|
1026 |
chrisfen |
597 |
else |
1027 |
chrisfen |
626 |
pref = pre14 * q_j / 3.0_dp |
1028 |
chrisfen |
597 |
vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + & |
1029 |
|
|
qyy_i * (3.0_dp*cy2 - 1.0_dp) + & |
1030 |
|
|
qzz_i * (3.0_dp*cz2 - 1.0_dp)) |
1031 |
chrisfen |
626 |
vpair = vpair + vterm |
1032 |
|
|
epot = epot + sw*vterm |
1033 |
chrisfen |
597 |
|
1034 |
chrisfen |
626 |
dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + sw*pref*ri4 * ( & |
1035 |
chrisfen |
597 |
qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + & |
1036 |
|
|
qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + & |
1037 |
|
|
qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) ) |
1038 |
chrisfen |
626 |
dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + sw*pref*ri4 * ( & |
1039 |
chrisfen |
597 |
qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + & |
1040 |
|
|
qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + & |
1041 |
|
|
qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) ) |
1042 |
chrisfen |
626 |
dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + sw*pref*ri4 * ( & |
1043 |
chrisfen |
597 |
qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + & |
1044 |
|
|
qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + & |
1045 |
|
|
qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) ) |
1046 |
|
|
|
1047 |
chrisfen |
626 |
dudux_i(1) = dudux_i(1) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*xhat) |
1048 |
|
|
dudux_i(2) = dudux_i(2) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*yhat) |
1049 |
|
|
dudux_i(3) = dudux_i(3) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*zhat) |
1050 |
chrisfen |
597 |
|
1051 |
chrisfen |
626 |
duduy_i(1) = duduy_i(1) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*xhat) |
1052 |
|
|
duduy_i(2) = duduy_i(2) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*yhat) |
1053 |
|
|
duduy_i(3) = duduy_i(3) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*zhat) |
1054 |
chrisfen |
597 |
|
1055 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*xhat) |
1056 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*yhat) |
1057 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*zhat) |
1058 |
chrisfen |
597 |
endif |
1059 |
gezelter |
439 |
endif |
1060 |
|
|
endif |
1061 |
gezelter |
507 |
|
1062 |
|
|
|
1063 |
gezelter |
411 |
if (do_pot) then |
1064 |
|
|
#ifdef IS_MPI |
1065 |
|
|
pot_row(atom1) = pot_row(atom1) + 0.5d0*epot |
1066 |
|
|
pot_col(atom2) = pot_col(atom2) + 0.5d0*epot |
1067 |
|
|
#else |
1068 |
|
|
pot = pot + epot |
1069 |
|
|
#endif |
1070 |
|
|
endif |
1071 |
gezelter |
507 |
|
1072 |
gezelter |
411 |
#ifdef IS_MPI |
1073 |
|
|
f_Row(1,atom1) = f_Row(1,atom1) + dudx |
1074 |
|
|
f_Row(2,atom1) = f_Row(2,atom1) + dudy |
1075 |
|
|
f_Row(3,atom1) = f_Row(3,atom1) + dudz |
1076 |
gezelter |
507 |
|
1077 |
gezelter |
411 |
f_Col(1,atom2) = f_Col(1,atom2) - dudx |
1078 |
|
|
f_Col(2,atom2) = f_Col(2,atom2) - dudy |
1079 |
|
|
f_Col(3,atom2) = f_Col(3,atom2) - dudz |
1080 |
gezelter |
507 |
|
1081 |
gezelter |
411 |
if (i_is_Dipole .or. i_is_Quadrupole) then |
1082 |
gezelter |
439 |
t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2) |
1083 |
|
|
t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3) |
1084 |
|
|
t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1) |
1085 |
gezelter |
411 |
endif |
1086 |
gezelter |
439 |
if (i_is_Quadrupole) then |
1087 |
|
|
t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2) |
1088 |
|
|
t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3) |
1089 |
|
|
t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1) |
1090 |
gezelter |
411 |
|
1091 |
gezelter |
439 |
t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2) |
1092 |
|
|
t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3) |
1093 |
|
|
t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1) |
1094 |
|
|
endif |
1095 |
|
|
|
1096 |
gezelter |
411 |
if (j_is_Dipole .or. j_is_Quadrupole) then |
1097 |
gezelter |
439 |
t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2) |
1098 |
|
|
t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3) |
1099 |
|
|
t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1) |
1100 |
gezelter |
411 |
endif |
1101 |
gezelter |
439 |
if (j_is_Quadrupole) then |
1102 |
|
|
t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2) |
1103 |
|
|
t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3) |
1104 |
|
|
t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1) |
1105 |
gezelter |
411 |
|
1106 |
gezelter |
439 |
t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2) |
1107 |
|
|
t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3) |
1108 |
|
|
t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1) |
1109 |
|
|
endif |
1110 |
|
|
|
1111 |
gezelter |
411 |
#else |
1112 |
|
|
f(1,atom1) = f(1,atom1) + dudx |
1113 |
|
|
f(2,atom1) = f(2,atom1) + dudy |
1114 |
|
|
f(3,atom1) = f(3,atom1) + dudz |
1115 |
gezelter |
507 |
|
1116 |
gezelter |
411 |
f(1,atom2) = f(1,atom2) - dudx |
1117 |
|
|
f(2,atom2) = f(2,atom2) - dudy |
1118 |
|
|
f(3,atom2) = f(3,atom2) - dudz |
1119 |
gezelter |
507 |
|
1120 |
gezelter |
411 |
if (i_is_Dipole .or. i_is_Quadrupole) then |
1121 |
gezelter |
439 |
t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2) |
1122 |
|
|
t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3) |
1123 |
|
|
t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1) |
1124 |
gezelter |
411 |
endif |
1125 |
gezelter |
439 |
if (i_is_Quadrupole) then |
1126 |
|
|
t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2) |
1127 |
|
|
t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3) |
1128 |
|
|
t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1) |
1129 |
|
|
|
1130 |
|
|
t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2) |
1131 |
|
|
t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3) |
1132 |
|
|
t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1) |
1133 |
|
|
endif |
1134 |
|
|
|
1135 |
gezelter |
411 |
if (j_is_Dipole .or. j_is_Quadrupole) then |
1136 |
gezelter |
439 |
t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2) |
1137 |
|
|
t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3) |
1138 |
|
|
t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1) |
1139 |
gezelter |
411 |
endif |
1140 |
gezelter |
439 |
if (j_is_Quadrupole) then |
1141 |
|
|
t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2) |
1142 |
|
|
t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3) |
1143 |
|
|
t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1) |
1144 |
|
|
|
1145 |
|
|
t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2) |
1146 |
|
|
t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3) |
1147 |
|
|
t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1) |
1148 |
|
|
endif |
1149 |
|
|
|
1150 |
gezelter |
411 |
#endif |
1151 |
gezelter |
507 |
|
1152 |
gezelter |
411 |
#ifdef IS_MPI |
1153 |
|
|
id1 = AtomRowToGlobal(atom1) |
1154 |
|
|
id2 = AtomColToGlobal(atom2) |
1155 |
|
|
#else |
1156 |
|
|
id1 = atom1 |
1157 |
|
|
id2 = atom2 |
1158 |
|
|
#endif |
1159 |
|
|
|
1160 |
|
|
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
1161 |
gezelter |
507 |
|
1162 |
gezelter |
411 |
fpair(1) = fpair(1) + dudx |
1163 |
|
|
fpair(2) = fpair(2) + dudy |
1164 |
|
|
fpair(3) = fpair(3) + dudz |
1165 |
|
|
|
1166 |
|
|
endif |
1167 |
|
|
|
1168 |
|
|
return |
1169 |
|
|
end subroutine doElectrostaticPair |
1170 |
chuckv |
492 |
|
1171 |
chrisfen |
532 |
!! calculates the switching functions and their derivatives for a given |
1172 |
|
|
subroutine calc_switch(r, mu, scale, dscale) |
1173 |
gezelter |
507 |
|
1174 |
chrisfen |
532 |
real (kind=dp), intent(in) :: r, mu |
1175 |
|
|
real (kind=dp), intent(inout) :: scale, dscale |
1176 |
|
|
real (kind=dp) :: rl, ru, mulow, minRatio, temp, scaleVal |
1177 |
|
|
|
1178 |
|
|
! distances must be in angstroms |
1179 |
|
|
rl = 2.75d0 |
1180 |
chrisfen |
534 |
ru = 3.75d0 |
1181 |
|
|
mulow = 0.0d0 !3.3856d0 ! 1.84 * 1.84 |
1182 |
chrisfen |
532 |
minRatio = mulow / (mu*mu) |
1183 |
|
|
scaleVal = 1.0d0 - minRatio |
1184 |
|
|
|
1185 |
|
|
if (r.lt.rl) then |
1186 |
|
|
scale = minRatio |
1187 |
|
|
dscale = 0.0d0 |
1188 |
|
|
elseif (r.gt.ru) then |
1189 |
|
|
scale = 1.0d0 |
1190 |
|
|
dscale = 0.0d0 |
1191 |
|
|
else |
1192 |
|
|
scale = 1.0d0 - scaleVal*((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) & |
1193 |
|
|
/ ((ru - rl)**3) |
1194 |
|
|
dscale = -scaleVal * 6.0d0 * (r-ru)*(r-rl)/((ru - rl)**3) |
1195 |
|
|
endif |
1196 |
|
|
|
1197 |
|
|
return |
1198 |
|
|
end subroutine calc_switch |
1199 |
|
|
|
1200 |
chuckv |
492 |
subroutine destroyElectrostaticTypes() |
1201 |
|
|
|
1202 |
gezelter |
507 |
if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap) |
1203 |
|
|
|
1204 |
chuckv |
492 |
end subroutine destroyElectrostaticTypes |
1205 |
|
|
|
1206 |
gezelter |
411 |
end module electrostatic_module |