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 |
chuckv |
656 |
|
58 |
gezelter |
602 |
#define __FORTRAN90 |
59 |
chuckv |
656 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
60 |
gezelter |
602 |
#include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" |
61 |
chrisfen |
710 |
#include "UseTheForce/DarkSide/fScreeningMethod.h" |
62 |
gezelter |
602 |
|
63 |
chuckv |
656 |
|
64 |
gezelter |
434 |
!! these prefactors convert the multipole interactions into kcal / mol |
65 |
|
|
!! all were computed assuming distances are measured in angstroms |
66 |
|
|
!! Charge-Charge, assuming charges are measured in electrons |
67 |
gezelter |
411 |
real(kind=dp), parameter :: pre11 = 332.0637778_dp |
68 |
gezelter |
434 |
!! Charge-Dipole, assuming charges are measured in electrons, and |
69 |
|
|
!! dipoles are measured in debyes |
70 |
|
|
real(kind=dp), parameter :: pre12 = 69.13373_dp |
71 |
|
|
!! Dipole-Dipole, assuming dipoles are measured in debyes |
72 |
|
|
real(kind=dp), parameter :: pre22 = 14.39325_dp |
73 |
|
|
!! Charge-Quadrupole, assuming charges are measured in electrons, and |
74 |
|
|
!! quadrupoles are measured in 10^-26 esu cm^2 |
75 |
|
|
!! This unit is also known affectionately as an esu centi-barn. |
76 |
|
|
real(kind=dp), parameter :: pre14 = 69.13373_dp |
77 |
gezelter |
411 |
|
78 |
chrisfen |
712 |
!! variables to handle different summation methods for long-range |
79 |
|
|
!! electrostatics: |
80 |
gezelter |
602 |
integer, save :: summationMethod = NONE |
81 |
chrisfen |
710 |
integer, save :: screeningMethod = UNDAMPED |
82 |
chrisfen |
603 |
logical, save :: summationMethodChecked = .false. |
83 |
gezelter |
602 |
real(kind=DP), save :: defaultCutoff = 0.0_DP |
84 |
chrisfen |
682 |
real(kind=DP), save :: defaultCutoff2 = 0.0_DP |
85 |
gezelter |
602 |
logical, save :: haveDefaultCutoff = .false. |
86 |
|
|
real(kind=DP), save :: dampingAlpha = 0.0_DP |
87 |
|
|
logical, save :: haveDampingAlpha = .false. |
88 |
chrisfen |
682 |
real(kind=DP), save :: dielectric = 1.0_DP |
89 |
gezelter |
602 |
logical, save :: haveDielectric = .false. |
90 |
|
|
real(kind=DP), save :: constERFC = 0.0_DP |
91 |
|
|
real(kind=DP), save :: constEXP = 0.0_DP |
92 |
chrisfen |
682 |
real(kind=dp), save :: rcuti = 0.0_DP |
93 |
|
|
real(kind=dp), save :: rcuti2 = 0.0_DP |
94 |
|
|
real(kind=dp), save :: rcuti3 = 0.0_DP |
95 |
|
|
real(kind=dp), save :: rcuti4 = 0.0_DP |
96 |
|
|
real(kind=dp), save :: alphaPi = 0.0_DP |
97 |
|
|
real(kind=dp), save :: invRootPi = 0.0_DP |
98 |
|
|
real(kind=dp), save :: rrf = 1.0_DP |
99 |
|
|
real(kind=dp), save :: rt = 1.0_DP |
100 |
|
|
real(kind=dp), save :: rrfsq = 1.0_DP |
101 |
|
|
real(kind=dp), save :: preRF = 0.0_DP |
102 |
chrisfen |
695 |
real(kind=dp), save :: preRF2 = 0.0_DP |
103 |
chrisfen |
682 |
|
104 |
chuckv |
632 |
#ifdef __IFC |
105 |
|
|
! error function for ifc version > 7. |
106 |
chuckv |
631 |
double precision, external :: derfc |
107 |
chuckv |
632 |
#endif |
108 |
|
|
|
109 |
gezelter |
602 |
public :: setElectrostaticSummationMethod |
110 |
chrisfen |
710 |
public :: setScreeningMethod |
111 |
gezelter |
602 |
public :: setElectrostaticCutoffRadius |
112 |
chrisfen |
710 |
public :: setDampingAlpha |
113 |
gezelter |
602 |
public :: setReactionFieldDielectric |
114 |
gezelter |
411 |
public :: newElectrostaticType |
115 |
|
|
public :: setCharge |
116 |
|
|
public :: setDipoleMoment |
117 |
|
|
public :: setSplitDipoleDistance |
118 |
|
|
public :: setQuadrupoleMoments |
119 |
|
|
public :: doElectrostaticPair |
120 |
|
|
public :: getCharge |
121 |
|
|
public :: getDipoleMoment |
122 |
chuckv |
492 |
public :: destroyElectrostaticTypes |
123 |
chrisfen |
703 |
public :: self_self |
124 |
chrisfen |
700 |
public :: rf_self_excludes |
125 |
gezelter |
411 |
|
126 |
|
|
type :: Electrostatic |
127 |
|
|
integer :: c_ident |
128 |
|
|
logical :: is_Charge = .false. |
129 |
|
|
logical :: is_Dipole = .false. |
130 |
|
|
logical :: is_SplitDipole = .false. |
131 |
|
|
logical :: is_Quadrupole = .false. |
132 |
chrisfen |
532 |
logical :: is_Tap = .false. |
133 |
gezelter |
411 |
real(kind=DP) :: charge = 0.0_DP |
134 |
|
|
real(kind=DP) :: dipole_moment = 0.0_DP |
135 |
|
|
real(kind=DP) :: split_dipole_distance = 0.0_DP |
136 |
|
|
real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP |
137 |
|
|
end type Electrostatic |
138 |
|
|
|
139 |
|
|
type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap |
140 |
|
|
|
141 |
|
|
contains |
142 |
|
|
|
143 |
gezelter |
602 |
subroutine setElectrostaticSummationMethod(the_ESM) |
144 |
|
|
integer, intent(in) :: the_ESM |
145 |
|
|
|
146 |
|
|
if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then |
147 |
|
|
call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method") |
148 |
|
|
endif |
149 |
|
|
|
150 |
chrisfen |
610 |
summationMethod = the_ESM |
151 |
chrisfen |
626 |
|
152 |
gezelter |
602 |
end subroutine setElectrostaticSummationMethod |
153 |
|
|
|
154 |
chrisfen |
710 |
subroutine setScreeningMethod(the_SM) |
155 |
|
|
integer, intent(in) :: the_SM |
156 |
|
|
screeningMethod = the_SM |
157 |
|
|
end subroutine setScreeningMethod |
158 |
|
|
|
159 |
chrisfen |
682 |
subroutine setElectrostaticCutoffRadius(thisRcut, thisRsw) |
160 |
gezelter |
602 |
real(kind=dp), intent(in) :: thisRcut |
161 |
chrisfen |
682 |
real(kind=dp), intent(in) :: thisRsw |
162 |
gezelter |
602 |
defaultCutoff = thisRcut |
163 |
chrisfen |
682 |
rrf = defaultCutoff |
164 |
|
|
rt = thisRsw |
165 |
gezelter |
602 |
haveDefaultCutoff = .true. |
166 |
|
|
end subroutine setElectrostaticCutoffRadius |
167 |
|
|
|
168 |
chrisfen |
710 |
subroutine setDampingAlpha(thisAlpha) |
169 |
gezelter |
602 |
real(kind=dp), intent(in) :: thisAlpha |
170 |
|
|
dampingAlpha = thisAlpha |
171 |
|
|
haveDampingAlpha = .true. |
172 |
chrisfen |
710 |
end subroutine setDampingAlpha |
173 |
gezelter |
602 |
|
174 |
|
|
subroutine setReactionFieldDielectric(thisDielectric) |
175 |
|
|
real(kind=dp), intent(in) :: thisDielectric |
176 |
|
|
dielectric = thisDielectric |
177 |
|
|
haveDielectric = .true. |
178 |
|
|
end subroutine setReactionFieldDielectric |
179 |
|
|
|
180 |
gezelter |
411 |
subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, & |
181 |
chrisfen |
532 |
is_SplitDipole, is_Quadrupole, is_Tap, status) |
182 |
gezelter |
507 |
|
183 |
gezelter |
411 |
integer, intent(in) :: c_ident |
184 |
|
|
logical, intent(in) :: is_Charge |
185 |
|
|
logical, intent(in) :: is_Dipole |
186 |
|
|
logical, intent(in) :: is_SplitDipole |
187 |
|
|
logical, intent(in) :: is_Quadrupole |
188 |
chrisfen |
532 |
logical, intent(in) :: is_Tap |
189 |
gezelter |
411 |
integer, intent(out) :: status |
190 |
|
|
integer :: nAtypes, myATID, i, j |
191 |
|
|
|
192 |
|
|
status = 0 |
193 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
194 |
gezelter |
507 |
|
195 |
gezelter |
411 |
!! Be simple-minded and assume that we need an ElectrostaticMap that |
196 |
|
|
!! is the same size as the total number of atom types |
197 |
|
|
|
198 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
199 |
gezelter |
507 |
|
200 |
gezelter |
411 |
nAtypes = getSize(atypes) |
201 |
gezelter |
507 |
|
202 |
gezelter |
411 |
if (nAtypes == 0) then |
203 |
|
|
status = -1 |
204 |
|
|
return |
205 |
|
|
end if |
206 |
gezelter |
507 |
|
207 |
gezelter |
411 |
if (.not. allocated(ElectrostaticMap)) then |
208 |
|
|
allocate(ElectrostaticMap(nAtypes)) |
209 |
|
|
endif |
210 |
gezelter |
507 |
|
211 |
gezelter |
411 |
end if |
212 |
|
|
|
213 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
214 |
|
|
status = -1 |
215 |
|
|
return |
216 |
|
|
endif |
217 |
gezelter |
507 |
|
218 |
gezelter |
411 |
! set the values for ElectrostaticMap for this atom type: |
219 |
|
|
|
220 |
|
|
ElectrostaticMap(myATID)%c_ident = c_ident |
221 |
|
|
ElectrostaticMap(myATID)%is_Charge = is_Charge |
222 |
|
|
ElectrostaticMap(myATID)%is_Dipole = is_Dipole |
223 |
|
|
ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole |
224 |
|
|
ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole |
225 |
chrisfen |
532 |
ElectrostaticMap(myATID)%is_Tap = is_Tap |
226 |
gezelter |
507 |
|
227 |
gezelter |
411 |
end subroutine newElectrostaticType |
228 |
|
|
|
229 |
|
|
subroutine setCharge(c_ident, charge, status) |
230 |
|
|
integer, intent(in) :: c_ident |
231 |
|
|
real(kind=dp), intent(in) :: charge |
232 |
|
|
integer, intent(out) :: status |
233 |
|
|
integer :: myATID |
234 |
|
|
|
235 |
|
|
status = 0 |
236 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
237 |
|
|
|
238 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
239 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!") |
240 |
|
|
status = -1 |
241 |
|
|
return |
242 |
|
|
end if |
243 |
|
|
|
244 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
245 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!") |
246 |
|
|
status = -1 |
247 |
|
|
return |
248 |
|
|
endif |
249 |
|
|
|
250 |
|
|
if (.not.ElectrostaticMap(myATID)%is_Charge) then |
251 |
|
|
call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!") |
252 |
|
|
status = -1 |
253 |
|
|
return |
254 |
gezelter |
507 |
endif |
255 |
gezelter |
411 |
|
256 |
|
|
ElectrostaticMap(myATID)%charge = charge |
257 |
|
|
end subroutine setCharge |
258 |
|
|
|
259 |
|
|
subroutine setDipoleMoment(c_ident, dipole_moment, status) |
260 |
|
|
integer, intent(in) :: c_ident |
261 |
|
|
real(kind=dp), intent(in) :: dipole_moment |
262 |
|
|
integer, intent(out) :: status |
263 |
|
|
integer :: myATID |
264 |
|
|
|
265 |
|
|
status = 0 |
266 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
267 |
|
|
|
268 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
269 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!") |
270 |
|
|
status = -1 |
271 |
|
|
return |
272 |
|
|
end if |
273 |
|
|
|
274 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
275 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!") |
276 |
|
|
status = -1 |
277 |
|
|
return |
278 |
|
|
endif |
279 |
|
|
|
280 |
|
|
if (.not.ElectrostaticMap(myATID)%is_Dipole) then |
281 |
|
|
call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!") |
282 |
|
|
status = -1 |
283 |
|
|
return |
284 |
|
|
endif |
285 |
|
|
|
286 |
|
|
ElectrostaticMap(myATID)%dipole_moment = dipole_moment |
287 |
|
|
end subroutine setDipoleMoment |
288 |
|
|
|
289 |
|
|
subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status) |
290 |
|
|
integer, intent(in) :: c_ident |
291 |
|
|
real(kind=dp), intent(in) :: split_dipole_distance |
292 |
|
|
integer, intent(out) :: status |
293 |
|
|
integer :: myATID |
294 |
|
|
|
295 |
|
|
status = 0 |
296 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
297 |
|
|
|
298 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
299 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!") |
300 |
|
|
status = -1 |
301 |
|
|
return |
302 |
|
|
end if |
303 |
|
|
|
304 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
305 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!") |
306 |
|
|
status = -1 |
307 |
|
|
return |
308 |
|
|
endif |
309 |
|
|
|
310 |
|
|
if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then |
311 |
|
|
call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!") |
312 |
|
|
status = -1 |
313 |
|
|
return |
314 |
|
|
endif |
315 |
|
|
|
316 |
|
|
ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance |
317 |
|
|
end subroutine setSplitDipoleDistance |
318 |
|
|
|
319 |
|
|
subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status) |
320 |
|
|
integer, intent(in) :: c_ident |
321 |
|
|
real(kind=dp), intent(in), dimension(3) :: quadrupole_moments |
322 |
|
|
integer, intent(out) :: status |
323 |
|
|
integer :: myATID, i, j |
324 |
|
|
|
325 |
|
|
status = 0 |
326 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
327 |
|
|
|
328 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
329 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!") |
330 |
|
|
status = -1 |
331 |
|
|
return |
332 |
|
|
end if |
333 |
|
|
|
334 |
|
|
if (myATID .gt. size(ElectrostaticMap)) then |
335 |
|
|
call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!") |
336 |
|
|
status = -1 |
337 |
|
|
return |
338 |
|
|
endif |
339 |
|
|
|
340 |
|
|
if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then |
341 |
|
|
call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!") |
342 |
|
|
status = -1 |
343 |
|
|
return |
344 |
|
|
endif |
345 |
gezelter |
507 |
|
346 |
gezelter |
411 |
do i = 1, 3 |
347 |
gezelter |
507 |
ElectrostaticMap(myATID)%quadrupole_moments(i) = & |
348 |
|
|
quadrupole_moments(i) |
349 |
|
|
enddo |
350 |
gezelter |
411 |
|
351 |
|
|
end subroutine setQuadrupoleMoments |
352 |
|
|
|
353 |
gezelter |
507 |
|
354 |
gezelter |
411 |
function getCharge(atid) result (c) |
355 |
|
|
integer, intent(in) :: atid |
356 |
|
|
integer :: localError |
357 |
|
|
real(kind=dp) :: c |
358 |
gezelter |
507 |
|
359 |
gezelter |
411 |
if (.not.allocated(ElectrostaticMap)) then |
360 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!") |
361 |
|
|
return |
362 |
|
|
end if |
363 |
gezelter |
507 |
|
364 |
gezelter |
411 |
if (.not.ElectrostaticMap(atid)%is_Charge) then |
365 |
|
|
call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!") |
366 |
|
|
return |
367 |
|
|
endif |
368 |
gezelter |
507 |
|
369 |
gezelter |
411 |
c = ElectrostaticMap(atid)%charge |
370 |
|
|
end function getCharge |
371 |
|
|
|
372 |
|
|
function getDipoleMoment(atid) result (dm) |
373 |
|
|
integer, intent(in) :: atid |
374 |
|
|
integer :: localError |
375 |
|
|
real(kind=dp) :: dm |
376 |
gezelter |
507 |
|
377 |
gezelter |
411 |
if (.not.allocated(ElectrostaticMap)) then |
378 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!") |
379 |
|
|
return |
380 |
|
|
end if |
381 |
gezelter |
507 |
|
382 |
gezelter |
411 |
if (.not.ElectrostaticMap(atid)%is_Dipole) then |
383 |
|
|
call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!") |
384 |
|
|
return |
385 |
|
|
endif |
386 |
gezelter |
507 |
|
387 |
gezelter |
411 |
dm = ElectrostaticMap(atid)%dipole_moment |
388 |
|
|
end function getDipoleMoment |
389 |
|
|
|
390 |
gezelter |
602 |
subroutine checkSummationMethod() |
391 |
|
|
|
392 |
chrisfen |
607 |
if (.not.haveDefaultCutoff) then |
393 |
|
|
call handleError("checkSummationMethod", "no Default Cutoff set!") |
394 |
|
|
endif |
395 |
|
|
|
396 |
|
|
rcuti = 1.0d0 / defaultCutoff |
397 |
|
|
rcuti2 = rcuti*rcuti |
398 |
|
|
rcuti3 = rcuti2*rcuti |
399 |
|
|
rcuti4 = rcuti2*rcuti2 |
400 |
|
|
|
401 |
chrisfen |
710 |
if (screeningMethod .eq. DAMPED) then |
402 |
chrisfen |
703 |
if (.not.haveDampingAlpha) then |
403 |
|
|
call handleError("checkSummationMethod", "no Damping Alpha set!") |
404 |
|
|
endif |
405 |
|
|
|
406 |
|
|
if (.not.haveDefaultCutoff) then |
407 |
|
|
call handleError("checkSummationMethod", "no Default Cutoff set!") |
408 |
|
|
endif |
409 |
chrisfen |
603 |
|
410 |
chrisfen |
703 |
constEXP = exp(-dampingAlpha*dampingAlpha*defaultCutoff*defaultCutoff) |
411 |
|
|
constERFC = derfc(dampingAlpha*defaultCutoff) |
412 |
|
|
invRootPi = 0.56418958354775628695d0 |
413 |
|
|
alphaPi = 2*dampingAlpha*invRootPi |
414 |
|
|
|
415 |
gezelter |
602 |
endif |
416 |
|
|
|
417 |
chrisfen |
603 |
if (summationMethod .eq. REACTION_FIELD) then |
418 |
chrisfen |
703 |
if (haveDielectric) then |
419 |
|
|
defaultCutoff2 = defaultCutoff*defaultCutoff |
420 |
|
|
preRF = (dielectric-1.0d0) / & |
421 |
|
|
((2.0d0*dielectric+1.0d0)*defaultCutoff2*defaultCutoff) |
422 |
|
|
preRF2 = 2.0d0*preRF |
423 |
|
|
else |
424 |
|
|
call handleError("checkSummationMethod", "Dielectric not set") |
425 |
chrisfen |
603 |
endif |
426 |
chrisfen |
703 |
|
427 |
chrisfen |
603 |
endif |
428 |
|
|
|
429 |
|
|
summationMethodChecked = .true. |
430 |
gezelter |
602 |
end subroutine checkSummationMethod |
431 |
|
|
|
432 |
chrisfen |
712 |
|
433 |
gezelter |
411 |
subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, & |
434 |
chrisfen |
712 |
vpair, fpair, pot, eFrame, f, t, do_pot) |
435 |
gezelter |
507 |
|
436 |
chrisfen |
703 |
logical, intent(in) :: do_pot |
437 |
gezelter |
507 |
|
438 |
gezelter |
411 |
integer, intent(in) :: atom1, atom2 |
439 |
|
|
integer :: localError |
440 |
|
|
|
441 |
chrisfen |
607 |
real(kind=dp), intent(in) :: rij, r2, sw |
442 |
gezelter |
411 |
real(kind=dp), intent(in), dimension(3) :: d |
443 |
|
|
real(kind=dp), intent(inout) :: vpair |
444 |
chrisfen |
703 |
real(kind=dp), intent(inout), dimension(3) :: fpair |
445 |
gezelter |
411 |
|
446 |
chrisfen |
626 |
real( kind = dp ) :: pot |
447 |
gezelter |
411 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
448 |
|
|
real( kind = dp ), dimension(3,nLocal) :: f |
449 |
chrisfen |
710 |
real( kind = dp ), dimension(3,nLocal) :: felec |
450 |
gezelter |
411 |
real( kind = dp ), dimension(3,nLocal) :: t |
451 |
gezelter |
507 |
|
452 |
gezelter |
439 |
real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i |
453 |
|
|
real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j |
454 |
|
|
real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i |
455 |
|
|
real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j |
456 |
gezelter |
411 |
|
457 |
|
|
logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole |
458 |
|
|
logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole |
459 |
chrisfen |
532 |
logical :: i_is_Tap, j_is_Tap |
460 |
gezelter |
411 |
integer :: me1, me2, id1, id2 |
461 |
|
|
real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j |
462 |
gezelter |
439 |
real (kind=dp) :: qxx_i, qyy_i, qzz_i |
463 |
|
|
real (kind=dp) :: qxx_j, qyy_j, qzz_j |
464 |
|
|
real (kind=dp) :: cx_i, cy_i, cz_i |
465 |
|
|
real (kind=dp) :: cx_j, cy_j, cz_j |
466 |
|
|
real (kind=dp) :: cx2, cy2, cz2 |
467 |
gezelter |
411 |
real (kind=dp) :: ct_i, ct_j, ct_ij, a1 |
468 |
gezelter |
421 |
real (kind=dp) :: riji, ri, ri2, ri3, ri4 |
469 |
chrisfen |
597 |
real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2 |
470 |
gezelter |
421 |
real (kind=dp) :: xhat, yhat, zhat |
471 |
gezelter |
411 |
real (kind=dp) :: dudx, dudy, dudz |
472 |
chrisfen |
626 |
real (kind=dp) :: scale, sc2, bigR |
473 |
chrisfen |
640 |
real (kind=dp) :: varERFC, varEXP |
474 |
chrisfen |
644 |
real (kind=dp) :: limScale |
475 |
chrisfen |
695 |
real (kind=dp) :: preVal, rfVal |
476 |
gezelter |
411 |
|
477 |
|
|
if (.not.allocated(ElectrostaticMap)) then |
478 |
|
|
call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!") |
479 |
|
|
return |
480 |
|
|
end if |
481 |
|
|
|
482 |
gezelter |
602 |
if (.not.summationMethodChecked) then |
483 |
|
|
call checkSummationMethod() |
484 |
|
|
endif |
485 |
|
|
|
486 |
gezelter |
411 |
#ifdef IS_MPI |
487 |
|
|
me1 = atid_Row(atom1) |
488 |
|
|
me2 = atid_Col(atom2) |
489 |
|
|
#else |
490 |
|
|
me1 = atid(atom1) |
491 |
|
|
me2 = atid(atom2) |
492 |
|
|
#endif |
493 |
|
|
|
494 |
chrisfen |
710 |
!!$ if (rij .ge. defaultCutoff) then |
495 |
|
|
!!$ write(*,*) 'warning: rij = ', rij, ' rcut = ', defaultCutoff, ' sw = ', sw |
496 |
|
|
!!$ endif |
497 |
|
|
|
498 |
gezelter |
411 |
!! some variables we'll need independent of electrostatic type: |
499 |
|
|
|
500 |
|
|
riji = 1.0d0 / rij |
501 |
chrisfen |
644 |
|
502 |
gezelter |
421 |
xhat = d(1) * riji |
503 |
|
|
yhat = d(2) * riji |
504 |
|
|
zhat = d(3) * riji |
505 |
gezelter |
411 |
|
506 |
|
|
!! logicals |
507 |
|
|
i_is_Charge = ElectrostaticMap(me1)%is_Charge |
508 |
|
|
i_is_Dipole = ElectrostaticMap(me1)%is_Dipole |
509 |
|
|
i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole |
510 |
|
|
i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole |
511 |
chrisfen |
532 |
i_is_Tap = ElectrostaticMap(me1)%is_Tap |
512 |
gezelter |
411 |
|
513 |
|
|
j_is_Charge = ElectrostaticMap(me2)%is_Charge |
514 |
|
|
j_is_Dipole = ElectrostaticMap(me2)%is_Dipole |
515 |
|
|
j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole |
516 |
|
|
j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole |
517 |
chrisfen |
532 |
j_is_Tap = ElectrostaticMap(me2)%is_Tap |
518 |
gezelter |
411 |
|
519 |
|
|
if (i_is_Charge) then |
520 |
|
|
q_i = ElectrostaticMap(me1)%charge |
521 |
|
|
endif |
522 |
gezelter |
507 |
|
523 |
gezelter |
411 |
if (i_is_Dipole) then |
524 |
|
|
mu_i = ElectrostaticMap(me1)%dipole_moment |
525 |
|
|
#ifdef IS_MPI |
526 |
gezelter |
439 |
uz_i(1) = eFrame_Row(3,atom1) |
527 |
|
|
uz_i(2) = eFrame_Row(6,atom1) |
528 |
|
|
uz_i(3) = eFrame_Row(9,atom1) |
529 |
gezelter |
411 |
#else |
530 |
gezelter |
439 |
uz_i(1) = eFrame(3,atom1) |
531 |
|
|
uz_i(2) = eFrame(6,atom1) |
532 |
|
|
uz_i(3) = eFrame(9,atom1) |
533 |
gezelter |
411 |
#endif |
534 |
gezelter |
439 |
ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat |
535 |
gezelter |
411 |
|
536 |
|
|
if (i_is_SplitDipole) then |
537 |
|
|
d_i = ElectrostaticMap(me1)%split_dipole_distance |
538 |
|
|
endif |
539 |
gezelter |
507 |
|
540 |
gezelter |
411 |
endif |
541 |
|
|
|
542 |
gezelter |
439 |
if (i_is_Quadrupole) then |
543 |
|
|
qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1) |
544 |
|
|
qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2) |
545 |
|
|
qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3) |
546 |
|
|
#ifdef IS_MPI |
547 |
|
|
ux_i(1) = eFrame_Row(1,atom1) |
548 |
|
|
ux_i(2) = eFrame_Row(4,atom1) |
549 |
|
|
ux_i(3) = eFrame_Row(7,atom1) |
550 |
|
|
uy_i(1) = eFrame_Row(2,atom1) |
551 |
|
|
uy_i(2) = eFrame_Row(5,atom1) |
552 |
|
|
uy_i(3) = eFrame_Row(8,atom1) |
553 |
|
|
uz_i(1) = eFrame_Row(3,atom1) |
554 |
|
|
uz_i(2) = eFrame_Row(6,atom1) |
555 |
|
|
uz_i(3) = eFrame_Row(9,atom1) |
556 |
|
|
#else |
557 |
|
|
ux_i(1) = eFrame(1,atom1) |
558 |
|
|
ux_i(2) = eFrame(4,atom1) |
559 |
|
|
ux_i(3) = eFrame(7,atom1) |
560 |
|
|
uy_i(1) = eFrame(2,atom1) |
561 |
|
|
uy_i(2) = eFrame(5,atom1) |
562 |
|
|
uy_i(3) = eFrame(8,atom1) |
563 |
|
|
uz_i(1) = eFrame(3,atom1) |
564 |
|
|
uz_i(2) = eFrame(6,atom1) |
565 |
|
|
uz_i(3) = eFrame(9,atom1) |
566 |
|
|
#endif |
567 |
|
|
cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat |
568 |
|
|
cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat |
569 |
|
|
cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat |
570 |
|
|
endif |
571 |
|
|
|
572 |
gezelter |
411 |
if (j_is_Charge) then |
573 |
|
|
q_j = ElectrostaticMap(me2)%charge |
574 |
|
|
endif |
575 |
gezelter |
507 |
|
576 |
gezelter |
411 |
if (j_is_Dipole) then |
577 |
|
|
mu_j = ElectrostaticMap(me2)%dipole_moment |
578 |
|
|
#ifdef IS_MPI |
579 |
gezelter |
439 |
uz_j(1) = eFrame_Col(3,atom2) |
580 |
|
|
uz_j(2) = eFrame_Col(6,atom2) |
581 |
|
|
uz_j(3) = eFrame_Col(9,atom2) |
582 |
gezelter |
411 |
#else |
583 |
gezelter |
439 |
uz_j(1) = eFrame(3,atom2) |
584 |
|
|
uz_j(2) = eFrame(6,atom2) |
585 |
|
|
uz_j(3) = eFrame(9,atom2) |
586 |
gezelter |
411 |
#endif |
587 |
chrisfen |
465 |
ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat |
588 |
gezelter |
411 |
|
589 |
|
|
if (j_is_SplitDipole) then |
590 |
|
|
d_j = ElectrostaticMap(me2)%split_dipole_distance |
591 |
|
|
endif |
592 |
|
|
endif |
593 |
|
|
|
594 |
gezelter |
439 |
if (j_is_Quadrupole) then |
595 |
|
|
qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1) |
596 |
|
|
qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2) |
597 |
|
|
qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3) |
598 |
|
|
#ifdef IS_MPI |
599 |
|
|
ux_j(1) = eFrame_Col(1,atom2) |
600 |
|
|
ux_j(2) = eFrame_Col(4,atom2) |
601 |
|
|
ux_j(3) = eFrame_Col(7,atom2) |
602 |
|
|
uy_j(1) = eFrame_Col(2,atom2) |
603 |
|
|
uy_j(2) = eFrame_Col(5,atom2) |
604 |
|
|
uy_j(3) = eFrame_Col(8,atom2) |
605 |
|
|
uz_j(1) = eFrame_Col(3,atom2) |
606 |
|
|
uz_j(2) = eFrame_Col(6,atom2) |
607 |
|
|
uz_j(3) = eFrame_Col(9,atom2) |
608 |
|
|
#else |
609 |
|
|
ux_j(1) = eFrame(1,atom2) |
610 |
|
|
ux_j(2) = eFrame(4,atom2) |
611 |
|
|
ux_j(3) = eFrame(7,atom2) |
612 |
|
|
uy_j(1) = eFrame(2,atom2) |
613 |
|
|
uy_j(2) = eFrame(5,atom2) |
614 |
|
|
uy_j(3) = eFrame(8,atom2) |
615 |
|
|
uz_j(1) = eFrame(3,atom2) |
616 |
|
|
uz_j(2) = eFrame(6,atom2) |
617 |
|
|
uz_j(3) = eFrame(9,atom2) |
618 |
|
|
#endif |
619 |
|
|
cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat |
620 |
|
|
cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat |
621 |
|
|
cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat |
622 |
|
|
endif |
623 |
chrisfen |
554 |
|
624 |
gezelter |
411 |
epot = 0.0_dp |
625 |
|
|
dudx = 0.0_dp |
626 |
|
|
dudy = 0.0_dp |
627 |
|
|
dudz = 0.0_dp |
628 |
|
|
|
629 |
gezelter |
439 |
dudux_i = 0.0_dp |
630 |
|
|
duduy_i = 0.0_dp |
631 |
|
|
duduz_i = 0.0_dp |
632 |
gezelter |
411 |
|
633 |
gezelter |
439 |
dudux_j = 0.0_dp |
634 |
|
|
duduy_j = 0.0_dp |
635 |
|
|
duduz_j = 0.0_dp |
636 |
gezelter |
411 |
|
637 |
|
|
if (i_is_Charge) then |
638 |
|
|
|
639 |
|
|
if (j_is_Charge) then |
640 |
gezelter |
507 |
|
641 |
chrisfen |
710 |
if (summationMethod .eq. SHIFTED_POTENTIAL) then |
642 |
|
|
|
643 |
chrisfen |
597 |
vterm = pre11 * q_i * q_j * (riji - rcuti) |
644 |
|
|
vpair = vpair + vterm |
645 |
chrisfen |
626 |
epot = epot + sw*vterm |
646 |
chrisfen |
597 |
|
647 |
chrisfen |
710 |
dudr = -sw*pre11*q_i*q_j * (riji*riji) |
648 |
chrisfen |
597 |
|
649 |
chrisfen |
703 |
dudx = dudx + dudr * xhat |
650 |
|
|
dudy = dudy + dudr * yhat |
651 |
|
|
dudz = dudz + dudr * zhat |
652 |
gezelter |
411 |
|
653 |
chrisfen |
710 |
!!$ elseif (summationMethod .eq. DAMPED_WOLF) then |
654 |
|
|
!!$ varERFC = derfc(dampingAlpha*rij) |
655 |
|
|
!!$ varEXP = exp(-dampingAlpha*dampingAlpha*rij*rij) |
656 |
|
|
!!$ vterm = pre11 * q_i * q_j * (varERFC*riji - constERFC*rcuti) |
657 |
|
|
!!$ vpair = vpair + vterm |
658 |
|
|
!!$ epot = epot + sw*vterm |
659 |
|
|
!!$ |
660 |
|
|
!!$! dudr = -sw*pre11*q_i*q_j * (((varERFC*riji*riji & |
661 |
|
|
!!$! + alphaPi*varEXP*riji) - (constERFC*rcuti2 & |
662 |
|
|
!!$! + alphaPi*constEXP*rcuti)) ) |
663 |
|
|
!!$ dudr = -sw*pre11*q_i*q_j * (varERFC*riji*riji & |
664 |
|
|
!!$ + alphaPi*varEXP*riji) |
665 |
|
|
!!$ |
666 |
|
|
!!$ dudx = dudx + dudr * xhat |
667 |
|
|
!!$ dudy = dudy + dudr * yhat |
668 |
|
|
!!$ dudz = dudz + dudr * zhat |
669 |
|
|
|
670 |
|
|
elseif (summationMethod .eq. SHIFTED_FORCE) then |
671 |
|
|
vterm = pre11 * q_i * q_j * (riji + rij*rcuti2 - 2.0d0*rcuti) |
672 |
chrisfen |
640 |
vpair = vpair + vterm |
673 |
|
|
epot = epot + sw*vterm |
674 |
|
|
|
675 |
chrisfen |
710 |
dudr = -sw*pre11*q_i*q_j * (riji*riji-rcuti2) |
676 |
|
|
|
677 |
chrisfen |
703 |
dudx = dudx + dudr * xhat |
678 |
|
|
dudy = dudy + dudr * yhat |
679 |
|
|
dudz = dudz + dudr * zhat |
680 |
chrisfen |
640 |
|
681 |
chrisfen |
710 |
!!$ elseif (summationMethod .eq. DAMPED_SHIFTED_FORCE) then |
682 |
|
|
!!$ varERFC = derfc(dampingAlpha*rij) |
683 |
|
|
!!$ varEXP = exp(-dampingAlpha*dampingAlpha*rij*rij) |
684 |
|
|
!!$ vterm = pre11 * q_i * q_j * (varERFC*riji - constERFC*rcuti) |
685 |
|
|
!!$ vpair = vpair + vterm |
686 |
|
|
!!$ epot = epot + sw*vterm |
687 |
|
|
!!$ |
688 |
|
|
!!$! dudr = -sw*pre11*q_i*q_j * (((varERFC*riji*riji & |
689 |
|
|
!!$! + alphaPi*varEXP*riji) - (constERFC*rcuti2 & |
690 |
|
|
!!$! + alphaPi*constEXP*rcuti)) ) |
691 |
|
|
!!$ dudr = -sw*pre11*q_i*q_j * (varERFC*riji*riji & |
692 |
|
|
!!$ + alphaPi*varEXP*riji) |
693 |
|
|
!!$ |
694 |
|
|
!!$ dudx = dudx + dudr * xhat |
695 |
|
|
!!$ dudy = dudy + dudr * yhat |
696 |
|
|
!!$ dudz = dudz + dudr * zhat |
697 |
|
|
|
698 |
chrisfen |
695 |
elseif (summationMethod .eq. REACTION_FIELD) then |
699 |
|
|
preVal = pre11 * q_i * q_j |
700 |
|
|
rfVal = preRF*rij*rij |
701 |
|
|
vterm = preVal * ( riji + rfVal ) |
702 |
chrisfen |
700 |
|
703 |
chrisfen |
695 |
vpair = vpair + vterm |
704 |
|
|
epot = epot + sw*vterm |
705 |
|
|
|
706 |
|
|
dudr = sw * preVal * ( 2.0d0*rfVal - riji )*riji |
707 |
|
|
|
708 |
|
|
dudx = dudx + dudr * xhat |
709 |
|
|
dudy = dudy + dudr * yhat |
710 |
|
|
dudz = dudz + dudr * zhat |
711 |
|
|
|
712 |
chrisfen |
597 |
else |
713 |
|
|
vterm = pre11 * q_i * q_j * riji |
714 |
|
|
vpair = vpair + vterm |
715 |
chrisfen |
626 |
epot = epot + sw*vterm |
716 |
chrisfen |
597 |
|
717 |
|
|
dudr = - sw * vterm * riji |
718 |
|
|
|
719 |
|
|
dudx = dudx + dudr * xhat |
720 |
|
|
dudy = dudy + dudr * yhat |
721 |
|
|
dudz = dudz + dudr * zhat |
722 |
|
|
|
723 |
|
|
endif |
724 |
|
|
|
725 |
gezelter |
411 |
endif |
726 |
|
|
|
727 |
|
|
if (j_is_Dipole) then |
728 |
|
|
|
729 |
chrisfen |
626 |
pref = pre12 * q_i * mu_j |
730 |
gezelter |
411 |
|
731 |
chrisfen |
710 |
!!$ if (summationMethod .eq. UNDAMPED_WOLF) then |
732 |
|
|
!!$ ri2 = riji * riji |
733 |
|
|
!!$ ri3 = ri2 * riji |
734 |
|
|
!!$ |
735 |
|
|
!!$ pref = pre12 * q_i * mu_j |
736 |
|
|
!!$ vterm = - pref * ct_j * (ri2 - rcuti2) |
737 |
|
|
!!$ vpair = vpair + vterm |
738 |
|
|
!!$ epot = epot + sw*vterm |
739 |
|
|
!!$ |
740 |
|
|
!!$ !! this has a + sign in the () because the rij vector is |
741 |
|
|
!!$ !! r_j - r_i and the charge-dipole potential takes the origin |
742 |
|
|
!!$ !! as the point dipole, which is atom j in this case. |
743 |
|
|
!!$ |
744 |
|
|
!!$ dudx = dudx - sw*pref * ( ri3*( uz_j(1) - 3.0d0*ct_j*xhat) & |
745 |
|
|
!!$ - rcuti3*( uz_j(1) - 3.0d0*ct_j*d(1)*rcuti ) ) |
746 |
|
|
!!$ dudy = dudy - sw*pref * ( ri3*( uz_j(2) - 3.0d0*ct_j*yhat) & |
747 |
|
|
!!$ - rcuti3*( uz_j(2) - 3.0d0*ct_j*d(2)*rcuti ) ) |
748 |
|
|
!!$ dudz = dudz - sw*pref * ( ri3*( uz_j(3) - 3.0d0*ct_j*zhat) & |
749 |
|
|
!!$ - rcuti3*( uz_j(3) - 3.0d0*ct_j*d(3)*rcuti ) ) |
750 |
|
|
!!$ |
751 |
|
|
!!$ duduz_j(1) = duduz_j(1) - sw*pref*( ri2*xhat - d(1)*rcuti3 ) |
752 |
|
|
!!$ duduz_j(2) = duduz_j(2) - sw*pref*( ri2*yhat - d(2)*rcuti3 ) |
753 |
|
|
!!$ duduz_j(3) = duduz_j(3) - sw*pref*( ri2*zhat - d(3)*rcuti3 ) |
754 |
|
|
!!$ |
755 |
|
|
!!$ elseif (summationMethod .eq. REACTION_FIELD) then |
756 |
gezelter |
507 |
|
757 |
chrisfen |
710 |
if (summationMethod .eq. REACTION_FIELD) then |
758 |
chrisfen |
700 |
ri2 = riji * riji |
759 |
|
|
ri3 = ri2 * riji |
760 |
chrisfen |
696 |
|
761 |
|
|
pref = pre12 * q_i * mu_j |
762 |
|
|
vterm = - pref * ct_j * ( ri2 - preRF2*rij ) |
763 |
|
|
vpair = vpair + vterm |
764 |
|
|
epot = epot + sw*vterm |
765 |
|
|
|
766 |
|
|
!! this has a + sign in the () because the rij vector is |
767 |
|
|
!! r_j - r_i and the charge-dipole potential takes the origin |
768 |
|
|
!! as the point dipole, which is atom j in this case. |
769 |
|
|
|
770 |
|
|
dudx = dudx - sw*pref*( ri3*(uz_j(1) - 3.0d0*ct_j*xhat) - & |
771 |
|
|
preRF2*uz_j(1) ) |
772 |
|
|
dudy = dudy - sw*pref*( ri3*(uz_j(2) - 3.0d0*ct_j*yhat) - & |
773 |
|
|
preRF2*uz_j(2) ) |
774 |
|
|
dudz = dudz - sw*pref*( ri3*(uz_j(3) - 3.0d0*ct_j*zhat) - & |
775 |
|
|
preRF2*uz_j(3) ) |
776 |
|
|
duduz_j(1) = duduz_j(1) - sw*pref * xhat * ( ri2 - preRF2*rij ) |
777 |
|
|
duduz_j(2) = duduz_j(2) - sw*pref * yhat * ( ri2 - preRF2*rij ) |
778 |
|
|
duduz_j(3) = duduz_j(3) - sw*pref * zhat * ( ri2 - preRF2*rij ) |
779 |
|
|
|
780 |
chrisfen |
597 |
else |
781 |
|
|
if (j_is_SplitDipole) then |
782 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_j * d_j) |
783 |
|
|
ri = 1.0_dp / BigR |
784 |
|
|
scale = rij * ri |
785 |
|
|
else |
786 |
|
|
ri = riji |
787 |
|
|
scale = 1.0_dp |
788 |
|
|
endif |
789 |
|
|
|
790 |
|
|
ri2 = ri * ri |
791 |
|
|
ri3 = ri2 * ri |
792 |
|
|
sc2 = scale * scale |
793 |
chrisfen |
626 |
|
794 |
|
|
pref = pre12 * q_i * mu_j |
795 |
chrisfen |
597 |
vterm = - pref * ct_j * ri2 * scale |
796 |
chrisfen |
626 |
vpair = vpair + vterm |
797 |
|
|
epot = epot + sw*vterm |
798 |
chrisfen |
597 |
|
799 |
|
|
!! this has a + sign in the () because the rij vector is |
800 |
|
|
!! r_j - r_i and the charge-dipole potential takes the origin |
801 |
|
|
!! as the point dipole, which is atom j in this case. |
802 |
|
|
|
803 |
chrisfen |
626 |
dudx = dudx - sw*pref * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2) |
804 |
|
|
dudy = dudy - sw*pref * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2) |
805 |
|
|
dudz = dudz - sw*pref * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2) |
806 |
chrisfen |
597 |
|
807 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) - sw*pref * ri2 * xhat * scale |
808 |
|
|
duduz_j(2) = duduz_j(2) - sw*pref * ri2 * yhat * scale |
809 |
|
|
duduz_j(3) = duduz_j(3) - sw*pref * ri2 * zhat * scale |
810 |
gezelter |
411 |
|
811 |
chrisfen |
597 |
endif |
812 |
gezelter |
411 |
endif |
813 |
gezelter |
421 |
|
814 |
gezelter |
439 |
if (j_is_Quadrupole) then |
815 |
|
|
ri2 = riji * riji |
816 |
|
|
ri3 = ri2 * riji |
817 |
gezelter |
440 |
ri4 = ri2 * ri2 |
818 |
gezelter |
439 |
cx2 = cx_j * cx_j |
819 |
|
|
cy2 = cy_j * cy_j |
820 |
|
|
cz2 = cz_j * cz_j |
821 |
|
|
|
822 |
chrisfen |
710 |
!!$ if (summationMethod .eq. UNDAMPED_WOLF) then |
823 |
|
|
!!$ pref = pre14 * q_i / 3.0_dp |
824 |
|
|
!!$ vterm1 = pref * ri3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + & |
825 |
|
|
!!$ qyy_j * (3.0_dp*cy2 - 1.0_dp) + & |
826 |
|
|
!!$ qzz_j * (3.0_dp*cz2 - 1.0_dp) ) |
827 |
|
|
!!$ vterm2 = pref * rcuti3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + & |
828 |
|
|
!!$ qyy_j * (3.0_dp*cy2 - 1.0_dp) + & |
829 |
|
|
!!$ qzz_j * (3.0_dp*cz2 - 1.0_dp) ) |
830 |
|
|
!!$ vpair = vpair + ( vterm1 - vterm2 ) |
831 |
|
|
!!$ epot = epot + sw*( vterm1 - vterm2 ) |
832 |
|
|
!!$ |
833 |
|
|
!!$ dudx = dudx - (5.0_dp * & |
834 |
|
|
!!$ (vterm1*riji*xhat - vterm2*rcuti2*d(1))) + sw*pref * ( & |
835 |
|
|
!!$ (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(1)) - & |
836 |
|
|
!!$ qxx_j*2.0_dp*(xhat - rcuti*d(1))) + & |
837 |
|
|
!!$ (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(1)) - & |
838 |
|
|
!!$ qyy_j*2.0_dp*(xhat - rcuti*d(1))) + & |
839 |
|
|
!!$ (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(1)) - & |
840 |
|
|
!!$ qzz_j*2.0_dp*(xhat - rcuti*d(1))) ) |
841 |
|
|
!!$ dudy = dudy - (5.0_dp * & |
842 |
|
|
!!$ (vterm1*riji*yhat - vterm2*rcuti2*d(2))) + sw*pref * ( & |
843 |
|
|
!!$ (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(2)) - & |
844 |
|
|
!!$ qxx_j*2.0_dp*(yhat - rcuti*d(2))) + & |
845 |
|
|
!!$ (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(2)) - & |
846 |
|
|
!!$ qyy_j*2.0_dp*(yhat - rcuti*d(2))) + & |
847 |
|
|
!!$ (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(2)) - & |
848 |
|
|
!!$ qzz_j*2.0_dp*(yhat - rcuti*d(2))) ) |
849 |
|
|
!!$ dudz = dudz - (5.0_dp * & |
850 |
|
|
!!$ (vterm1*riji*zhat - vterm2*rcuti2*d(3))) + sw*pref * ( & |
851 |
|
|
!!$ (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(3)) - & |
852 |
|
|
!!$ qxx_j*2.0_dp*(zhat - rcuti*d(3))) + & |
853 |
|
|
!!$ (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(3)) - & |
854 |
|
|
!!$ qyy_j*2.0_dp*(zhat - rcuti*d(3))) + & |
855 |
|
|
!!$ (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(3)) - & |
856 |
|
|
!!$ qzz_j*2.0_dp*(zhat - rcuti*d(3))) ) |
857 |
|
|
!!$ |
858 |
|
|
!!$ dudux_j(1) = dudux_j(1) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*xhat) -& |
859 |
|
|
!!$ rcuti4*(qxx_j*6.0_dp*cx_j*d(1))) |
860 |
|
|
!!$ dudux_j(2) = dudux_j(2) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*yhat) -& |
861 |
|
|
!!$ rcuti4*(qxx_j*6.0_dp*cx_j*d(2))) |
862 |
|
|
!!$ dudux_j(3) = dudux_j(3) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*zhat) -& |
863 |
|
|
!!$ rcuti4*(qxx_j*6.0_dp*cx_j*d(3))) |
864 |
|
|
!!$ |
865 |
|
|
!!$ duduy_j(1) = duduy_j(1) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*xhat) -& |
866 |
|
|
!!$ rcuti4*(qyy_j*6.0_dp*cx_j*d(1))) |
867 |
|
|
!!$ duduy_j(2) = duduy_j(2) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*yhat) -& |
868 |
|
|
!!$ rcuti4*(qyy_j*6.0_dp*cx_j*d(2))) |
869 |
|
|
!!$ duduy_j(3) = duduy_j(3) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*zhat) -& |
870 |
|
|
!!$ rcuti4*(qyy_j*6.0_dp*cx_j*d(3))) |
871 |
|
|
!!$ |
872 |
|
|
!!$ duduz_j(1) = duduz_j(1) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*xhat) -& |
873 |
|
|
!!$ rcuti4*(qzz_j*6.0_dp*cx_j*d(1))) |
874 |
|
|
!!$ duduz_j(2) = duduz_j(2) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*yhat) -& |
875 |
|
|
!!$ rcuti4*(qzz_j*6.0_dp*cx_j*d(2))) |
876 |
|
|
!!$ duduz_j(3) = duduz_j(3) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*zhat) -& |
877 |
|
|
!!$ rcuti4*(qzz_j*6.0_dp*cx_j*d(3))) |
878 |
|
|
!!$ |
879 |
|
|
!!$ else |
880 |
chrisfen |
626 |
pref = pre14 * q_i / 3.0_dp |
881 |
chrisfen |
597 |
vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + & |
882 |
|
|
qyy_j * (3.0_dp*cy2 - 1.0_dp) + & |
883 |
|
|
qzz_j * (3.0_dp*cz2 - 1.0_dp)) |
884 |
chrisfen |
626 |
vpair = vpair + vterm |
885 |
|
|
epot = epot + sw*vterm |
886 |
chrisfen |
597 |
|
887 |
chrisfen |
626 |
dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + sw*pref * ri4 * ( & |
888 |
chrisfen |
597 |
qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + & |
889 |
|
|
qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + & |
890 |
|
|
qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) ) |
891 |
chrisfen |
626 |
dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + sw*pref * ri4 * ( & |
892 |
chrisfen |
597 |
qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + & |
893 |
|
|
qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + & |
894 |
|
|
qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) ) |
895 |
chrisfen |
626 |
dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + sw*pref * ri4 * ( & |
896 |
chrisfen |
597 |
qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + & |
897 |
|
|
qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + & |
898 |
|
|
qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) ) |
899 |
|
|
|
900 |
chrisfen |
626 |
dudux_j(1) = dudux_j(1) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*xhat) |
901 |
|
|
dudux_j(2) = dudux_j(2) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*yhat) |
902 |
|
|
dudux_j(3) = dudux_j(3) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*zhat) |
903 |
chrisfen |
597 |
|
904 |
chrisfen |
626 |
duduy_j(1) = duduy_j(1) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*xhat) |
905 |
|
|
duduy_j(2) = duduy_j(2) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*yhat) |
906 |
|
|
duduy_j(3) = duduy_j(3) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*zhat) |
907 |
chrisfen |
597 |
|
908 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*xhat) |
909 |
|
|
duduz_j(2) = duduz_j(2) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*yhat) |
910 |
|
|
duduz_j(3) = duduz_j(3) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*zhat) |
911 |
chrisfen |
597 |
|
912 |
chrisfen |
710 |
!!$ endif |
913 |
gezelter |
439 |
endif |
914 |
gezelter |
411 |
endif |
915 |
gezelter |
507 |
|
916 |
gezelter |
411 |
if (i_is_Dipole) then |
917 |
gezelter |
507 |
|
918 |
gezelter |
411 |
if (j_is_Charge) then |
919 |
chrisfen |
626 |
|
920 |
|
|
pref = pre12 * q_j * mu_i |
921 |
|
|
|
922 |
chrisfen |
710 |
!!$ if (summationMethod .eq. UNDAMPED_WOLF) then |
923 |
|
|
!!$ ri2 = riji * riji |
924 |
|
|
!!$ ri3 = ri2 * riji |
925 |
|
|
!!$ |
926 |
|
|
!!$ pref = pre12 * q_j * mu_i |
927 |
|
|
!!$ vterm = pref * ct_i * (ri2 - rcuti2) |
928 |
|
|
!!$ vpair = vpair + vterm |
929 |
|
|
!!$ epot = epot + sw*vterm |
930 |
|
|
!!$ |
931 |
|
|
!!$ dudx = dudx + sw*pref * ( ri3*( uz_i(1) - 3.0d0*ct_i*xhat) & |
932 |
|
|
!!$ - rcuti3*( uz_i(1) - 3.0d0*ct_i*d(1)*rcuti ) ) |
933 |
|
|
!!$ dudy = dudy + sw*pref * ( ri3*( uz_i(2) - 3.0d0*ct_i*yhat) & |
934 |
|
|
!!$ - rcuti3*( uz_i(2) - 3.0d0*ct_i*d(2)*rcuti ) ) |
935 |
|
|
!!$ dudz = dudz + sw*pref * ( ri3*( uz_i(3) - 3.0d0*ct_i*zhat) & |
936 |
|
|
!!$ - rcuti3*( uz_i(3) - 3.0d0*ct_i*d(3)*rcuti ) ) |
937 |
|
|
!!$ |
938 |
|
|
!!$ duduz_i(1) = duduz_i(1) + sw*pref*( ri2*xhat - d(1)*rcuti3 ) |
939 |
|
|
!!$ duduz_i(2) = duduz_i(2) + sw*pref*( ri2*yhat - d(2)*rcuti3 ) |
940 |
|
|
!!$ duduz_i(3) = duduz_i(3) + sw*pref*( ri2*zhat - d(3)*rcuti3 ) |
941 |
|
|
!!$ |
942 |
|
|
!!$ elseif (summationMethod .eq. REACTION_FIELD) then |
943 |
|
|
if (summationMethod .eq. REACTION_FIELD) then |
944 |
chrisfen |
597 |
ri2 = riji * riji |
945 |
|
|
ri3 = ri2 * riji |
946 |
gezelter |
507 |
|
947 |
chrisfen |
626 |
pref = pre12 * q_j * mu_i |
948 |
chrisfen |
700 |
vterm = pref * ct_i * ( ri2 - preRF2*rij ) |
949 |
chrisfen |
696 |
vpair = vpair + vterm |
950 |
|
|
epot = epot + sw*vterm |
951 |
|
|
|
952 |
chrisfen |
700 |
dudx = dudx + sw*pref * ( ri3*(uz_i(1) - 3.0d0*ct_i*xhat) - & |
953 |
|
|
preRF2*uz_i(1) ) |
954 |
|
|
dudy = dudy + sw*pref * ( ri3*(uz_i(2) - 3.0d0*ct_i*yhat) - & |
955 |
|
|
preRF2*uz_i(2) ) |
956 |
|
|
dudz = dudz + sw*pref * ( ri3*(uz_i(3) - 3.0d0*ct_i*zhat) - & |
957 |
|
|
preRF2*uz_i(3) ) |
958 |
chrisfen |
696 |
|
959 |
chrisfen |
700 |
duduz_i(1) = duduz_i(1) + sw*pref * xhat * ( ri2 - preRF2*rij ) |
960 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref * yhat * ( ri2 - preRF2*rij ) |
961 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref * zhat * ( ri2 - preRF2*rij ) |
962 |
chrisfen |
696 |
|
963 |
chrisfen |
597 |
else |
964 |
|
|
if (i_is_SplitDipole) then |
965 |
gezelter |
421 |
BigR = sqrt(r2 + 0.25_dp * d_i * d_i) |
966 |
|
|
ri = 1.0_dp / BigR |
967 |
chrisfen |
597 |
scale = rij * ri |
968 |
|
|
else |
969 |
gezelter |
421 |
ri = riji |
970 |
|
|
scale = 1.0_dp |
971 |
|
|
endif |
972 |
chrisfen |
597 |
|
973 |
|
|
ri2 = ri * ri |
974 |
|
|
ri3 = ri2 * ri |
975 |
|
|
sc2 = scale * scale |
976 |
chrisfen |
626 |
|
977 |
|
|
pref = pre12 * q_j * mu_i |
978 |
chrisfen |
597 |
vterm = pref * ct_i * ri2 * scale |
979 |
chrisfen |
626 |
vpair = vpair + vterm |
980 |
|
|
epot = epot + sw*vterm |
981 |
chrisfen |
597 |
|
982 |
chrisfen |
626 |
dudx = dudx + sw*pref * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2) |
983 |
|
|
dudy = dudy + sw*pref * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2) |
984 |
|
|
dudz = dudz + sw*pref * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2) |
985 |
chrisfen |
597 |
|
986 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref * ri2 * xhat * scale |
987 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref * ri2 * yhat * scale |
988 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref * ri2 * zhat * scale |
989 |
gezelter |
421 |
endif |
990 |
chrisfen |
597 |
endif |
991 |
chrisfen |
626 |
|
992 |
chrisfen |
597 |
if (j_is_Dipole) then |
993 |
gezelter |
421 |
|
994 |
chrisfen |
710 |
!!$ if (summationMethod .eq. UNDAMPED_WOLF) then |
995 |
chrisfen |
703 |
!!$ ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3) |
996 |
|
|
!!$ |
997 |
|
|
!!$ ri2 = riji * riji |
998 |
|
|
!!$ ri3 = ri2 * riji |
999 |
|
|
!!$ ri4 = ri2 * ri2 |
1000 |
|
|
!!$ |
1001 |
|
|
!!$ pref = pre22 * mu_i * mu_j |
1002 |
|
|
!!$ vterm = pref * (ri3 - rcuti3) * (ct_ij - 3.0d0 * ct_i * ct_j) |
1003 |
|
|
!!$ vpair = vpair + vterm |
1004 |
|
|
!!$ epot = epot + sw*vterm |
1005 |
|
|
!!$ |
1006 |
|
|
!!$ a1 = 5.0d0 * ct_i * ct_j - ct_ij |
1007 |
|
|
!!$ |
1008 |
|
|
!!$ dudx = dudx + sw*pref*3.0d0*( & |
1009 |
|
|
!!$ ri4*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) & |
1010 |
|
|
!!$ - rcuti4*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) ) |
1011 |
|
|
!!$ dudy = dudy + sw*pref*3.0d0*( & |
1012 |
|
|
!!$ ri4*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) & |
1013 |
|
|
!!$ - rcuti4*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) ) |
1014 |
|
|
!!$ dudz = dudz + sw*pref*3.0d0*( & |
1015 |
|
|
!!$ ri4*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) & |
1016 |
|
|
!!$ - rcuti4*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) ) |
1017 |
|
|
!!$ |
1018 |
|
|
!!$ duduz_i(1) = duduz_i(1) + sw*pref*(ri3*(uz_j(1)-3.0d0*ct_j*xhat) & |
1019 |
|
|
!!$ - rcuti3*(uz_j(1) - 3.0d0*ct_j*xhat)) |
1020 |
|
|
!!$ duduz_i(2) = duduz_i(2) + sw*pref*(ri3*(uz_j(2)-3.0d0*ct_j*yhat) & |
1021 |
|
|
!!$ - rcuti3*(uz_j(2) - 3.0d0*ct_j*yhat)) |
1022 |
|
|
!!$ duduz_i(3) = duduz_i(3) + sw*pref*(ri3*(uz_j(3)-3.0d0*ct_j*zhat) & |
1023 |
|
|
!!$ - rcuti3*(uz_j(3) - 3.0d0*ct_j*zhat)) |
1024 |
|
|
!!$ duduz_j(1) = duduz_j(1) + sw*pref*(ri3*(uz_i(1)-3.0d0*ct_i*xhat) & |
1025 |
|
|
!!$ - rcuti3*(uz_i(1) - 3.0d0*ct_i*xhat)) |
1026 |
|
|
!!$ duduz_j(2) = duduz_j(2) + sw*pref*(ri3*(uz_i(2)-3.0d0*ct_i*yhat) & |
1027 |
|
|
!!$ - rcuti3*(uz_i(2) - 3.0d0*ct_i*yhat)) |
1028 |
|
|
!!$ duduz_j(3) = duduz_j(3) + sw*pref*(ri3*(uz_i(3)-3.0d0*ct_i*zhat) & |
1029 |
|
|
!!$ - rcuti3*(uz_i(3) - 3.0d0*ct_i*zhat)) |
1030 |
chrisfen |
710 |
!!$ |
1031 |
|
|
!!$ elseif (summationMethod .eq. DAMPED_WOLF) then |
1032 |
|
|
!!$ ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3) |
1033 |
|
|
!!$ |
1034 |
|
|
!!$ ri2 = riji * riji |
1035 |
|
|
!!$ ri3 = ri2 * riji |
1036 |
|
|
!!$ ri4 = ri2 * ri2 |
1037 |
|
|
!!$ sc2 = scale * scale |
1038 |
|
|
!!$ |
1039 |
|
|
!!$ pref = pre22 * mu_i * mu_j |
1040 |
|
|
!!$ vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j) |
1041 |
|
|
!!$ vpair = vpair + vterm |
1042 |
|
|
!!$ epot = epot + sw*vterm |
1043 |
|
|
!!$ |
1044 |
|
|
!!$ a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij |
1045 |
|
|
!!$ |
1046 |
|
|
!!$ dudx = dudx + sw*pref*3.0d0*ri4*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) |
1047 |
|
|
!!$ dudy = dudy + sw*pref*3.0d0*ri4*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) |
1048 |
|
|
!!$ dudz = dudz + sw*pref*3.0d0*ri4*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) |
1049 |
|
|
!!$ |
1050 |
|
|
!!$ duduz_i(1) = duduz_i(1) + sw*pref*ri3 *(uz_j(1) - 3.0d0*ct_j*xhat) |
1051 |
|
|
!!$ duduz_i(2) = duduz_i(2) + sw*pref*ri3 *(uz_j(2) - 3.0d0*ct_j*yhat) |
1052 |
|
|
!!$ duduz_i(3) = duduz_i(3) + sw*pref*ri3 *(uz_j(3) - 3.0d0*ct_j*zhat) |
1053 |
|
|
!!$ |
1054 |
|
|
!!$ duduz_j(1) = duduz_j(1) + sw*pref*ri3 *(uz_i(1) - 3.0d0*ct_i*xhat) |
1055 |
|
|
!!$ duduz_j(2) = duduz_j(2) + sw*pref*ri3 *(uz_i(2) - 3.0d0*ct_i*yhat) |
1056 |
|
|
!!$ duduz_j(3) = duduz_j(3) + sw*pref*ri3 *(uz_i(3) - 3.0d0*ct_i*zhat) |
1057 |
|
|
!!$ |
1058 |
|
|
!!$ elseif (summationMethod .eq. REACTION_FIELD) then |
1059 |
|
|
if (summationMethod .eq. REACTION_FIELD) then |
1060 |
chrisfen |
703 |
ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3) |
1061 |
chrisfen |
695 |
|
1062 |
|
|
ri2 = riji * riji |
1063 |
|
|
ri3 = ri2 * riji |
1064 |
|
|
ri4 = ri2 * ri2 |
1065 |
|
|
|
1066 |
|
|
pref = pre22 * mu_i * mu_j |
1067 |
|
|
|
1068 |
|
|
vterm = pref*( ri3*(ct_ij - 3.0d0 * ct_i * ct_j) - & |
1069 |
|
|
preRF2*ct_ij ) |
1070 |
|
|
vpair = vpair + vterm |
1071 |
|
|
epot = epot + sw*vterm |
1072 |
|
|
|
1073 |
|
|
a1 = 5.0d0 * ct_i * ct_j - ct_ij |
1074 |
|
|
|
1075 |
|
|
dudx = dudx + sw*pref*3.0d0*ri4 & |
1076 |
|
|
* (a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) |
1077 |
|
|
dudy = dudy + sw*pref*3.0d0*ri4 & |
1078 |
|
|
* (a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) |
1079 |
|
|
dudz = dudz + sw*pref*3.0d0*ri4 & |
1080 |
|
|
* (a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) |
1081 |
|
|
|
1082 |
|
|
duduz_i(1) = duduz_i(1) + sw*pref*(ri3*(uz_j(1)-3.0d0*ct_j*xhat) & |
1083 |
|
|
- preRF2*uz_j(1)) |
1084 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref*(ri3*(uz_j(2)-3.0d0*ct_j*yhat) & |
1085 |
|
|
- preRF2*uz_j(2)) |
1086 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref*(ri3*(uz_j(3)-3.0d0*ct_j*zhat) & |
1087 |
|
|
- preRF2*uz_j(3)) |
1088 |
|
|
duduz_j(1) = duduz_j(1) + sw*pref*(ri3*(uz_i(1)-3.0d0*ct_i*xhat) & |
1089 |
|
|
- preRF2*uz_i(1)) |
1090 |
|
|
duduz_j(2) = duduz_j(2) + sw*pref*(ri3*(uz_i(2)-3.0d0*ct_i*yhat) & |
1091 |
|
|
- preRF2*uz_i(2)) |
1092 |
|
|
duduz_j(3) = duduz_j(3) + sw*pref*(ri3*(uz_i(3)-3.0d0*ct_i*zhat) & |
1093 |
|
|
- preRF2*uz_i(3)) |
1094 |
|
|
|
1095 |
chrisfen |
597 |
else |
1096 |
|
|
if (i_is_SplitDipole) then |
1097 |
|
|
if (j_is_SplitDipole) then |
1098 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j) |
1099 |
|
|
else |
1100 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_i * d_i) |
1101 |
|
|
endif |
1102 |
|
|
ri = 1.0_dp / BigR |
1103 |
|
|
scale = rij * ri |
1104 |
|
|
else |
1105 |
|
|
if (j_is_SplitDipole) then |
1106 |
|
|
BigR = sqrt(r2 + 0.25_dp * d_j * d_j) |
1107 |
|
|
ri = 1.0_dp / BigR |
1108 |
|
|
scale = rij * ri |
1109 |
|
|
else |
1110 |
|
|
ri = riji |
1111 |
|
|
scale = 1.0_dp |
1112 |
|
|
endif |
1113 |
|
|
endif |
1114 |
|
|
|
1115 |
|
|
ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3) |
1116 |
|
|
|
1117 |
|
|
ri2 = ri * ri |
1118 |
|
|
ri3 = ri2 * ri |
1119 |
|
|
ri4 = ri2 * ri2 |
1120 |
|
|
sc2 = scale * scale |
1121 |
|
|
|
1122 |
chrisfen |
626 |
pref = pre22 * mu_i * mu_j |
1123 |
chrisfen |
597 |
vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2) |
1124 |
chrisfen |
626 |
vpair = vpair + vterm |
1125 |
|
|
epot = epot + sw*vterm |
1126 |
chrisfen |
597 |
|
1127 |
|
|
a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij |
1128 |
|
|
|
1129 |
chrisfen |
626 |
dudx = dudx + sw*pref*3.0d0*ri4*scale & |
1130 |
|
|
*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) |
1131 |
|
|
dudy = dudy + sw*pref*3.0d0*ri4*scale & |
1132 |
|
|
*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) |
1133 |
|
|
dudz = dudz + sw*pref*3.0d0*ri4*scale & |
1134 |
|
|
*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) |
1135 |
chrisfen |
597 |
|
1136 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref*ri3 & |
1137 |
|
|
*(uz_j(1) - 3.0d0*ct_j*xhat*sc2) |
1138 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref*ri3 & |
1139 |
|
|
*(uz_j(2) - 3.0d0*ct_j*yhat*sc2) |
1140 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref*ri3 & |
1141 |
|
|
*(uz_j(3) - 3.0d0*ct_j*zhat*sc2) |
1142 |
chrisfen |
597 |
|
1143 |
chrisfen |
626 |
duduz_j(1) = duduz_j(1) + sw*pref*ri3 & |
1144 |
|
|
*(uz_i(1) - 3.0d0*ct_i*xhat*sc2) |
1145 |
|
|
duduz_j(2) = duduz_j(2) + sw*pref*ri3 & |
1146 |
|
|
*(uz_i(2) - 3.0d0*ct_i*yhat*sc2) |
1147 |
|
|
duduz_j(3) = duduz_j(3) + sw*pref*ri3 & |
1148 |
|
|
*(uz_i(3) - 3.0d0*ct_i*zhat*sc2) |
1149 |
chrisfen |
597 |
endif |
1150 |
gezelter |
411 |
endif |
1151 |
|
|
endif |
1152 |
gezelter |
439 |
|
1153 |
|
|
if (i_is_Quadrupole) then |
1154 |
|
|
if (j_is_Charge) then |
1155 |
gezelter |
507 |
|
1156 |
gezelter |
439 |
ri2 = riji * riji |
1157 |
|
|
ri3 = ri2 * riji |
1158 |
gezelter |
440 |
ri4 = ri2 * ri2 |
1159 |
gezelter |
439 |
cx2 = cx_i * cx_i |
1160 |
|
|
cy2 = cy_i * cy_i |
1161 |
|
|
cz2 = cz_i * cz_i |
1162 |
gezelter |
507 |
|
1163 |
chrisfen |
710 |
!!$ if (summationMethod .eq. UNDAMPED_WOLF) then |
1164 |
|
|
!!$ pref = pre14 * q_j / 3.0_dp |
1165 |
|
|
!!$ vterm1 = pref * ri3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + & |
1166 |
|
|
!!$ qyy_i * (3.0_dp*cy2 - 1.0_dp) + & |
1167 |
|
|
!!$ qzz_i * (3.0_dp*cz2 - 1.0_dp) ) |
1168 |
|
|
!!$ vterm2 = pref * rcuti3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + & |
1169 |
|
|
!!$ qyy_i * (3.0_dp*cy2 - 1.0_dp) + & |
1170 |
|
|
!!$ qzz_i * (3.0_dp*cz2 - 1.0_dp) ) |
1171 |
|
|
!!$ vpair = vpair + ( vterm1 - vterm2 ) |
1172 |
|
|
!!$ epot = epot + sw*( vterm1 - vterm2 ) |
1173 |
|
|
!!$ |
1174 |
|
|
!!$ dudx = dudx - sw*(5.0_dp*(vterm1*riji*xhat-vterm2*rcuti2*d(1))) +& |
1175 |
|
|
!!$ sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(1)) - & |
1176 |
|
|
!!$ qxx_i*2.0_dp*(xhat - rcuti*d(1))) + & |
1177 |
|
|
!!$ (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(1)) - & |
1178 |
|
|
!!$ qyy_i*2.0_dp*(xhat - rcuti*d(1))) + & |
1179 |
|
|
!!$ (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(1)) - & |
1180 |
|
|
!!$ qzz_i*2.0_dp*(xhat - rcuti*d(1))) ) |
1181 |
|
|
!!$ dudy = dudy - sw*(5.0_dp*(vterm1*riji*yhat-vterm2*rcuti2*d(2))) +& |
1182 |
|
|
!!$ sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(2)) - & |
1183 |
|
|
!!$ qxx_i*2.0_dp*(yhat - rcuti*d(2))) + & |
1184 |
|
|
!!$ (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(2)) - & |
1185 |
|
|
!!$ qyy_i*2.0_dp*(yhat - rcuti*d(2))) + & |
1186 |
|
|
!!$ (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(2)) - & |
1187 |
|
|
!!$ qzz_i*2.0_dp*(yhat - rcuti*d(2))) ) |
1188 |
|
|
!!$ dudz = dudz - sw*(5.0_dp*(vterm1*riji*zhat-vterm2*rcuti2*d(3))) +& |
1189 |
|
|
!!$ sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(3)) - & |
1190 |
|
|
!!$ qxx_i*2.0_dp*(zhat - rcuti*d(3))) + & |
1191 |
|
|
!!$ (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(3)) - & |
1192 |
|
|
!!$ qyy_i*2.0_dp*(zhat - rcuti*d(3))) + & |
1193 |
|
|
!!$ (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(3)) - & |
1194 |
|
|
!!$ qzz_i*2.0_dp*(zhat - rcuti*d(3))) ) |
1195 |
|
|
!!$ |
1196 |
|
|
!!$ dudux_i(1) = dudux_i(1) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*xhat) -& |
1197 |
|
|
!!$ rcuti4*(qxx_i*6.0_dp*cx_i*d(1))) |
1198 |
|
|
!!$ dudux_i(2) = dudux_i(2) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*yhat) -& |
1199 |
|
|
!!$ rcuti4*(qxx_i*6.0_dp*cx_i*d(2))) |
1200 |
|
|
!!$ dudux_i(3) = dudux_i(3) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*zhat) -& |
1201 |
|
|
!!$ rcuti4*(qxx_i*6.0_dp*cx_i*d(3))) |
1202 |
|
|
!!$ |
1203 |
|
|
!!$ duduy_i(1) = duduy_i(1) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*xhat) -& |
1204 |
|
|
!!$ rcuti4*(qyy_i*6.0_dp*cx_i*d(1))) |
1205 |
|
|
!!$ duduy_i(2) = duduy_i(2) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*yhat) -& |
1206 |
|
|
!!$ rcuti4*(qyy_i*6.0_dp*cx_i*d(2))) |
1207 |
|
|
!!$ duduy_i(3) = duduy_i(3) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*zhat) -& |
1208 |
|
|
!!$ rcuti4*(qyy_i*6.0_dp*cx_i*d(3))) |
1209 |
|
|
!!$ |
1210 |
|
|
!!$ duduz_i(1) = duduz_i(1) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*xhat) -& |
1211 |
|
|
!!$ rcuti4*(qzz_i*6.0_dp*cx_i*d(1))) |
1212 |
|
|
!!$ duduz_i(2) = duduz_i(2) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*yhat) -& |
1213 |
|
|
!!$ rcuti4*(qzz_i*6.0_dp*cx_i*d(2))) |
1214 |
|
|
!!$ duduz_i(3) = duduz_i(3) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*zhat) -& |
1215 |
|
|
!!$ rcuti4*(qzz_i*6.0_dp*cx_i*d(3))) |
1216 |
|
|
!!$ |
1217 |
|
|
!!$ else |
1218 |
chrisfen |
626 |
pref = pre14 * q_j / 3.0_dp |
1219 |
chrisfen |
597 |
vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + & |
1220 |
|
|
qyy_i * (3.0_dp*cy2 - 1.0_dp) + & |
1221 |
|
|
qzz_i * (3.0_dp*cz2 - 1.0_dp)) |
1222 |
chrisfen |
626 |
vpair = vpair + vterm |
1223 |
|
|
epot = epot + sw*vterm |
1224 |
chrisfen |
597 |
|
1225 |
chrisfen |
626 |
dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + sw*pref*ri4 * ( & |
1226 |
chrisfen |
597 |
qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + & |
1227 |
|
|
qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + & |
1228 |
|
|
qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) ) |
1229 |
chrisfen |
626 |
dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + sw*pref*ri4 * ( & |
1230 |
chrisfen |
597 |
qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + & |
1231 |
|
|
qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + & |
1232 |
|
|
qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) ) |
1233 |
chrisfen |
626 |
dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + sw*pref*ri4 * ( & |
1234 |
chrisfen |
597 |
qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + & |
1235 |
|
|
qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + & |
1236 |
|
|
qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) ) |
1237 |
|
|
|
1238 |
chrisfen |
626 |
dudux_i(1) = dudux_i(1) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*xhat) |
1239 |
|
|
dudux_i(2) = dudux_i(2) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*yhat) |
1240 |
|
|
dudux_i(3) = dudux_i(3) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*zhat) |
1241 |
chrisfen |
597 |
|
1242 |
chrisfen |
626 |
duduy_i(1) = duduy_i(1) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*xhat) |
1243 |
|
|
duduy_i(2) = duduy_i(2) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*yhat) |
1244 |
|
|
duduy_i(3) = duduy_i(3) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*zhat) |
1245 |
chrisfen |
597 |
|
1246 |
chrisfen |
626 |
duduz_i(1) = duduz_i(1) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*xhat) |
1247 |
|
|
duduz_i(2) = duduz_i(2) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*yhat) |
1248 |
|
|
duduz_i(3) = duduz_i(3) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*zhat) |
1249 |
chrisfen |
710 |
!!$ endif |
1250 |
gezelter |
439 |
endif |
1251 |
|
|
endif |
1252 |
gezelter |
507 |
|
1253 |
|
|
|
1254 |
gezelter |
411 |
if (do_pot) then |
1255 |
|
|
#ifdef IS_MPI |
1256 |
chuckv |
656 |
pot_row(ELECTROSTATIC_POT,atom1) = pot_row(ELECTROSTATIC_POT,atom1) + 0.5d0*epot |
1257 |
|
|
pot_col(ELECTROSTATIC_POT,atom2) = pot_col(ELECTROSTATIC_POT,atom2) + 0.5d0*epot |
1258 |
gezelter |
411 |
#else |
1259 |
|
|
pot = pot + epot |
1260 |
|
|
#endif |
1261 |
|
|
endif |
1262 |
gezelter |
507 |
|
1263 |
gezelter |
411 |
#ifdef IS_MPI |
1264 |
|
|
f_Row(1,atom1) = f_Row(1,atom1) + dudx |
1265 |
|
|
f_Row(2,atom1) = f_Row(2,atom1) + dudy |
1266 |
|
|
f_Row(3,atom1) = f_Row(3,atom1) + dudz |
1267 |
gezelter |
507 |
|
1268 |
gezelter |
411 |
f_Col(1,atom2) = f_Col(1,atom2) - dudx |
1269 |
|
|
f_Col(2,atom2) = f_Col(2,atom2) - dudy |
1270 |
|
|
f_Col(3,atom2) = f_Col(3,atom2) - dudz |
1271 |
gezelter |
507 |
|
1272 |
gezelter |
411 |
if (i_is_Dipole .or. i_is_Quadrupole) then |
1273 |
gezelter |
439 |
t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2) |
1274 |
|
|
t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3) |
1275 |
|
|
t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1) |
1276 |
gezelter |
411 |
endif |
1277 |
gezelter |
439 |
if (i_is_Quadrupole) then |
1278 |
|
|
t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2) |
1279 |
|
|
t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3) |
1280 |
|
|
t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1) |
1281 |
gezelter |
411 |
|
1282 |
gezelter |
439 |
t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2) |
1283 |
|
|
t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3) |
1284 |
|
|
t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1) |
1285 |
|
|
endif |
1286 |
|
|
|
1287 |
gezelter |
411 |
if (j_is_Dipole .or. j_is_Quadrupole) then |
1288 |
gezelter |
439 |
t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2) |
1289 |
|
|
t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3) |
1290 |
|
|
t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1) |
1291 |
gezelter |
411 |
endif |
1292 |
gezelter |
439 |
if (j_is_Quadrupole) then |
1293 |
|
|
t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2) |
1294 |
|
|
t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3) |
1295 |
|
|
t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1) |
1296 |
gezelter |
411 |
|
1297 |
gezelter |
439 |
t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2) |
1298 |
|
|
t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3) |
1299 |
|
|
t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1) |
1300 |
|
|
endif |
1301 |
|
|
|
1302 |
gezelter |
411 |
#else |
1303 |
|
|
f(1,atom1) = f(1,atom1) + dudx |
1304 |
|
|
f(2,atom1) = f(2,atom1) + dudy |
1305 |
|
|
f(3,atom1) = f(3,atom1) + dudz |
1306 |
gezelter |
507 |
|
1307 |
gezelter |
411 |
f(1,atom2) = f(1,atom2) - dudx |
1308 |
|
|
f(2,atom2) = f(2,atom2) - dudy |
1309 |
|
|
f(3,atom2) = f(3,atom2) - dudz |
1310 |
gezelter |
507 |
|
1311 |
gezelter |
411 |
if (i_is_Dipole .or. i_is_Quadrupole) then |
1312 |
gezelter |
439 |
t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2) |
1313 |
|
|
t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3) |
1314 |
|
|
t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1) |
1315 |
gezelter |
411 |
endif |
1316 |
gezelter |
439 |
if (i_is_Quadrupole) then |
1317 |
|
|
t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2) |
1318 |
|
|
t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3) |
1319 |
|
|
t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1) |
1320 |
|
|
|
1321 |
|
|
t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2) |
1322 |
|
|
t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3) |
1323 |
|
|
t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1) |
1324 |
|
|
endif |
1325 |
|
|
|
1326 |
gezelter |
411 |
if (j_is_Dipole .or. j_is_Quadrupole) then |
1327 |
gezelter |
439 |
t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2) |
1328 |
|
|
t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3) |
1329 |
|
|
t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1) |
1330 |
gezelter |
411 |
endif |
1331 |
gezelter |
439 |
if (j_is_Quadrupole) then |
1332 |
|
|
t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2) |
1333 |
|
|
t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3) |
1334 |
|
|
t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1) |
1335 |
|
|
|
1336 |
|
|
t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2) |
1337 |
|
|
t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3) |
1338 |
|
|
t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1) |
1339 |
|
|
endif |
1340 |
|
|
|
1341 |
gezelter |
411 |
#endif |
1342 |
gezelter |
507 |
|
1343 |
gezelter |
411 |
#ifdef IS_MPI |
1344 |
|
|
id1 = AtomRowToGlobal(atom1) |
1345 |
|
|
id2 = AtomColToGlobal(atom2) |
1346 |
|
|
#else |
1347 |
|
|
id1 = atom1 |
1348 |
|
|
id2 = atom2 |
1349 |
|
|
#endif |
1350 |
|
|
|
1351 |
|
|
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
1352 |
gezelter |
507 |
|
1353 |
gezelter |
411 |
fpair(1) = fpair(1) + dudx |
1354 |
|
|
fpair(2) = fpair(2) + dudy |
1355 |
|
|
fpair(3) = fpair(3) + dudz |
1356 |
|
|
|
1357 |
|
|
endif |
1358 |
|
|
|
1359 |
|
|
return |
1360 |
|
|
end subroutine doElectrostaticPair |
1361 |
chuckv |
492 |
|
1362 |
|
|
subroutine destroyElectrostaticTypes() |
1363 |
|
|
|
1364 |
gezelter |
507 |
if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap) |
1365 |
|
|
|
1366 |
chuckv |
492 |
end subroutine destroyElectrostaticTypes |
1367 |
|
|
|
1368 |
chrisfen |
703 |
subroutine self_self(atom1, eFrame, mypot, t, do_pot) |
1369 |
chrisfen |
695 |
logical, intent(in) :: do_pot |
1370 |
chrisfen |
682 |
integer, intent(in) :: atom1 |
1371 |
chrisfen |
695 |
integer :: atid1 |
1372 |
chrisfen |
682 |
real(kind=dp), dimension(9,nLocal) :: eFrame |
1373 |
chrisfen |
695 |
real(kind=dp), dimension(3,nLocal) :: t |
1374 |
chrisfen |
703 |
real(kind=dp) :: mu1, c1 |
1375 |
|
|
real(kind=dp) :: preVal, epot, mypot |
1376 |
chrisfen |
695 |
real(kind=dp) :: eix, eiy, eiz |
1377 |
chrisfen |
682 |
|
1378 |
chrisfen |
695 |
! this is a local only array, so we use the local atom type id's: |
1379 |
|
|
atid1 = atid(atom1) |
1380 |
chrisfen |
703 |
|
1381 |
|
|
if (.not.summationMethodChecked) then |
1382 |
|
|
call checkSummationMethod() |
1383 |
|
|
endif |
1384 |
chrisfen |
695 |
|
1385 |
chrisfen |
703 |
if (summationMethod .eq. REACTION_FIELD) then |
1386 |
|
|
if (ElectrostaticMap(atid1)%is_Dipole) then |
1387 |
|
|
mu1 = getDipoleMoment(atid1) |
1388 |
|
|
|
1389 |
|
|
preVal = pre22 * preRF2 * mu1*mu1 |
1390 |
|
|
mypot = mypot - 0.5d0*preVal |
1391 |
|
|
|
1392 |
|
|
! The self-correction term adds into the reaction field vector |
1393 |
|
|
|
1394 |
|
|
eix = preVal * eFrame(3,atom1) |
1395 |
|
|
eiy = preVal * eFrame(6,atom1) |
1396 |
|
|
eiz = preVal * eFrame(9,atom1) |
1397 |
|
|
|
1398 |
|
|
! once again, this is self-self, so only the local arrays are needed |
1399 |
|
|
! even for MPI jobs: |
1400 |
|
|
|
1401 |
|
|
t(1,atom1)=t(1,atom1) - eFrame(6,atom1)*eiz + & |
1402 |
|
|
eFrame(9,atom1)*eiy |
1403 |
|
|
t(2,atom1)=t(2,atom1) - eFrame(9,atom1)*eix + & |
1404 |
|
|
eFrame(3,atom1)*eiz |
1405 |
|
|
t(3,atom1)=t(3,atom1) - eFrame(3,atom1)*eiy + & |
1406 |
|
|
eFrame(6,atom1)*eix |
1407 |
|
|
|
1408 |
|
|
endif |
1409 |
|
|
|
1410 |
chrisfen |
710 |
!!$ elseif (summationMethod .eq. UNDAMPED_WOLF) then |
1411 |
|
|
!!$ if (ElectrostaticMap(atid1)%is_Charge) then |
1412 |
|
|
!!$ c1 = getCharge(atid1) |
1413 |
|
|
!!$ |
1414 |
|
|
!!$ mypot = mypot - (rcuti * 0.5_dp * c1 * c1) |
1415 |
|
|
!!$ endif |
1416 |
|
|
!!$ |
1417 |
|
|
!!$ elseif (summationMethod .eq. DAMPED_WOLF) then |
1418 |
|
|
!!$ if (ElectrostaticMap(atid1)%is_Charge) then |
1419 |
|
|
!!$ c1 = getCharge(atid1) |
1420 |
|
|
!!$ |
1421 |
|
|
!!$ mypot = mypot - (constERFC * rcuti * 0.5_dp + & |
1422 |
|
|
!!$ dampingAlpha*invRootPi) * c1 * c1 |
1423 |
|
|
!!$ endif |
1424 |
chrisfen |
682 |
endif |
1425 |
chrisfen |
695 |
|
1426 |
chrisfen |
682 |
return |
1427 |
chrisfen |
703 |
end subroutine self_self |
1428 |
chrisfen |
682 |
|
1429 |
chrisfen |
703 |
subroutine rf_self_excludes(atom1, atom2, sw, eFrame, d, rij, vpair, myPot, & |
1430 |
chrisfen |
700 |
f, t, do_pot) |
1431 |
|
|
logical, intent(in) :: do_pot |
1432 |
|
|
integer, intent(in) :: atom1 |
1433 |
|
|
integer, intent(in) :: atom2 |
1434 |
|
|
logical :: i_is_Charge, j_is_Charge |
1435 |
|
|
logical :: i_is_Dipole, j_is_Dipole |
1436 |
|
|
integer :: atid1 |
1437 |
|
|
integer :: atid2 |
1438 |
|
|
real(kind=dp), intent(in) :: rij |
1439 |
|
|
real(kind=dp), intent(in) :: sw |
1440 |
|
|
real(kind=dp), intent(in), dimension(3) :: d |
1441 |
|
|
real(kind=dp), intent(inout) :: vpair |
1442 |
|
|
real(kind=dp), dimension(9,nLocal) :: eFrame |
1443 |
|
|
real(kind=dp), dimension(3,nLocal) :: f |
1444 |
|
|
real(kind=dp), dimension(3,nLocal) :: t |
1445 |
|
|
real (kind = dp), dimension(3) :: duduz_i |
1446 |
|
|
real (kind = dp), dimension(3) :: duduz_j |
1447 |
|
|
real (kind = dp), dimension(3) :: uz_i |
1448 |
|
|
real (kind = dp), dimension(3) :: uz_j |
1449 |
|
|
real(kind=dp) :: q_i, q_j, mu_i, mu_j |
1450 |
|
|
real(kind=dp) :: xhat, yhat, zhat |
1451 |
|
|
real(kind=dp) :: ct_i, ct_j |
1452 |
|
|
real(kind=dp) :: ri2, ri3, riji, vterm |
1453 |
chrisfen |
703 |
real(kind=dp) :: pref, preVal, rfVal, myPot |
1454 |
chrisfen |
700 |
real(kind=dp) :: dudx, dudy, dudz, dudr |
1455 |
|
|
|
1456 |
chrisfen |
703 |
if (.not.summationMethodChecked) then |
1457 |
|
|
call checkSummationMethod() |
1458 |
chrisfen |
700 |
endif |
1459 |
|
|
|
1460 |
|
|
dudx = 0.0d0 |
1461 |
|
|
dudy = 0.0d0 |
1462 |
|
|
dudz = 0.0d0 |
1463 |
|
|
|
1464 |
|
|
riji = 1.0d0/rij |
1465 |
|
|
|
1466 |
|
|
xhat = d(1) * riji |
1467 |
|
|
yhat = d(2) * riji |
1468 |
|
|
zhat = d(3) * riji |
1469 |
|
|
|
1470 |
|
|
! this is a local only array, so we use the local atom type id's: |
1471 |
|
|
atid1 = atid(atom1) |
1472 |
|
|
atid2 = atid(atom2) |
1473 |
|
|
i_is_Charge = ElectrostaticMap(atid1)%is_Charge |
1474 |
|
|
j_is_Charge = ElectrostaticMap(atid2)%is_Charge |
1475 |
|
|
i_is_Dipole = ElectrostaticMap(atid1)%is_Dipole |
1476 |
|
|
j_is_Dipole = ElectrostaticMap(atid2)%is_Dipole |
1477 |
|
|
|
1478 |
|
|
if (i_is_Charge.and.j_is_Charge) then |
1479 |
|
|
q_i = ElectrostaticMap(atid1)%charge |
1480 |
|
|
q_j = ElectrostaticMap(atid2)%charge |
1481 |
|
|
|
1482 |
|
|
preVal = pre11 * q_i * q_j |
1483 |
|
|
rfVal = preRF*rij*rij |
1484 |
|
|
vterm = preVal * rfVal |
1485 |
|
|
|
1486 |
chrisfen |
703 |
myPot = myPot + sw*vterm |
1487 |
|
|
|
1488 |
chrisfen |
700 |
dudr = sw*preVal * 2.0d0*rfVal*riji |
1489 |
chrisfen |
703 |
|
1490 |
chrisfen |
700 |
dudx = dudx + dudr * xhat |
1491 |
|
|
dudy = dudy + dudr * yhat |
1492 |
|
|
dudz = dudz + dudr * zhat |
1493 |
chrisfen |
703 |
|
1494 |
chrisfen |
700 |
elseif (i_is_Charge.and.j_is_Dipole) then |
1495 |
|
|
q_i = ElectrostaticMap(atid1)%charge |
1496 |
|
|
mu_j = ElectrostaticMap(atid2)%dipole_moment |
1497 |
|
|
uz_j(1) = eFrame(3,atom2) |
1498 |
|
|
uz_j(2) = eFrame(6,atom2) |
1499 |
|
|
uz_j(3) = eFrame(9,atom2) |
1500 |
|
|
ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat |
1501 |
chrisfen |
703 |
|
1502 |
chrisfen |
700 |
ri2 = riji * riji |
1503 |
|
|
ri3 = ri2 * riji |
1504 |
|
|
|
1505 |
|
|
pref = pre12 * q_i * mu_j |
1506 |
|
|
vterm = - pref * ct_j * ( ri2 - preRF2*rij ) |
1507 |
chrisfen |
703 |
myPot = myPot + sw*vterm |
1508 |
|
|
|
1509 |
|
|
dudx = dudx - sw*pref*( ri3*(uz_j(1)-3.0d0*ct_j*xhat) & |
1510 |
|
|
- preRF2*uz_j(1) ) |
1511 |
|
|
dudy = dudy - sw*pref*( ri3*(uz_j(2)-3.0d0*ct_j*yhat) & |
1512 |
|
|
- preRF2*uz_j(2) ) |
1513 |
|
|
dudz = dudz - sw*pref*( ri3*(uz_j(3)-3.0d0*ct_j*zhat) & |
1514 |
|
|
- preRF2*uz_j(3) ) |
1515 |
|
|
|
1516 |
chrisfen |
700 |
duduz_j(1) = duduz_j(1) - sw * pref * xhat * ( ri2 - preRF2*rij ) |
1517 |
|
|
duduz_j(2) = duduz_j(2) - sw * pref * yhat * ( ri2 - preRF2*rij ) |
1518 |
|
|
duduz_j(3) = duduz_j(3) - sw * pref * zhat * ( ri2 - preRF2*rij ) |
1519 |
chrisfen |
703 |
|
1520 |
chrisfen |
700 |
elseif (i_is_Dipole.and.j_is_Charge) then |
1521 |
|
|
mu_i = ElectrostaticMap(atid1)%dipole_moment |
1522 |
|
|
q_j = ElectrostaticMap(atid2)%charge |
1523 |
|
|
uz_i(1) = eFrame(3,atom1) |
1524 |
|
|
uz_i(2) = eFrame(6,atom1) |
1525 |
|
|
uz_i(3) = eFrame(9,atom1) |
1526 |
|
|
ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat |
1527 |
chrisfen |
703 |
|
1528 |
chrisfen |
700 |
ri2 = riji * riji |
1529 |
|
|
ri3 = ri2 * riji |
1530 |
|
|
|
1531 |
|
|
pref = pre12 * q_j * mu_i |
1532 |
|
|
vterm = pref * ct_i * ( ri2 - preRF2*rij ) |
1533 |
chrisfen |
703 |
myPot = myPot + sw*vterm |
1534 |
chrisfen |
700 |
|
1535 |
chrisfen |
703 |
dudx = dudx + sw*pref*( ri3*(uz_i(1)-3.0d0*ct_i*xhat) & |
1536 |
|
|
- preRF2*uz_i(1) ) |
1537 |
|
|
dudy = dudy + sw*pref*( ri3*(uz_i(2)-3.0d0*ct_i*yhat) & |
1538 |
|
|
- preRF2*uz_i(2) ) |
1539 |
|
|
dudz = dudz + sw*pref*( ri3*(uz_i(3)-3.0d0*ct_i*zhat) & |
1540 |
|
|
- preRF2*uz_i(3) ) |
1541 |
chrisfen |
700 |
|
1542 |
|
|
duduz_i(1) = duduz_i(1) + sw * pref * xhat * ( ri2 - preRF2*rij ) |
1543 |
|
|
duduz_i(2) = duduz_i(2) + sw * pref * yhat * ( ri2 - preRF2*rij ) |
1544 |
|
|
duduz_i(3) = duduz_i(3) + sw * pref * zhat * ( ri2 - preRF2*rij ) |
1545 |
|
|
|
1546 |
|
|
endif |
1547 |
chrisfen |
703 |
|
1548 |
|
|
|
1549 |
|
|
! accumulate the forces and torques resulting from the self term |
1550 |
chrisfen |
700 |
f(1,atom1) = f(1,atom1) + dudx |
1551 |
|
|
f(2,atom1) = f(2,atom1) + dudy |
1552 |
|
|
f(3,atom1) = f(3,atom1) + dudz |
1553 |
|
|
|
1554 |
|
|
f(1,atom2) = f(1,atom2) - dudx |
1555 |
|
|
f(2,atom2) = f(2,atom2) - dudy |
1556 |
|
|
f(3,atom2) = f(3,atom2) - dudz |
1557 |
|
|
|
1558 |
|
|
if (i_is_Dipole) then |
1559 |
|
|
t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2) |
1560 |
|
|
t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3) |
1561 |
|
|
t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1) |
1562 |
|
|
elseif (j_is_Dipole) then |
1563 |
|
|
t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2) |
1564 |
|
|
t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3) |
1565 |
|
|
t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1) |
1566 |
|
|
endif |
1567 |
|
|
|
1568 |
|
|
return |
1569 |
|
|
end subroutine rf_self_excludes |
1570 |
|
|
|
1571 |
gezelter |
411 |
end module electrostatic_module |