ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 1467
Committed: Sat Jul 17 15:33:03 2010 UTC (14 years, 10 months ago) by gezelter
File size: 48954 byte(s)
Log Message:
well, it compiles, but still segfaults

File Contents

# Content
1 !!
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. Redistributions of source code must retain the above copyright
10 !! notice, this list of conditions and the following disclaimer.
11 !!
12 !! 2. Redistributions in binary form must reproduce the above copyright
13 !! notice, this list of conditions and the following disclaimer in the
14 !! documentation and/or other materials provided with the
15 !! distribution.
16 !!
17 !! This software is provided "AS IS," without a warranty of any
18 !! kind. All express or implied conditions, representations and
19 !! warranties, including any implied warranty of merchantability,
20 !! fitness for a particular purpose or non-infringement, are hereby
21 !! excluded. The University of Notre Dame and its licensors shall not
22 !! be liable for any damages suffered by licensee as a result of
23 !! using, modifying or distributing the software or its
24 !! derivatives. In no event will the University of Notre Dame or its
25 !! licensors be liable for any lost revenue, profit or data, or for
26 !! direct, indirect, special, consequential, incidental or punitive
27 !! damages, however caused and regardless of the theory of liability,
28 !! arising out of the use of or inability to use software, even if the
29 !! University of Notre Dame has been advised of the possibility of
30 !! such damages.
31 !!
32 !! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 !! research, please cite the appropriate papers when you publish your
34 !! work. Good starting points are:
35 !!
36 !! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 !! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 !! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 !! [4] Vardeman & Gezelter, in progress (2009).
40 !!
41
42 module electrostatic_module
43
44 use force_globals
45 use definitions
46 use atype_module
47 use vector_class
48 use simulation
49 use status
50 use interpolation
51 implicit none
52
53 PRIVATE
54
55
56 #define __FORTRAN90
57 #include "UseTheForce/DarkSide/fInteractionMap.h"
58 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
60
61
62 !! these prefactors convert the multipole interactions into kcal / mol
63 !! all were computed assuming distances are measured in angstroms
64 !! Charge-Charge, assuming charges are measured in electrons
65 real(kind=dp), parameter :: pre11 = 332.0637778_dp
66 !! Charge-Dipole, assuming charges are measured in electrons, and
67 !! dipoles are measured in debyes
68 real(kind=dp), parameter :: pre12 = 69.13373_dp
69 !! Dipole-Dipole, assuming dipoles are measured in debyes
70 real(kind=dp), parameter :: pre22 = 14.39325_dp
71 !! Charge-Quadrupole, assuming charges are measured in electrons, and
72 !! quadrupoles are measured in 10^-26 esu cm^2
73 !! This unit is also known affectionately as an esu centi-barn.
74 real(kind=dp), parameter :: pre14 = 69.13373_dp
75
76 real(kind=dp), parameter :: zero = 0.0_dp
77
78 !! conversions for the simulation box dipole moment
79 real(kind=dp), parameter :: chargeToC = 1.60217733e-19_dp
80 real(kind=dp), parameter :: angstromToM = 1.0e-10_dp
81 real(kind=dp), parameter :: debyeToCm = 3.33564095198e-30_dp
82
83 !! number of points for electrostatic splines
84 integer, parameter :: np = 100
85
86 !! variables to handle different summation methods for long-range
87 !! electrostatics:
88 integer, save :: summationMethod = NONE
89 integer, save :: screeningMethod = UNDAMPED
90 logical, save :: summationMethodChecked = .false.
91 real(kind=DP), save :: defaultCutoff = 0.0_DP
92 real(kind=DP), save :: defaultCutoff2 = 0.0_DP
93 logical, save :: haveDefaultCutoff = .false.
94 real(kind=DP), save :: dampingAlpha = 0.0_DP
95 real(kind=DP), save :: alpha2 = 0.0_DP
96 real(kind=DP), save :: alpha4 = 0.0_DP
97 real(kind=DP), save :: alpha6 = 0.0_DP
98 real(kind=DP), save :: alpha8 = 0.0_DP
99 logical, save :: haveDampingAlpha = .false.
100 real(kind=DP), save :: dielectric = 1.0_DP
101 logical, save :: haveDielectric = .false.
102 real(kind=DP), save :: constEXP = 0.0_DP
103 real(kind=dp), save :: rcuti = 0.0_DP
104 real(kind=dp), save :: rcuti2 = 0.0_DP
105 real(kind=dp), save :: rcuti3 = 0.0_DP
106 real(kind=dp), save :: rcuti4 = 0.0_DP
107 real(kind=dp), save :: alphaPi = 0.0_DP
108 real(kind=dp), save :: invRootPi = 0.0_DP
109 real(kind=dp), save :: rrf = 1.0_DP
110 real(kind=dp), save :: rt = 1.0_DP
111 real(kind=dp), save :: rrfsq = 1.0_DP
112 real(kind=dp), save :: preRF = 0.0_DP
113 real(kind=dp), save :: preRF2 = 0.0_DP
114 real(kind=dp), save :: erfcVal = 1.0_DP
115 real(kind=dp), save :: derfcVal = 0.0_DP
116 type(cubicSpline), save :: erfcSpline
117 logical, save :: haveElectroSpline = .false.
118 real(kind=dp), save :: c1 = 1.0_DP
119 real(kind=dp), save :: c2 = 1.0_DP
120 real(kind=dp), save :: c3 = 0.0_DP
121 real(kind=dp), save :: c4 = 0.0_DP
122 real(kind=dp), save :: c5 = 0.0_DP
123 real(kind=dp), save :: c6 = 0.0_DP
124 real(kind=dp), save :: c1c = 1.0_DP
125 real(kind=dp), save :: c2c = 1.0_DP
126 real(kind=dp), save :: c3c = 0.0_DP
127 real(kind=dp), save :: c4c = 0.0_DP
128 real(kind=dp), save :: c5c = 0.0_DP
129 real(kind=dp), save :: c6c = 0.0_DP
130 real(kind=dp), save :: one_third = 1.0_DP / 3.0_DP
131
132 #if defined(__IFC) || defined(__PGI)
133 ! error function for ifc version > 7.
134 real(kind=dp), external :: erfc
135 #endif
136
137 public :: setElectrostaticSummationMethod
138 public :: setScreeningMethod
139 public :: setElectrostaticCutoffRadius
140 public :: setDampingAlpha
141 public :: setReactionFieldDielectric
142 public :: buildElectroSpline
143 public :: newElectrostaticType
144 public :: setCharge
145 public :: setDipoleMoment
146 public :: setSplitDipoleDistance
147 public :: setQuadrupoleMoments
148 public :: doElectrostaticPair
149 public :: getCharge
150 public :: getDipoleMoment
151 public :: destroyElectrostaticTypes
152 public :: self_self
153 public :: rf_self_excludes
154 public :: accumulate_box_dipole
155
156 type :: Electrostatic
157 integer :: c_ident
158 logical :: is_Charge = .false.
159 logical :: is_Dipole = .false.
160 logical :: is_SplitDipole = .false.
161 logical :: is_Quadrupole = .false.
162 logical :: is_Tap = .false.
163 real(kind=DP) :: charge = 0.0_DP
164 real(kind=DP) :: dipole_moment = 0.0_DP
165 real(kind=DP) :: split_dipole_distance = 0.0_DP
166 real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP
167 end type Electrostatic
168
169 type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
170
171 logical, save :: hasElectrostaticMap
172
173 contains
174
175 subroutine setElectrostaticSummationMethod(the_ESM)
176 integer, intent(in) :: the_ESM
177
178 if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then
179 call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method")
180 endif
181
182 summationMethod = the_ESM
183
184 end subroutine setElectrostaticSummationMethod
185
186 subroutine setScreeningMethod(the_SM)
187 integer, intent(in) :: the_SM
188 screeningMethod = the_SM
189 end subroutine setScreeningMethod
190
191 subroutine setElectrostaticCutoffRadius(thisRcut, thisRsw)
192 real(kind=dp), intent(in) :: thisRcut
193 real(kind=dp), intent(in) :: thisRsw
194 defaultCutoff = thisRcut
195 defaultCutoff2 = defaultCutoff*defaultCutoff
196 rrf = defaultCutoff
197 rt = thisRsw
198 haveDefaultCutoff = .true.
199 end subroutine setElectrostaticCutoffRadius
200
201 subroutine setDampingAlpha(thisAlpha)
202 real(kind=dp), intent(in) :: thisAlpha
203 dampingAlpha = thisAlpha
204 alpha2 = dampingAlpha*dampingAlpha
205 alpha4 = alpha2*alpha2
206 alpha6 = alpha4*alpha2
207 alpha8 = alpha4*alpha4
208 haveDampingAlpha = .true.
209 end subroutine setDampingAlpha
210
211 subroutine setReactionFieldDielectric(thisDielectric)
212 real(kind=dp), intent(in) :: thisDielectric
213 dielectric = thisDielectric
214 haveDielectric = .true.
215 end subroutine setReactionFieldDielectric
216
217 subroutine buildElectroSpline()
218 real( kind = dp ), dimension(np) :: xvals, yvals
219 real( kind = dp ) :: dx, rmin, rval
220 integer :: i
221
222 rmin = 0.0_dp
223
224 dx = (defaultCutoff-rmin) / dble(np-1)
225
226 do i = 1, np
227 rval = rmin + dble(i-1)*dx
228 xvals(i) = rval
229 yvals(i) = erfc(dampingAlpha*rval)
230 enddo
231
232 call newSpline(erfcSpline, xvals, yvals, .true.)
233
234 haveElectroSpline = .true.
235 end subroutine buildElectroSpline
236
237 subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
238 is_SplitDipole, is_Quadrupole, is_Tap, status)
239
240 integer, intent(in) :: c_ident
241 logical, intent(in) :: is_Charge
242 logical, intent(in) :: is_Dipole
243 logical, intent(in) :: is_SplitDipole
244 logical, intent(in) :: is_Quadrupole
245 logical, intent(in) :: is_Tap
246 integer, intent(out) :: status
247 integer :: nAtypes, myATID, i, j
248
249 status = 0
250 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
251
252 !! Be simple-minded and assume that we need an ElectrostaticMap that
253 !! is the same size as the total number of atom types
254
255 if (.not.allocated(ElectrostaticMap)) then
256
257 nAtypes = getSize(atypes)
258
259 if (nAtypes == 0) then
260 status = -1
261 return
262 end if
263
264 allocate(ElectrostaticMap(nAtypes))
265
266 end if
267
268 if (myATID .gt. size(ElectrostaticMap)) then
269 status = -1
270 return
271 endif
272
273 ! set the values for ElectrostaticMap for this atom type:
274
275 ElectrostaticMap(myATID)%c_ident = c_ident
276 ElectrostaticMap(myATID)%is_Charge = is_Charge
277 ElectrostaticMap(myATID)%is_Dipole = is_Dipole
278 ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole
279 ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole
280 ElectrostaticMap(myATID)%is_Tap = is_Tap
281
282 hasElectrostaticMap = .true.
283
284 end subroutine newElectrostaticType
285
286 subroutine setCharge(c_ident, charge, status)
287 integer, intent(in) :: c_ident
288 real(kind=dp), intent(in) :: charge
289 integer, intent(out) :: status
290 integer :: myATID
291
292 status = 0
293 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
294
295 if (.not.hasElectrostaticMap) then
296 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!")
297 status = -1
298 return
299 end if
300
301 if (myATID .gt. size(ElectrostaticMap)) then
302 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!")
303 status = -1
304 return
305 endif
306
307 if (.not.ElectrostaticMap(myATID)%is_Charge) then
308 call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!")
309 status = -1
310 return
311 endif
312
313 ElectrostaticMap(myATID)%charge = charge
314 end subroutine setCharge
315
316 subroutine setDipoleMoment(c_ident, dipole_moment, status)
317 integer, intent(in) :: c_ident
318 real(kind=dp), intent(in) :: dipole_moment
319 integer, intent(out) :: status
320 integer :: myATID
321
322 status = 0
323 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
324
325 if (.not.hasElectrostaticMap) then
326 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!")
327 status = -1
328 return
329 end if
330
331 if (myATID .gt. size(ElectrostaticMap)) then
332 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!")
333 status = -1
334 return
335 endif
336
337 if (.not.ElectrostaticMap(myATID)%is_Dipole) then
338 call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!")
339 status = -1
340 return
341 endif
342
343 ElectrostaticMap(myATID)%dipole_moment = dipole_moment
344 end subroutine setDipoleMoment
345
346 subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status)
347 integer, intent(in) :: c_ident
348 real(kind=dp), intent(in) :: split_dipole_distance
349 integer, intent(out) :: status
350 integer :: myATID
351
352 status = 0
353 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
354
355 if (.not.hasElectrostaticMap) then
356 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!")
357 status = -1
358 return
359 end if
360
361 if (myATID .gt. size(ElectrostaticMap)) then
362 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!")
363 status = -1
364 return
365 endif
366
367 if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then
368 call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!")
369 status = -1
370 return
371 endif
372
373 ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance
374 end subroutine setSplitDipoleDistance
375
376 subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status)
377 integer, intent(in) :: c_ident
378 real(kind=dp), intent(in), dimension(3) :: quadrupole_moments
379 integer, intent(out) :: status
380 integer :: myATID, i, j
381
382 status = 0
383 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
384
385 if (.not.hasElectrostaticMap) then
386 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!")
387 status = -1
388 return
389 end if
390
391 if (myATID .gt. size(ElectrostaticMap)) then
392 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!")
393 status = -1
394 return
395 endif
396
397 if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then
398 call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!")
399 status = -1
400 return
401 endif
402
403 do i = 1, 3
404 ElectrostaticMap(myATID)%quadrupole_moments(i) = &
405 quadrupole_moments(i)
406 enddo
407
408 end subroutine setQuadrupoleMoments
409
410
411 function getCharge(atid) result (c)
412 integer, intent(in) :: atid
413 integer :: localError
414 real(kind=dp) :: c
415
416 if (.not.hasElectrostaticMap) then
417 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
418 return
419 end if
420
421 if (.not.ElectrostaticMap(atid)%is_Charge) then
422 call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!")
423 return
424 endif
425
426 c = ElectrostaticMap(atid)%charge
427 end function getCharge
428
429 function getDipoleMoment(atid) result (dm)
430 integer, intent(in) :: atid
431 integer :: localError
432 real(kind=dp) :: dm
433
434 if (.not.hasElectrostaticMap) then
435 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
436 return
437 end if
438
439 if (.not.ElectrostaticMap(atid)%is_Dipole) then
440 call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!")
441 return
442 endif
443
444 dm = ElectrostaticMap(atid)%dipole_moment
445 end function getDipoleMoment
446
447 subroutine checkSummationMethod()
448
449 if (.not.haveDefaultCutoff) then
450 call handleError("checkSummationMethod", "no Default Cutoff set!")
451 endif
452
453 rcuti = 1.0_dp / defaultCutoff
454 rcuti2 = rcuti*rcuti
455 rcuti3 = rcuti2*rcuti
456 rcuti4 = rcuti2*rcuti2
457
458 if (screeningMethod .eq. DAMPED) then
459 if (.not.haveDampingAlpha) then
460 call handleError("checkSummationMethod", "no Damping Alpha set!")
461 endif
462
463 if (.not.haveDefaultCutoff) then
464 call handleError("checkSummationMethod", "no Default Cutoff set!")
465 endif
466
467 constEXP = exp(-alpha2*defaultCutoff2)
468 invRootPi = 0.56418958354775628695_dp
469 alphaPi = 2.0_dp*dampingAlpha*invRootPi
470
471 c1c = erfc(dampingAlpha*defaultCutoff) * rcuti
472 c2c = alphaPi*constEXP*rcuti + c1c*rcuti
473 c3c = 2.0_dp*alphaPi*alpha2 + 3.0_dp*c2c*rcuti
474 c4c = 4.0_dp*alphaPi*alpha4 + 5.0_dp*c3c*rcuti2
475 c5c = 8.0_dp*alphaPi*alpha6 + 7.0_dp*c4c*rcuti2
476 c6c = 16.0_dp*alphaPi*alpha8 + 9.0_dp*c5c*rcuti2
477 else
478 c1c = rcuti
479 c2c = c1c*rcuti
480 c3c = 3.0_dp*c2c*rcuti
481 c4c = 5.0_dp*c3c*rcuti2
482 c5c = 7.0_dp*c4c*rcuti2
483 c6c = 9.0_dp*c5c*rcuti2
484 endif
485
486 if (summationMethod .eq. REACTION_FIELD) then
487 if (haveDielectric) then
488 defaultCutoff2 = defaultCutoff*defaultCutoff
489 preRF = (dielectric-1.0_dp) / &
490 ((2.0_dp*dielectric+1.0_dp)*defaultCutoff2*defaultCutoff)
491 preRF2 = 2.0_dp*preRF
492 else
493 call handleError("checkSummationMethod", "Dielectric not set")
494 endif
495
496 endif
497
498 if (.not.haveElectroSpline) then
499 call buildElectroSpline()
500 end if
501
502 summationMethodChecked = .true.
503 end subroutine checkSummationMethod
504
505
506 subroutine doElectrostaticPair(me1, me2, d, rij, r2, rcut, sw, &
507 electroMult, vpair, pot, eF1, eF2, f1, t1, t2)
508
509 integer, intent(in) :: me1, me2
510 integer :: localError
511
512 real(kind=dp), intent(in) :: rij, r2, sw, rcut, electroMult
513 real(kind=dp), intent(in), dimension(3) :: d
514 real(kind=dp), intent(inout) :: vpair
515
516 real( kind = dp ) :: pot
517 real( kind = dp ), dimension(9) :: eF1, eF2 ! eFrame = electroFrame
518 real( kind = dp ), dimension(3) :: f1
519 real( kind = dp ), dimension(3,nLocal) :: felec
520 real( kind = dp ), dimension(3) :: t1, t2
521
522 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
523 real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
524 real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
525 real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
526
527 logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
528 logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
529 logical :: i_is_Tap, j_is_Tap
530 integer :: id1, id2
531 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
532 real (kind=dp) :: qxx_i, qyy_i, qzz_i
533 real (kind=dp) :: qxx_j, qyy_j, qzz_j
534 real (kind=dp) :: cx_i, cy_i, cz_i
535 real (kind=dp) :: cx_j, cy_j, cz_j
536 real (kind=dp) :: cx2, cy2, cz2
537 real (kind=dp) :: ct_i, ct_j, ct_ij, a0, a1
538 real (kind=dp) :: riji, ri, ri2, ri3, ri4
539 real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2
540 real (kind=dp) :: xhat, yhat, zhat
541 real (kind=dp) :: dudx, dudy, dudz
542 real (kind=dp) :: scale, sc2, bigR
543 real (kind=dp) :: varEXP
544 real (kind=dp) :: pot_term
545 real (kind=dp) :: preVal, rfVal
546 real (kind=dp) :: c2ri, c3ri, c4rij
547 real (kind=dp) :: cti3, ctj3, ctidotj
548 real (kind=dp) :: preSw, preSwSc
549 real (kind=dp) :: xhatdot2, yhatdot2, zhatdot2
550 real (kind=dp) :: xhatc4, yhatc4, zhatc4
551
552 if (.not.summationMethodChecked) then
553 call checkSummationMethod()
554 endif
555
556 !! some variables we'll need independent of electrostatic type:
557
558 riji = 1.0_dp / rij
559
560 xhat = d(1) * riji
561 yhat = d(2) * riji
562 zhat = d(3) * riji
563
564 !! logicals
565 i_is_Charge = ElectrostaticMap(me1)%is_Charge
566 i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
567 i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
568 i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
569 i_is_Tap = ElectrostaticMap(me1)%is_Tap
570
571 j_is_Charge = ElectrostaticMap(me2)%is_Charge
572 j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
573 j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
574 j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
575 j_is_Tap = ElectrostaticMap(me2)%is_Tap
576
577 if (i_is_Charge) then
578 q_i = ElectrostaticMap(me1)%charge
579 endif
580
581 if (i_is_Dipole) then
582 mu_i = ElectrostaticMap(me1)%dipole_moment
583
584 uz_i(1) = eF1(3)
585 uz_i(2) = eF1(6)
586 uz_i(3) = eF1(9)
587
588 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
589
590 if (i_is_SplitDipole) then
591 d_i = ElectrostaticMap(me1)%split_dipole_distance
592 endif
593 duduz_i = zero
594 endif
595
596 if (i_is_Quadrupole) then
597 qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
598 qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
599 qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
600
601 ux_i(1) = eF1(1)
602 ux_i(2) = eF1(4)
603 ux_i(3) = eF1(7)
604 uy_i(1) = eF1(2)
605 uy_i(2) = eF1(5)
606 uy_i(3) = eF1(8)
607 uz_i(1) = eF1(3)
608 uz_i(2) = eF1(6)
609 uz_i(3) = eF1(9)
610
611 cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
612 cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
613 cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
614 dudux_i = zero
615 duduy_i = zero
616 duduz_i = zero
617 endif
618
619 if (j_is_Charge) then
620 q_j = ElectrostaticMap(me2)%charge
621 endif
622
623 if (j_is_Dipole) then
624 mu_j = ElectrostaticMap(me2)%dipole_moment
625
626 uz_j(1) = eF2(3)
627 uz_j(2) = eF2(6)
628 uz_j(3) = eF2(9)
629
630 ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
631
632 if (j_is_SplitDipole) then
633 d_j = ElectrostaticMap(me2)%split_dipole_distance
634 endif
635 duduz_j = zero
636 endif
637
638 if (j_is_Quadrupole) then
639 qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
640 qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
641 qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
642
643 ux_j(1) = eF2(1)
644 ux_j(2) = eF2(4)
645 ux_j(3) = eF2(7)
646 uy_j(1) = eF2(2)
647 uy_j(2) = eF2(5)
648 uy_j(3) = eF2(8)
649 uz_j(1) = eF2(3)
650 uz_j(2) = eF2(6)
651 uz_j(3) = eF2(9)
652
653 cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
654 cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
655 cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
656 dudux_j = zero
657 duduy_j = zero
658 duduz_j = zero
659 endif
660
661 epot = zero
662 dudx = zero
663 dudy = zero
664 dudz = zero
665
666 if (i_is_Charge) then
667
668 if (j_is_Charge) then
669 if (screeningMethod .eq. DAMPED) then
670 ! assemble the damping variables
671 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
672 c1 = erfcVal*riji
673 c2 = (-derfcVal + c1)*riji
674 else
675 c1 = riji
676 c2 = c1*riji
677 endif
678
679 preVal = electroMult * pre11 * q_i * q_j
680
681 if (summationMethod .eq. SHIFTED_POTENTIAL) then
682 vterm = preVal * (c1 - c1c)
683
684 dudr = -sw * preVal * c2
685
686 elseif (summationMethod .eq. SHIFTED_FORCE) then
687 vterm = preVal * ( c1 - c1c + c2c*(rij - defaultCutoff) )
688
689 dudr = sw * preVal * (c2c - c2)
690
691 elseif (summationMethod .eq. REACTION_FIELD) then
692 rfVal = electroMult * preRF*rij*rij
693 vterm = preVal * ( riji + rfVal )
694
695 dudr = sw * preVal * ( 2.0_dp*rfVal - riji )*riji
696
697 else
698 vterm = preVal * riji*erfcVal
699
700 dudr = - sw * preVal * c2
701
702 endif
703
704 vpair = vpair + vterm
705 epot = epot + sw*vterm
706
707 dudx = dudx + dudr * xhat
708 dudy = dudy + dudr * yhat
709 dudz = dudz + dudr * zhat
710
711 endif
712
713 if (j_is_Dipole) then
714 ! pref is used by all the possible methods
715 pref = electroMult * pre12 * q_i * mu_j
716 preSw = sw*pref
717
718 if (summationMethod .eq. REACTION_FIELD) then
719 ri2 = riji * riji
720 ri3 = ri2 * riji
721
722 vterm = - pref * ct_j * ( ri2 - preRF2*rij )
723 vpair = vpair + vterm
724 epot = epot + sw*vterm
725
726 dudx = dudx - preSw*( ri3*(uz_j(1) - 3.0_dp*ct_j*xhat) - &
727 preRF2*uz_j(1) )
728 dudy = dudy - preSw*( ri3*(uz_j(2) - 3.0_dp*ct_j*yhat) - &
729 preRF2*uz_j(2) )
730 dudz = dudz - preSw*( ri3*(uz_j(3) - 3.0_dp*ct_j*zhat) - &
731 preRF2*uz_j(3) )
732 duduz_j(1) = duduz_j(1) - preSw * xhat * ( ri2 - preRF2*rij )
733 duduz_j(2) = duduz_j(2) - preSw * yhat * ( ri2 - preRF2*rij )
734 duduz_j(3) = duduz_j(3) - preSw * zhat * ( ri2 - preRF2*rij )
735
736 else
737 ! determine the inverse r used if we have split dipoles
738 if (j_is_SplitDipole) then
739 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
740 ri = 1.0_dp / BigR
741 scale = rij * ri
742 else
743 ri = riji
744 scale = 1.0_dp
745 endif
746
747 sc2 = scale * scale
748
749 if (screeningMethod .eq. DAMPED) then
750 ! assemble the damping variables
751 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
752 c1 = erfcVal*ri
753 c2 = (-derfcVal + c1)*ri
754 c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
755 else
756 c1 = ri
757 c2 = c1*ri
758 c3 = 3.0_dp*c2*ri
759 endif
760
761 c2ri = c2*ri
762
763 ! calculate the potential
764 pot_term = scale * c2
765 vterm = -pref * ct_j * pot_term
766 vpair = vpair + vterm
767 epot = epot + sw*vterm
768
769 ! calculate derivatives for forces and torques
770 dudx = dudx - preSw*( uz_j(1)*c2ri - ct_j*xhat*sc2*c3 )
771 dudy = dudy - preSw*( uz_j(2)*c2ri - ct_j*yhat*sc2*c3 )
772 dudz = dudz - preSw*( uz_j(3)*c2ri - ct_j*zhat*sc2*c3 )
773
774 duduz_j(1) = duduz_j(1) - preSw * pot_term * xhat
775 duduz_j(2) = duduz_j(2) - preSw * pot_term * yhat
776 duduz_j(3) = duduz_j(3) - preSw * pot_term * zhat
777
778 endif
779 endif
780
781 if (j_is_Quadrupole) then
782 ! first precalculate some necessary variables
783 cx2 = cx_j * cx_j
784 cy2 = cy_j * cy_j
785 cz2 = cz_j * cz_j
786 pref = electroMult * pre14 * q_i * one_third
787
788 if (screeningMethod .eq. DAMPED) then
789 ! assemble the damping variables
790 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
791 c1 = erfcVal*riji
792 c2 = (-derfcVal + c1)*riji
793 c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*riji
794 c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*riji*riji
795 else
796 c1 = riji
797 c2 = c1*riji
798 c3 = 3.0_dp*c2*riji
799 c4 = 5.0_dp*c3*riji*riji
800 endif
801
802 ! precompute variables for convenience
803 preSw = sw*pref
804 c2ri = c2*riji
805 c3ri = c3*riji
806 c4rij = c4*rij
807 xhatdot2 = 2.0_dp*xhat*c3
808 yhatdot2 = 2.0_dp*yhat*c3
809 zhatdot2 = 2.0_dp*zhat*c3
810 xhatc4 = xhat*c4rij
811 yhatc4 = yhat*c4rij
812 zhatc4 = zhat*c4rij
813
814 ! calculate the potential
815 pot_term = ( qxx_j*(cx2*c3 - c2ri) + qyy_j*(cy2*c3 - c2ri) + &
816 qzz_j*(cz2*c3 - c2ri) )
817 vterm = pref * pot_term
818 vpair = vpair + vterm
819 epot = epot + sw*vterm
820
821 ! calculate derivatives for the forces and torques
822 dudx = dudx - preSw * ( &
823 qxx_j*(cx2*xhatc4 - (2.0_dp*cx_j*ux_j(1) + xhat)*c3ri) + &
824 qyy_j*(cy2*xhatc4 - (2.0_dp*cy_j*uy_j(1) + xhat)*c3ri) + &
825 qzz_j*(cz2*xhatc4 - (2.0_dp*cz_j*uz_j(1) + xhat)*c3ri) )
826 dudy = dudy - preSw * ( &
827 qxx_j*(cx2*yhatc4 - (2.0_dp*cx_j*ux_j(2) + yhat)*c3ri) + &
828 qyy_j*(cy2*yhatc4 - (2.0_dp*cy_j*uy_j(2) + yhat)*c3ri) + &
829 qzz_j*(cz2*yhatc4 - (2.0_dp*cz_j*uz_j(2) + yhat)*c3ri) )
830 dudz = dudz - preSw * ( &
831 qxx_j*(cx2*zhatc4 - (2.0_dp*cx_j*ux_j(3) + zhat)*c3ri) + &
832 qyy_j*(cy2*zhatc4 - (2.0_dp*cy_j*uy_j(3) + zhat)*c3ri) + &
833 qzz_j*(cz2*zhatc4 - (2.0_dp*cz_j*uz_j(3) + zhat)*c3ri) )
834
835 dudux_j(1) = dudux_j(1) + preSw*(qxx_j*cx_j*xhatdot2)
836 dudux_j(2) = dudux_j(2) + preSw*(qxx_j*cx_j*yhatdot2)
837 dudux_j(3) = dudux_j(3) + preSw*(qxx_j*cx_j*zhatdot2)
838
839 duduy_j(1) = duduy_j(1) + preSw*(qyy_j*cy_j*xhatdot2)
840 duduy_j(2) = duduy_j(2) + preSw*(qyy_j*cy_j*yhatdot2)
841 duduy_j(3) = duduy_j(3) + preSw*(qyy_j*cy_j*zhatdot2)
842
843 duduz_j(1) = duduz_j(1) + preSw*(qzz_j*cz_j*xhatdot2)
844 duduz_j(2) = duduz_j(2) + preSw*(qzz_j*cz_j*yhatdot2)
845 duduz_j(3) = duduz_j(3) + preSw*(qzz_j*cz_j*zhatdot2)
846
847
848 endif
849 endif
850
851 if (i_is_Dipole) then
852
853 if (j_is_Charge) then
854 ! variables used by all the methods
855 pref = electroMult * pre12 * q_j * mu_i
856 preSw = sw*pref
857
858 if (summationMethod .eq. REACTION_FIELD) then
859
860 ri2 = riji * riji
861 ri3 = ri2 * riji
862
863 vterm = pref * ct_i * ( ri2 - preRF2*rij )
864 vpair = vpair + vterm
865 epot = epot + sw*vterm
866
867 dudx = dudx + preSw * ( ri3*(uz_i(1) - 3.0_dp*ct_i*xhat) - &
868 preRF2*uz_i(1) )
869 dudy = dudy + preSw * ( ri3*(uz_i(2) - 3.0_dp*ct_i*yhat) - &
870 preRF2*uz_i(2) )
871 dudz = dudz + preSw * ( ri3*(uz_i(3) - 3.0_dp*ct_i*zhat) - &
872 preRF2*uz_i(3) )
873
874 duduz_i(1) = duduz_i(1) + preSw * xhat * ( ri2 - preRF2*rij )
875 duduz_i(2) = duduz_i(2) + preSw * yhat * ( ri2 - preRF2*rij )
876 duduz_i(3) = duduz_i(3) + preSw * zhat * ( ri2 - preRF2*rij )
877
878 else
879 ! determine inverse r if we are using split dipoles
880 if (i_is_SplitDipole) then
881 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
882 ri = 1.0_dp / BigR
883 scale = rij * ri
884 else
885 ri = riji
886 scale = 1.0_dp
887 endif
888
889 sc2 = scale * scale
890
891 if (screeningMethod .eq. DAMPED) then
892 ! assemble the damping variables
893 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
894 c1 = erfcVal*ri
895 c2 = (-derfcVal + c1)*ri
896 c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
897 else
898 c1 = ri
899 c2 = c1*ri
900 c3 = 3.0_dp*c2*ri
901 endif
902
903 c2ri = c2*ri
904
905 ! calculate the potential
906 pot_term = c2 * scale
907 vterm = pref * ct_i * pot_term
908 vpair = vpair + vterm
909 epot = epot + sw*vterm
910
911 ! calculate derivatives for the forces and torques
912 dudx = dudx + preSw * ( uz_i(1)*c2ri - ct_i*xhat*sc2*c3 )
913 dudy = dudy + preSw * ( uz_i(2)*c2ri - ct_i*yhat*sc2*c3 )
914 dudz = dudz + preSw * ( uz_i(3)*c2ri - ct_i*zhat*sc2*c3 )
915
916 duduz_i(1) = duduz_i(1) + preSw * pot_term * xhat
917 duduz_i(2) = duduz_i(2) + preSw * pot_term * yhat
918 duduz_i(3) = duduz_i(3) + preSw * pot_term * zhat
919
920 endif
921 endif
922
923 if (j_is_Dipole) then
924 ! variables used by all methods
925 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
926 pref = electroMult * pre22 * mu_i * mu_j
927 preSw = sw*pref
928
929 if (summationMethod .eq. REACTION_FIELD) then
930 ri2 = riji * riji
931 ri3 = ri2 * riji
932 ri4 = ri2 * ri2
933
934 vterm = pref*( ri3*(ct_ij - 3.0_dp * ct_i * ct_j) - &
935 preRF2*ct_ij )
936 vpair = vpair + vterm
937 epot = epot + sw*vterm
938
939 a1 = 5.0_dp * ct_i * ct_j - ct_ij
940
941 dudx = dudx + preSw*3.0_dp*ri4*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
942 dudy = dudy + preSw*3.0_dp*ri4*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
943 dudz = dudz + preSw*3.0_dp*ri4*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
944
945 duduz_i(1) = duduz_i(1) + preSw*(ri3*(uz_j(1)-3.0_dp*ct_j*xhat) &
946 - preRF2*uz_j(1))
947 duduz_i(2) = duduz_i(2) + preSw*(ri3*(uz_j(2)-3.0_dp*ct_j*yhat) &
948 - preRF2*uz_j(2))
949 duduz_i(3) = duduz_i(3) + preSw*(ri3*(uz_j(3)-3.0_dp*ct_j*zhat) &
950 - preRF2*uz_j(3))
951 duduz_j(1) = duduz_j(1) + preSw*(ri3*(uz_i(1)-3.0_dp*ct_i*xhat) &
952 - preRF2*uz_i(1))
953 duduz_j(2) = duduz_j(2) + preSw*(ri3*(uz_i(2)-3.0_dp*ct_i*yhat) &
954 - preRF2*uz_i(2))
955 duduz_j(3) = duduz_j(3) + preSw*(ri3*(uz_i(3)-3.0_dp*ct_i*zhat) &
956 - preRF2*uz_i(3))
957
958 else
959 if (i_is_SplitDipole) then
960 if (j_is_SplitDipole) then
961 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
962 else
963 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
964 endif
965 ri = 1.0_dp / BigR
966 scale = rij * ri
967 else
968 if (j_is_SplitDipole) then
969 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
970 ri = 1.0_dp / BigR
971 scale = rij * ri
972 else
973 ri = riji
974 scale = 1.0_dp
975 endif
976 endif
977
978 if (screeningMethod .eq. DAMPED) then
979 ! assemble the damping variables
980 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
981 c1 = erfcVal*ri
982 c2 = (-derfcVal + c1)*ri
983 c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
984 c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*ri*ri
985 else
986 c1 = ri
987 c2 = c1*ri
988 c3 = 3.0_dp*c2*ri
989 c4 = 5.0_dp*c3*ri*ri
990 endif
991
992 ! precompute variables for convenience
993 sc2 = scale * scale
994 cti3 = ct_i*sc2*c3
995 ctj3 = ct_j*sc2*c3
996 ctidotj = ct_i * ct_j * sc2
997 preSwSc = preSw*scale
998 c2ri = c2*ri
999 c3ri = c3*ri
1000 c4rij = c4*rij
1001
1002
1003 ! calculate the potential
1004 pot_term = (ct_ij*c2ri - ctidotj*c3)
1005 vterm = pref * pot_term
1006 vpair = vpair + vterm
1007 epot = epot + sw*vterm
1008
1009 ! calculate derivatives for the forces and torques
1010 dudx = dudx + preSwSc * ( ctidotj*xhat*c4rij - &
1011 (ct_i*uz_j(1) + ct_j*uz_i(1) + ct_ij*xhat)*c3ri )
1012 dudy = dudy + preSwSc * ( ctidotj*yhat*c4rij - &
1013 (ct_i*uz_j(2) + ct_j*uz_i(2) + ct_ij*yhat)*c3ri )
1014 dudz = dudz + preSwSc * ( ctidotj*zhat*c4rij - &
1015 (ct_i*uz_j(3) + ct_j*uz_i(3) + ct_ij*zhat)*c3ri )
1016
1017 duduz_i(1) = duduz_i(1) + preSw * ( uz_j(1)*c2ri - ctj3*xhat )
1018 duduz_i(2) = duduz_i(2) + preSw * ( uz_j(2)*c2ri - ctj3*yhat )
1019 duduz_i(3) = duduz_i(3) + preSw * ( uz_j(3)*c2ri - ctj3*zhat )
1020
1021 duduz_j(1) = duduz_j(1) + preSw * ( uz_i(1)*c2ri - cti3*xhat )
1022 duduz_j(2) = duduz_j(2) + preSw * ( uz_i(2)*c2ri - cti3*yhat )
1023 duduz_j(3) = duduz_j(3) + preSw * ( uz_i(3)*c2ri - cti3*zhat )
1024
1025 endif
1026 endif
1027 endif
1028
1029 if (i_is_Quadrupole) then
1030 if (j_is_Charge) then
1031 ! precompute some necessary variables
1032 cx2 = cx_i * cx_i
1033 cy2 = cy_i * cy_i
1034 cz2 = cz_i * cz_i
1035 pref = electroMult * pre14 * q_j * one_third
1036
1037 if (screeningMethod .eq. DAMPED) then
1038 ! assemble the damping variables
1039 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
1040 c1 = erfcVal*riji
1041 c2 = (-derfcVal + c1)*riji
1042 c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*riji
1043 c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*riji*riji
1044 else
1045 c1 = riji
1046 c2 = c1*riji
1047 c3 = 3.0_dp*c2*riji
1048 c4 = 5.0_dp*c3*riji*riji
1049 endif
1050
1051 ! precompute some variables for convenience
1052 preSw = sw*pref
1053 c2ri = c2*riji
1054 c3ri = c3*riji
1055 c4rij = c4*rij
1056 xhatdot2 = 2.0_dp*xhat*c3
1057 yhatdot2 = 2.0_dp*yhat*c3
1058 zhatdot2 = 2.0_dp*zhat*c3
1059 xhatc4 = xhat*c4rij
1060 yhatc4 = yhat*c4rij
1061 zhatc4 = zhat*c4rij
1062
1063 ! calculate the potential
1064 pot_term = ( qxx_i * (cx2*c3 - c2ri) + qyy_i * (cy2*c3 - c2ri) + &
1065 qzz_i * (cz2*c3 - c2ri) )
1066
1067 vterm = pref * pot_term
1068 vpair = vpair + vterm
1069 epot = epot + sw*vterm
1070
1071 ! calculate the derivatives for the forces and torques
1072 dudx = dudx - preSw * ( &
1073 qxx_i*(cx2*xhatc4 - (2.0_dp*cx_i*ux_i(1) + xhat)*c3ri) + &
1074 qyy_i*(cy2*xhatc4 - (2.0_dp*cy_i*uy_i(1) + xhat)*c3ri) + &
1075 qzz_i*(cz2*xhatc4 - (2.0_dp*cz_i*uz_i(1) + xhat)*c3ri) )
1076 dudy = dudy - preSw * ( &
1077 qxx_i*(cx2*yhatc4 - (2.0_dp*cx_i*ux_i(2) + yhat)*c3ri) + &
1078 qyy_i*(cy2*yhatc4 - (2.0_dp*cy_i*uy_i(2) + yhat)*c3ri) + &
1079 qzz_i*(cz2*yhatc4 - (2.0_dp*cz_i*uz_i(2) + yhat)*c3ri) )
1080 dudz = dudz - preSw * ( &
1081 qxx_i*(cx2*zhatc4 - (2.0_dp*cx_i*ux_i(3) + zhat)*c3ri) + &
1082 qyy_i*(cy2*zhatc4 - (2.0_dp*cy_i*uy_i(3) + zhat)*c3ri) + &
1083 qzz_i*(cz2*zhatc4 - (2.0_dp*cz_i*uz_i(3) + zhat)*c3ri) )
1084
1085 dudux_i(1) = dudux_i(1) + preSw*(qxx_i*cx_i*xhatdot2)
1086 dudux_i(2) = dudux_i(2) + preSw*(qxx_i*cx_i*yhatdot2)
1087 dudux_i(3) = dudux_i(3) + preSw*(qxx_i*cx_i*zhatdot2)
1088
1089 duduy_i(1) = duduy_i(1) + preSw*(qyy_i*cy_i*xhatdot2)
1090 duduy_i(2) = duduy_i(2) + preSw*(qyy_i*cy_i*yhatdot2)
1091 duduy_i(3) = duduy_i(3) + preSw*(qyy_i*cy_i*zhatdot2)
1092
1093 duduz_i(1) = duduz_i(1) + preSw*(qzz_i*cz_i*xhatdot2)
1094 duduz_i(2) = duduz_i(2) + preSw*(qzz_i*cz_i*yhatdot2)
1095 duduz_i(3) = duduz_i(3) + preSw*(qzz_i*cz_i*zhatdot2)
1096 endif
1097 endif
1098
1099 pot = pot + epot
1100
1101 f1(1) = f1(1) + dudx
1102 f1(2) = f1(2) + dudy
1103 f1(3) = f1(3) + dudz
1104
1105 if (i_is_Dipole .or. i_is_Quadrupole) then
1106 t1(1) = t1(1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1107 t1(2) = t1(2) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1108 t1(3) = t1(3) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1109 endif
1110 if (i_is_Quadrupole) then
1111 t1(1) = t1(1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
1112 t1(2) = t1(2) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
1113 t1(3) = t1(3) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
1114
1115 t1(1) = t1(1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
1116 t1(2) = t1(2) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
1117 t1(3) = t1(3) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
1118 endif
1119
1120 if (j_is_Dipole .or. j_is_Quadrupole) then
1121 t2(1) = t2(1) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1122 t2(2) = t2(2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1123 t2(3) = t2(3) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1124 endif
1125 if (j_is_Quadrupole) then
1126 t2(1) = t2(1) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
1127 t2(2) = t2(2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
1128 t2(3) = t2(3) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
1129
1130 t2(1) = t2(1) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
1131 t2(2) = t2(2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
1132 t2(3) = t2(3) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
1133 endif
1134
1135 return
1136 end subroutine doElectrostaticPair
1137
1138 subroutine destroyElectrostaticTypes()
1139
1140 if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap)
1141
1142 end subroutine destroyElectrostaticTypes
1143
1144 subroutine self_self(atom1, eFrame, skch, mypot, t)
1145 integer, intent(in) :: atom1
1146 integer :: atid1
1147 real(kind=dp), dimension(9,nLocal) :: eFrame
1148 real(kind=dp), dimension(3,nLocal) :: t
1149 real(kind=dp) :: mu1, chg1, c1e
1150 real(kind=dp) :: preVal, epot, mypot, skch
1151 real(kind=dp) :: eix, eiy, eiz, self
1152
1153 ! this is a local only array, so we use the local atom type id's:
1154 atid1 = atid(atom1)
1155
1156 if (.not.summationMethodChecked) then
1157 call checkSummationMethod()
1158 endif
1159
1160 if (summationMethod .eq. REACTION_FIELD) then
1161 if (ElectrostaticMap(atid1)%is_Dipole) then
1162 mu1 = getDipoleMoment(atid1)
1163
1164 preVal = pre22 * preRF2 * mu1*mu1
1165 mypot = mypot - 0.5_dp*preVal
1166
1167 ! The self-correction term adds into the reaction field vector
1168
1169 eix = preVal * eFrame(3,atom1)
1170 eiy = preVal * eFrame(6,atom1)
1171 eiz = preVal * eFrame(9,atom1)
1172
1173 ! once again, this is self-self, so only the local arrays are needed
1174 ! even for MPI jobs:
1175
1176 t(1,atom1)=t(1,atom1) - eFrame(6,atom1)*eiz + &
1177 eFrame(9,atom1)*eiy
1178 t(2,atom1)=t(2,atom1) - eFrame(9,atom1)*eix + &
1179 eFrame(3,atom1)*eiz
1180 t(3,atom1)=t(3,atom1) - eFrame(3,atom1)*eiy + &
1181 eFrame(6,atom1)*eix
1182
1183 endif
1184
1185 elseif ( (summationMethod .eq. SHIFTED_FORCE) .or. &
1186 (summationMethod .eq. SHIFTED_POTENTIAL) ) then
1187 if (ElectrostaticMap(atid1)%is_Charge) then
1188 chg1 = getCharge(atid1)
1189 if (screeningMethod .eq. DAMPED) then
1190 self = - 0.5_dp * (c1c+alphaPi) * chg1 * (chg1+skch) * pre11
1191 else
1192 self = - 0.5_dp * rcuti * chg1 * (chg1+skch) * pre11
1193 endif
1194
1195 mypot = mypot + self
1196 endif
1197 endif
1198
1199
1200
1201 return
1202 end subroutine self_self
1203
1204 subroutine rf_self_excludes(atom1, atom2, sw, electroMult, eFrame, d, &
1205 rij, vpair, myPot, f, t)
1206 integer, intent(in) :: atom1
1207 integer, intent(in) :: atom2
1208 logical :: i_is_Charge, j_is_Charge
1209 logical :: i_is_Dipole, j_is_Dipole
1210 integer :: atid1
1211 integer :: atid2
1212 real(kind=dp), intent(in) :: rij
1213 real(kind=dp), intent(in) :: sw, electroMult
1214 real(kind=dp), intent(in), dimension(3) :: d
1215 real(kind=dp), intent(inout) :: vpair
1216 real(kind=dp), dimension(9,nLocal) :: eFrame
1217 real(kind=dp), dimension(3,nLocal) :: f
1218 real(kind=dp), dimension(3,nLocal) :: t
1219 real (kind = dp), dimension(3) :: duduz_i
1220 real (kind = dp), dimension(3) :: duduz_j
1221 real (kind = dp), dimension(3) :: uz_i
1222 real (kind = dp), dimension(3) :: uz_j
1223 real(kind=dp) :: q_i, q_j, mu_i, mu_j
1224 real(kind=dp) :: xhat, yhat, zhat
1225 real(kind=dp) :: ct_i, ct_j
1226 real(kind=dp) :: ri2, ri3, riji, vterm
1227 real(kind=dp) :: pref, preVal, rfVal, myPot
1228 real(kind=dp) :: dudx, dudy, dudz, dudr
1229
1230 if (.not.summationMethodChecked) then
1231 call checkSummationMethod()
1232 endif
1233
1234 dudx = zero
1235 dudy = zero
1236 dudz = zero
1237
1238 riji = 1.0_dp/rij
1239
1240 xhat = d(1) * riji
1241 yhat = d(2) * riji
1242 zhat = d(3) * riji
1243
1244 ! this is a local only array, so we use the local atom type id's:
1245 atid1 = atid(atom1)
1246 atid2 = atid(atom2)
1247 i_is_Charge = ElectrostaticMap(atid1)%is_Charge
1248 j_is_Charge = ElectrostaticMap(atid2)%is_Charge
1249 i_is_Dipole = ElectrostaticMap(atid1)%is_Dipole
1250 j_is_Dipole = ElectrostaticMap(atid2)%is_Dipole
1251
1252 if (i_is_Charge.and.j_is_Charge) then
1253 q_i = ElectrostaticMap(atid1)%charge
1254 q_j = ElectrostaticMap(atid2)%charge
1255
1256 preVal = electroMult * pre11 * q_i * q_j
1257 rfVal = preRF*rij*rij
1258 vterm = preVal * rfVal
1259
1260 myPot = myPot + sw*vterm
1261
1262 dudr = sw*preVal * 2.0_dp*rfVal*riji
1263
1264 dudx = dudx + dudr * xhat
1265 dudy = dudy + dudr * yhat
1266 dudz = dudz + dudr * zhat
1267
1268 elseif (i_is_Charge.and.j_is_Dipole) then
1269 q_i = ElectrostaticMap(atid1)%charge
1270 mu_j = ElectrostaticMap(atid2)%dipole_moment
1271 uz_j(1) = eFrame(3,atom2)
1272 uz_j(2) = eFrame(6,atom2)
1273 uz_j(3) = eFrame(9,atom2)
1274 ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
1275
1276 ri2 = riji * riji
1277 ri3 = ri2 * riji
1278
1279 pref = electroMult * pre12 * q_i * mu_j
1280 vterm = - pref * ct_j * ( ri2 - preRF2*rij )
1281 myPot = myPot + sw*vterm
1282
1283 dudx = dudx - sw*pref*( ri3*(uz_j(1)-3.0_dp*ct_j*xhat) &
1284 - preRF2*uz_j(1) )
1285 dudy = dudy - sw*pref*( ri3*(uz_j(2)-3.0_dp*ct_j*yhat) &
1286 - preRF2*uz_j(2) )
1287 dudz = dudz - sw*pref*( ri3*(uz_j(3)-3.0_dp*ct_j*zhat) &
1288 - preRF2*uz_j(3) )
1289
1290 duduz_j(1) = duduz_j(1) - sw * pref * xhat * ( ri2 - preRF2*rij )
1291 duduz_j(2) = duduz_j(2) - sw * pref * yhat * ( ri2 - preRF2*rij )
1292 duduz_j(3) = duduz_j(3) - sw * pref * zhat * ( ri2 - preRF2*rij )
1293
1294 elseif (i_is_Dipole.and.j_is_Charge) then
1295 mu_i = ElectrostaticMap(atid1)%dipole_moment
1296 q_j = ElectrostaticMap(atid2)%charge
1297 uz_i(1) = eFrame(3,atom1)
1298 uz_i(2) = eFrame(6,atom1)
1299 uz_i(3) = eFrame(9,atom1)
1300 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
1301
1302 ri2 = riji * riji
1303 ri3 = ri2 * riji
1304
1305 pref = electroMult * pre12 * q_j * mu_i
1306 vterm = pref * ct_i * ( ri2 - preRF2*rij )
1307 myPot = myPot + sw*vterm
1308
1309 dudx = dudx + sw*pref*( ri3*(uz_i(1)-3.0_dp*ct_i*xhat) &
1310 - preRF2*uz_i(1) )
1311 dudy = dudy + sw*pref*( ri3*(uz_i(2)-3.0_dp*ct_i*yhat) &
1312 - preRF2*uz_i(2) )
1313 dudz = dudz + sw*pref*( ri3*(uz_i(3)-3.0_dp*ct_i*zhat) &
1314 - preRF2*uz_i(3) )
1315
1316 duduz_i(1) = duduz_i(1) + sw * pref * xhat * ( ri2 - preRF2*rij )
1317 duduz_i(2) = duduz_i(2) + sw * pref * yhat * ( ri2 - preRF2*rij )
1318 duduz_i(3) = duduz_i(3) + sw * pref * zhat * ( ri2 - preRF2*rij )
1319
1320 endif
1321
1322
1323 ! accumulate the forces and torques resulting from the self term
1324 f(1,atom1) = f(1,atom1) + dudx
1325 f(2,atom1) = f(2,atom1) + dudy
1326 f(3,atom1) = f(3,atom1) + dudz
1327
1328 f(1,atom2) = f(1,atom2) - dudx
1329 f(2,atom2) = f(2,atom2) - dudy
1330 f(3,atom2) = f(3,atom2) - dudz
1331
1332 if (i_is_Dipole) then
1333 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1334 t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1335 t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1336 elseif (j_is_Dipole) then
1337 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1338 t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1339 t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1340 endif
1341
1342 return
1343 end subroutine rf_self_excludes
1344
1345 subroutine accumulate_box_dipole(atom1, eFrame, d, pChg, nChg, pChgPos, &
1346 nChgPos, dipVec, pChgCount, nChgCount)
1347 integer, intent(in) :: atom1
1348 logical :: i_is_Charge
1349 logical :: i_is_Dipole
1350 integer :: atid1
1351 integer :: pChgCount
1352 integer :: nChgCount
1353 real(kind=dp), intent(in), dimension(3) :: d
1354 real(kind=dp), dimension(9,nLocal) :: eFrame
1355 real(kind=dp) :: pChg
1356 real(kind=dp) :: nChg
1357 real(kind=dp), dimension(3) :: pChgPos
1358 real(kind=dp), dimension(3) :: nChgPos
1359 real(kind=dp), dimension(3) :: dipVec
1360 real(kind=dp), dimension(3) :: uz_i
1361 real(kind=dp), dimension(3) :: pos
1362 real(kind=dp) :: q_i, mu_i
1363 real(kind=dp) :: pref, preVal
1364
1365 if (.not.summationMethodChecked) then
1366 call checkSummationMethod()
1367 endif
1368
1369 ! this is a local only array, so we use the local atom type id's:
1370 atid1 = atid(atom1)
1371 i_is_Charge = ElectrostaticMap(atid1)%is_Charge
1372 i_is_Dipole = ElectrostaticMap(atid1)%is_Dipole
1373
1374 if (i_is_Charge) then
1375 q_i = ElectrostaticMap(atid1)%charge
1376 ! convert to the proper units
1377 q_i = q_i * chargeToC
1378 pos = d * angstromToM
1379
1380 if (q_i.le.0.0_dp) then
1381 nChg = nChg - q_i
1382 nChgPos(1) = nChgPos(1) + pos(1)
1383 nChgPos(2) = nChgPos(2) + pos(2)
1384 nChgPos(3) = nChgPos(3) + pos(3)
1385 nChgCount = nChgCount + 1
1386
1387 else
1388 pChg = pChg + q_i
1389 pChgPos(1) = pChgPos(1) + pos(1)
1390 pChgPos(2) = pChgPos(2) + pos(2)
1391 pChgPos(3) = pChgPos(3) + pos(3)
1392 pChgCount = pChgCount + 1
1393
1394 endif
1395
1396 endif
1397
1398 if (i_is_Dipole) then
1399 mu_i = ElectrostaticMap(atid1)%dipole_moment
1400 uz_i(1) = eFrame(3,atom1)
1401 uz_i(2) = eFrame(6,atom1)
1402 uz_i(3) = eFrame(9,atom1)
1403 ! convert to the proper units
1404 mu_i = mu_i * debyeToCm
1405
1406 dipVec(1) = dipVec(1) + uz_i(1)*mu_i
1407 dipVec(2) = dipVec(2) + uz_i(2)*mu_i
1408 dipVec(3) = dipVec(3) + uz_i(3)*mu_i
1409
1410 endif
1411
1412 return
1413 end subroutine accumulate_box_dipole
1414
1415 end module electrostatic_module

Properties

Name Value
svn:keywords Author Id Revision Date