ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 1022
Committed: Mon Aug 14 20:32:20 2006 UTC (18 years, 9 months ago) by chrisfen
File size: 53030 byte(s)
Log Message:
fixed an ordering issue in quadrupole-charge interactions

File Contents

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