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