ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2717
Committed: Mon Apr 17 21:49:12 2006 UTC (18 years, 5 months ago) by gezelter
File size: 52362 byte(s)
Log Message:
Many performance improvements

File Contents

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