ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2325
Committed: Sat Sep 24 17:39:36 2005 UTC (18 years, 11 months ago) by chrisfen
File size: 45753 byte(s)
Log Message:
slowdown fixed - now roughly the same speed as the old version when using dipoles

energies are now exactly the same between the old version of OOPSE and this version

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     #ifdef IS_MPI
51     use mpiSimulation
52     #endif
53     implicit none
54    
55     PRIVATE
56    
57 gezelter 2301 #define __FORTRAN90
58     #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59    
60 gezelter 2118 !! these prefactors convert the multipole interactions into kcal / mol
61     !! all were computed assuming distances are measured in angstroms
62     !! Charge-Charge, assuming charges are measured in electrons
63 gezelter 2095 real(kind=dp), parameter :: pre11 = 332.0637778_dp
64 gezelter 2118 !! Charge-Dipole, assuming charges are measured in electrons, and
65     !! dipoles are measured in debyes
66     real(kind=dp), parameter :: pre12 = 69.13373_dp
67     !! Dipole-Dipole, assuming dipoles are measured in debyes
68     real(kind=dp), parameter :: pre22 = 14.39325_dp
69     !! Charge-Quadrupole, assuming charges are measured in electrons, and
70     !! quadrupoles are measured in 10^-26 esu cm^2
71     !! This unit is also known affectionately as an esu centi-barn.
72     real(kind=dp), parameter :: pre14 = 69.13373_dp
73 gezelter 2095
74 gezelter 2301 !! variables to handle different summation methods for long-range electrostatics:
75     integer, save :: summationMethod = NONE
76 chrisfen 2302 logical, save :: summationMethodChecked = .false.
77 gezelter 2301 real(kind=DP), save :: defaultCutoff = 0.0_DP
78     logical, save :: haveDefaultCutoff = .false.
79     real(kind=DP), save :: dampingAlpha = 0.0_DP
80     logical, save :: haveDampingAlpha = .false.
81     real(kind=DP), save :: dielectric = 0.0_DP
82     logical, save :: haveDielectric = .false.
83     real(kind=DP), save :: constERFC = 0.0_DP
84     real(kind=DP), save :: constEXP = 0.0_DP
85     logical, save :: haveDWAconstants = .false.
86 chrisfen 2306 real(kind=dp), save :: rcuti = 0.0_dp
87     real(kind=dp), save :: rcuti2 = 0.0_dp
88     real(kind=dp), save :: rcuti3 = 0.0_dp
89     real(kind=dp), save :: rcuti4 = 0.0_dp
90 chrisfen 2325 logical, save :: is_Undamped_Wolf = .false.
91     logical, save :: is_Damped_Wolf = .false.
92 gezelter 2301
93    
94     public :: setElectrostaticSummationMethod
95     public :: setElectrostaticCutoffRadius
96     public :: setDampedWolfAlpha
97     public :: setReactionFieldDielectric
98 gezelter 2095 public :: newElectrostaticType
99     public :: setCharge
100     public :: setDipoleMoment
101     public :: setSplitDipoleDistance
102     public :: setQuadrupoleMoments
103     public :: doElectrostaticPair
104     public :: getCharge
105     public :: getDipoleMoment
106 chrisfen 2129 public :: pre22
107 chuckv 2189 public :: destroyElectrostaticTypes
108 gezelter 2095
109     type :: Electrostatic
110     integer :: c_ident
111     logical :: is_Charge = .false.
112     logical :: is_Dipole = .false.
113     logical :: is_SplitDipole = .false.
114     logical :: is_Quadrupole = .false.
115 chrisfen 2229 logical :: is_Tap = .false.
116 gezelter 2095 real(kind=DP) :: charge = 0.0_DP
117     real(kind=DP) :: dipole_moment = 0.0_DP
118     real(kind=DP) :: split_dipole_distance = 0.0_DP
119     real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP
120     end type Electrostatic
121    
122     type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
123    
124     contains
125    
126 gezelter 2301 subroutine setElectrostaticSummationMethod(the_ESM)
127     integer, intent(in) :: the_ESM
128    
129     if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then
130     call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method")
131     endif
132    
133 chrisfen 2309 summationMethod = the_ESM
134 chrisfen 2325
135     if (summationMethod .eq. UNDAMPED_WOLF) is_Undamped_Wolf = .true.
136     if (summationMethod .eq. DAMPED_WOLF) is_Damped_Wolf = .true.
137 gezelter 2301 end subroutine setElectrostaticSummationMethod
138    
139     subroutine setElectrostaticCutoffRadius(thisRcut)
140     real(kind=dp), intent(in) :: thisRcut
141     defaultCutoff = thisRcut
142     haveDefaultCutoff = .true.
143     end subroutine setElectrostaticCutoffRadius
144    
145     subroutine setDampedWolfAlpha(thisAlpha)
146     real(kind=dp), intent(in) :: thisAlpha
147     dampingAlpha = thisAlpha
148     haveDampingAlpha = .true.
149     end subroutine setDampedWolfAlpha
150    
151     subroutine setReactionFieldDielectric(thisDielectric)
152     real(kind=dp), intent(in) :: thisDielectric
153     dielectric = thisDielectric
154     haveDielectric = .true.
155     end subroutine setReactionFieldDielectric
156    
157 gezelter 2095 subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
158 chrisfen 2229 is_SplitDipole, is_Quadrupole, is_Tap, status)
159 gezelter 2204
160 gezelter 2095 integer, intent(in) :: c_ident
161     logical, intent(in) :: is_Charge
162     logical, intent(in) :: is_Dipole
163     logical, intent(in) :: is_SplitDipole
164     logical, intent(in) :: is_Quadrupole
165 chrisfen 2229 logical, intent(in) :: is_Tap
166 gezelter 2095 integer, intent(out) :: status
167     integer :: nAtypes, myATID, i, j
168    
169     status = 0
170     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
171 gezelter 2204
172 gezelter 2095 !! Be simple-minded and assume that we need an ElectrostaticMap that
173     !! is the same size as the total number of atom types
174    
175     if (.not.allocated(ElectrostaticMap)) then
176 gezelter 2204
177 gezelter 2095 nAtypes = getSize(atypes)
178 gezelter 2204
179 gezelter 2095 if (nAtypes == 0) then
180     status = -1
181     return
182     end if
183 gezelter 2204
184 gezelter 2095 if (.not. allocated(ElectrostaticMap)) then
185     allocate(ElectrostaticMap(nAtypes))
186     endif
187 gezelter 2204
188 gezelter 2095 end if
189    
190     if (myATID .gt. size(ElectrostaticMap)) then
191     status = -1
192     return
193     endif
194 gezelter 2204
195 gezelter 2095 ! set the values for ElectrostaticMap for this atom type:
196    
197     ElectrostaticMap(myATID)%c_ident = c_ident
198     ElectrostaticMap(myATID)%is_Charge = is_Charge
199     ElectrostaticMap(myATID)%is_Dipole = is_Dipole
200     ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole
201     ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole
202 chrisfen 2229 ElectrostaticMap(myATID)%is_Tap = is_Tap
203 gezelter 2204
204 gezelter 2095 end subroutine newElectrostaticType
205    
206     subroutine setCharge(c_ident, charge, status)
207     integer, intent(in) :: c_ident
208     real(kind=dp), intent(in) :: charge
209     integer, intent(out) :: status
210     integer :: myATID
211    
212     status = 0
213     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
214    
215     if (.not.allocated(ElectrostaticMap)) then
216     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!")
217     status = -1
218     return
219     end if
220    
221     if (myATID .gt. size(ElectrostaticMap)) then
222     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!")
223     status = -1
224     return
225     endif
226    
227     if (.not.ElectrostaticMap(myATID)%is_Charge) then
228     call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!")
229     status = -1
230     return
231 gezelter 2204 endif
232 gezelter 2095
233     ElectrostaticMap(myATID)%charge = charge
234     end subroutine setCharge
235    
236     subroutine setDipoleMoment(c_ident, dipole_moment, status)
237     integer, intent(in) :: c_ident
238     real(kind=dp), intent(in) :: dipole_moment
239     integer, intent(out) :: status
240     integer :: myATID
241    
242     status = 0
243     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
244    
245     if (.not.allocated(ElectrostaticMap)) then
246     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!")
247     status = -1
248     return
249     end if
250    
251     if (myATID .gt. size(ElectrostaticMap)) then
252     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!")
253     status = -1
254     return
255     endif
256    
257     if (.not.ElectrostaticMap(myATID)%is_Dipole) then
258     call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!")
259     status = -1
260     return
261     endif
262    
263     ElectrostaticMap(myATID)%dipole_moment = dipole_moment
264     end subroutine setDipoleMoment
265    
266     subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status)
267     integer, intent(in) :: c_ident
268     real(kind=dp), intent(in) :: split_dipole_distance
269     integer, intent(out) :: status
270     integer :: myATID
271    
272     status = 0
273     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
274    
275     if (.not.allocated(ElectrostaticMap)) then
276     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!")
277     status = -1
278     return
279     end if
280    
281     if (myATID .gt. size(ElectrostaticMap)) then
282     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!")
283     status = -1
284     return
285     endif
286    
287     if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then
288     call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!")
289     status = -1
290     return
291     endif
292    
293     ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance
294     end subroutine setSplitDipoleDistance
295    
296     subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status)
297     integer, intent(in) :: c_ident
298     real(kind=dp), intent(in), dimension(3) :: quadrupole_moments
299     integer, intent(out) :: status
300     integer :: myATID, i, j
301    
302     status = 0
303     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
304    
305     if (.not.allocated(ElectrostaticMap)) then
306     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!")
307     status = -1
308     return
309     end if
310    
311     if (myATID .gt. size(ElectrostaticMap)) then
312     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!")
313     status = -1
314     return
315     endif
316    
317     if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then
318     call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!")
319     status = -1
320     return
321     endif
322 gezelter 2204
323 gezelter 2095 do i = 1, 3
324 gezelter 2204 ElectrostaticMap(myATID)%quadrupole_moments(i) = &
325     quadrupole_moments(i)
326     enddo
327 gezelter 2095
328     end subroutine setQuadrupoleMoments
329    
330 gezelter 2204
331 gezelter 2095 function getCharge(atid) result (c)
332     integer, intent(in) :: atid
333     integer :: localError
334     real(kind=dp) :: c
335 gezelter 2204
336 gezelter 2095 if (.not.allocated(ElectrostaticMap)) then
337     call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
338     return
339     end if
340 gezelter 2204
341 gezelter 2095 if (.not.ElectrostaticMap(atid)%is_Charge) then
342     call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!")
343     return
344     endif
345 gezelter 2204
346 gezelter 2095 c = ElectrostaticMap(atid)%charge
347     end function getCharge
348    
349     function getDipoleMoment(atid) result (dm)
350     integer, intent(in) :: atid
351     integer :: localError
352     real(kind=dp) :: dm
353 gezelter 2204
354 gezelter 2095 if (.not.allocated(ElectrostaticMap)) then
355     call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
356     return
357     end if
358 gezelter 2204
359 gezelter 2095 if (.not.ElectrostaticMap(atid)%is_Dipole) then
360     call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!")
361     return
362     endif
363 gezelter 2204
364 gezelter 2095 dm = ElectrostaticMap(atid)%dipole_moment
365     end function getDipoleMoment
366    
367 gezelter 2301 subroutine checkSummationMethod()
368    
369 chrisfen 2306 if (.not.haveDefaultCutoff) then
370     call handleError("checkSummationMethod", "no Default Cutoff set!")
371     endif
372    
373     rcuti = 1.0d0 / defaultCutoff
374     rcuti2 = rcuti*rcuti
375     rcuti3 = rcuti2*rcuti
376     rcuti4 = rcuti2*rcuti2
377    
378 gezelter 2301 if (summationMethod .eq. DAMPED_WOLF) then
379     if (.not.haveDWAconstants) then
380    
381     if (.not.haveDampingAlpha) then
382     call handleError("checkSummationMethod", "no Damping Alpha set!")
383     endif
384    
385 chrisfen 2302 if (.not.haveDefaultCutoff) then
386     call handleError("checkSummationMethod", "no Default Cutoff set!")
387     endif
388    
389     constEXP = exp(-dampingAlpha*dampingAlpha*defaultCutoff*defaultCutoff)
390     constERFC = erfc(dampingAlpha*defaultCutoff)
391 gezelter 2301
392     haveDWAconstants = .true.
393     endif
394     endif
395    
396 chrisfen 2302 if (summationMethod .eq. REACTION_FIELD) then
397     if (.not.haveDielectric) then
398     call handleError("checkSummationMethod", "no reaction field Dielectric set!")
399     endif
400     endif
401    
402     summationMethodChecked = .true.
403 gezelter 2301 end subroutine checkSummationMethod
404    
405    
406    
407 gezelter 2095 subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, &
408 chrisfen 2306 vpair, fpair, pot, eFrame, f, t, do_pot)
409 gezelter 2204
410 gezelter 2095 logical, intent(in) :: do_pot
411 gezelter 2204
412 gezelter 2095 integer, intent(in) :: atom1, atom2
413     integer :: localError
414    
415 chrisfen 2306 real(kind=dp), intent(in) :: rij, r2, sw
416 gezelter 2095 real(kind=dp), intent(in), dimension(3) :: d
417     real(kind=dp), intent(inout) :: vpair
418     real(kind=dp), intent(inout), dimension(3) :: fpair
419    
420 chrisfen 2325 real( kind = dp ) :: pot
421 gezelter 2095 real( kind = dp ), dimension(9,nLocal) :: eFrame
422     real( kind = dp ), dimension(3,nLocal) :: f
423     real( kind = dp ), dimension(3,nLocal) :: t
424 gezelter 2204
425 gezelter 2123 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
426     real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
427     real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
428     real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
429 gezelter 2095
430     logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
431     logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
432 chrisfen 2229 logical :: i_is_Tap, j_is_Tap
433 gezelter 2095 integer :: me1, me2, id1, id2
434     real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
435 gezelter 2123 real (kind=dp) :: qxx_i, qyy_i, qzz_i
436     real (kind=dp) :: qxx_j, qyy_j, qzz_j
437     real (kind=dp) :: cx_i, cy_i, cz_i
438     real (kind=dp) :: cx_j, cy_j, cz_j
439     real (kind=dp) :: cx2, cy2, cz2
440 gezelter 2095 real (kind=dp) :: ct_i, ct_j, ct_ij, a1
441 gezelter 2105 real (kind=dp) :: riji, ri, ri2, ri3, ri4
442 chrisfen 2296 real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2
443 gezelter 2105 real (kind=dp) :: xhat, yhat, zhat
444 gezelter 2095 real (kind=dp) :: dudx, dudy, dudz
445 chrisfen 2325 real (kind=dp) :: scale, sc2, bigR
446 gezelter 2095
447     if (.not.allocated(ElectrostaticMap)) then
448     call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
449     return
450     end if
451    
452 gezelter 2301 if (.not.summationMethodChecked) then
453     call checkSummationMethod()
454 chrisfen 2310
455 gezelter 2301 endif
456    
457    
458 gezelter 2095 #ifdef IS_MPI
459     me1 = atid_Row(atom1)
460     me2 = atid_Col(atom2)
461     #else
462     me1 = atid(atom1)
463     me2 = atid(atom2)
464     #endif
465    
466     !! some variables we'll need independent of electrostatic type:
467    
468     riji = 1.0d0 / rij
469    
470 gezelter 2105 xhat = d(1) * riji
471     yhat = d(2) * riji
472     zhat = d(3) * riji
473 gezelter 2095
474     !! logicals
475     i_is_Charge = ElectrostaticMap(me1)%is_Charge
476     i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
477     i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
478     i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
479 chrisfen 2229 i_is_Tap = ElectrostaticMap(me1)%is_Tap
480 gezelter 2095
481     j_is_Charge = ElectrostaticMap(me2)%is_Charge
482     j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
483     j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
484     j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
485 chrisfen 2229 j_is_Tap = ElectrostaticMap(me2)%is_Tap
486 gezelter 2095
487     if (i_is_Charge) then
488     q_i = ElectrostaticMap(me1)%charge
489     endif
490 gezelter 2204
491 gezelter 2095 if (i_is_Dipole) then
492     mu_i = ElectrostaticMap(me1)%dipole_moment
493     #ifdef IS_MPI
494 gezelter 2123 uz_i(1) = eFrame_Row(3,atom1)
495     uz_i(2) = eFrame_Row(6,atom1)
496     uz_i(3) = eFrame_Row(9,atom1)
497 gezelter 2095 #else
498 gezelter 2123 uz_i(1) = eFrame(3,atom1)
499     uz_i(2) = eFrame(6,atom1)
500     uz_i(3) = eFrame(9,atom1)
501 gezelter 2095 #endif
502 gezelter 2123 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
503 gezelter 2095
504     if (i_is_SplitDipole) then
505     d_i = ElectrostaticMap(me1)%split_dipole_distance
506     endif
507 gezelter 2204
508 gezelter 2095 endif
509    
510 gezelter 2123 if (i_is_Quadrupole) then
511     qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
512     qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
513     qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
514     #ifdef IS_MPI
515     ux_i(1) = eFrame_Row(1,atom1)
516     ux_i(2) = eFrame_Row(4,atom1)
517     ux_i(3) = eFrame_Row(7,atom1)
518     uy_i(1) = eFrame_Row(2,atom1)
519     uy_i(2) = eFrame_Row(5,atom1)
520     uy_i(3) = eFrame_Row(8,atom1)
521     uz_i(1) = eFrame_Row(3,atom1)
522     uz_i(2) = eFrame_Row(6,atom1)
523     uz_i(3) = eFrame_Row(9,atom1)
524     #else
525     ux_i(1) = eFrame(1,atom1)
526     ux_i(2) = eFrame(4,atom1)
527     ux_i(3) = eFrame(7,atom1)
528     uy_i(1) = eFrame(2,atom1)
529     uy_i(2) = eFrame(5,atom1)
530     uy_i(3) = eFrame(8,atom1)
531     uz_i(1) = eFrame(3,atom1)
532     uz_i(2) = eFrame(6,atom1)
533     uz_i(3) = eFrame(9,atom1)
534     #endif
535     cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
536     cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
537     cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
538     endif
539    
540 gezelter 2095 if (j_is_Charge) then
541     q_j = ElectrostaticMap(me2)%charge
542     endif
543 gezelter 2204
544 gezelter 2095 if (j_is_Dipole) then
545     mu_j = ElectrostaticMap(me2)%dipole_moment
546     #ifdef IS_MPI
547 gezelter 2123 uz_j(1) = eFrame_Col(3,atom2)
548     uz_j(2) = eFrame_Col(6,atom2)
549     uz_j(3) = eFrame_Col(9,atom2)
550 gezelter 2095 #else
551 gezelter 2123 uz_j(1) = eFrame(3,atom2)
552     uz_j(2) = eFrame(6,atom2)
553     uz_j(3) = eFrame(9,atom2)
554 gezelter 2095 #endif
555 chrisfen 2162 ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
556 gezelter 2095
557     if (j_is_SplitDipole) then
558     d_j = ElectrostaticMap(me2)%split_dipole_distance
559     endif
560     endif
561    
562 gezelter 2123 if (j_is_Quadrupole) then
563     qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
564     qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
565     qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
566     #ifdef IS_MPI
567     ux_j(1) = eFrame_Col(1,atom2)
568     ux_j(2) = eFrame_Col(4,atom2)
569     ux_j(3) = eFrame_Col(7,atom2)
570     uy_j(1) = eFrame_Col(2,atom2)
571     uy_j(2) = eFrame_Col(5,atom2)
572     uy_j(3) = eFrame_Col(8,atom2)
573     uz_j(1) = eFrame_Col(3,atom2)
574     uz_j(2) = eFrame_Col(6,atom2)
575     uz_j(3) = eFrame_Col(9,atom2)
576     #else
577     ux_j(1) = eFrame(1,atom2)
578     ux_j(2) = eFrame(4,atom2)
579     ux_j(3) = eFrame(7,atom2)
580     uy_j(1) = eFrame(2,atom2)
581     uy_j(2) = eFrame(5,atom2)
582     uy_j(3) = eFrame(8,atom2)
583     uz_j(1) = eFrame(3,atom2)
584     uz_j(2) = eFrame(6,atom2)
585     uz_j(3) = eFrame(9,atom2)
586     #endif
587     cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
588     cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
589     cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
590     endif
591 chrisfen 2251
592 gezelter 2095 epot = 0.0_dp
593     dudx = 0.0_dp
594     dudy = 0.0_dp
595     dudz = 0.0_dp
596    
597 gezelter 2123 dudux_i = 0.0_dp
598     duduy_i = 0.0_dp
599     duduz_i = 0.0_dp
600 gezelter 2095
601 gezelter 2123 dudux_j = 0.0_dp
602     duduy_j = 0.0_dp
603     duduz_j = 0.0_dp
604 gezelter 2095
605     if (i_is_Charge) then
606    
607     if (j_is_Charge) then
608 gezelter 2204
609 chrisfen 2310 if (summationMethod .eq. UNDAMPED_WOLF) then
610 chrisfen 2325
611 chrisfen 2296 vterm = pre11 * q_i * q_j * (riji - rcuti)
612     vpair = vpair + vterm
613 chrisfen 2325 epot = epot + sw*vterm
614 chrisfen 2296
615     dudr = - sw * pre11 * q_i * q_j * (riji*riji*riji - rcuti2*rcuti)
616    
617     dudx = dudx + dudr * d(1)
618     dudy = dudy + dudr * d(2)
619     dudz = dudz + dudr * d(3)
620 gezelter 2095
621 chrisfen 2296 else
622 chrisfen 2325
623 chrisfen 2296 vterm = pre11 * q_i * q_j * riji
624     vpair = vpair + vterm
625 chrisfen 2325 epot = epot + sw*vterm
626 chrisfen 2296
627     dudr = - sw * vterm * riji
628    
629     dudx = dudx + dudr * xhat
630     dudy = dudy + dudr * yhat
631     dudz = dudz + dudr * zhat
632    
633     endif
634    
635 gezelter 2095 endif
636    
637     if (j_is_Dipole) then
638    
639 chrisfen 2325 pref = pre12 * q_i * mu_j
640 gezelter 2095
641 chrisfen 2310 if (summationMethod .eq. UNDAMPED_WOLF) then
642 chrisfen 2296 ri2 = riji * riji
643     ri3 = ri2 * riji
644 gezelter 2204
645 chrisfen 2325 pref = pre12 * q_i * mu_j
646 chrisfen 2296 vterm = - pref * ct_j * (ri2 - rcuti2)
647 chrisfen 2325 vpair = vpair + vterm
648     epot = epot + sw*vterm
649 chrisfen 2296
650     !! this has a + sign in the () because the rij vector is
651     !! r_j - r_i and the charge-dipole potential takes the origin
652     !! as the point dipole, which is atom j in this case.
653    
654 chrisfen 2325 dudx = dudx - sw*pref * ( ri3*( uz_j(1) - 3.0d0*ct_j*xhat) &
655 chrisfen 2296 - rcuti3*( uz_j(1) - 3.0d0*ct_j*d(1)*rcuti ) )
656 chrisfen 2325 dudy = dudy - sw*pref * ( ri3*( uz_j(2) - 3.0d0*ct_j*yhat) &
657 chrisfen 2296 - rcuti3*( uz_j(2) - 3.0d0*ct_j*d(2)*rcuti ) )
658 chrisfen 2325 dudz = dudz - sw*pref * ( ri3*( uz_j(3) - 3.0d0*ct_j*zhat) &
659 chrisfen 2296 - rcuti3*( uz_j(3) - 3.0d0*ct_j*d(3)*rcuti ) )
660    
661 chrisfen 2325 duduz_j(1) = duduz_j(1) - sw*pref*( ri2*xhat - d(1)*rcuti3 )
662     duduz_j(2) = duduz_j(2) - sw*pref*( ri2*yhat - d(2)*rcuti3 )
663     duduz_j(3) = duduz_j(3) - sw*pref*( ri2*zhat - d(3)*rcuti3 )
664 gezelter 2095
665 chrisfen 2296 else
666     if (j_is_SplitDipole) then
667     BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
668     ri = 1.0_dp / BigR
669     scale = rij * ri
670     else
671     ri = riji
672     scale = 1.0_dp
673     endif
674    
675     ri2 = ri * ri
676     ri3 = ri2 * ri
677     sc2 = scale * scale
678 chrisfen 2325
679     pref = pre12 * q_i * mu_j
680 chrisfen 2296 vterm = - pref * ct_j * ri2 * scale
681 chrisfen 2325 vpair = vpair + vterm
682     epot = epot + sw*vterm
683 chrisfen 2296
684     !! this has a + sign in the () because the rij vector is
685     !! r_j - r_i and the charge-dipole potential takes the origin
686     !! as the point dipole, which is atom j in this case.
687    
688 chrisfen 2325 dudx = dudx - sw*pref * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2)
689     dudy = dudy - sw*pref * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2)
690     dudz = dudz - sw*pref * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2)
691 chrisfen 2296
692 chrisfen 2325 duduz_j(1) = duduz_j(1) - sw*pref * ri2 * xhat * scale
693     duduz_j(2) = duduz_j(2) - sw*pref * ri2 * yhat * scale
694     duduz_j(3) = duduz_j(3) - sw*pref * ri2 * zhat * scale
695 gezelter 2095
696 chrisfen 2296 endif
697 gezelter 2095 endif
698 gezelter 2105
699 gezelter 2123 if (j_is_Quadrupole) then
700     ri2 = riji * riji
701     ri3 = ri2 * riji
702 gezelter 2124 ri4 = ri2 * ri2
703 gezelter 2123 cx2 = cx_j * cx_j
704     cy2 = cy_j * cy_j
705     cz2 = cz_j * cz_j
706    
707 chrisfen 2310 if (summationMethod .eq. UNDAMPED_WOLF) then
708 chrisfen 2325 pref = pre14 * q_i / 3.0_dp
709 chrisfen 2296 vterm1 = pref * ri3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
710     qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
711     qzz_j * (3.0_dp*cz2 - 1.0_dp) )
712     vterm2 = pref * rcuti3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
713     qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
714     qzz_j * (3.0_dp*cz2 - 1.0_dp) )
715 chrisfen 2325 vpair = vpair + ( vterm1 - vterm2 )
716     epot = epot + sw*( vterm1 - vterm2 )
717 chrisfen 2296
718     dudx = dudx - (5.0_dp * &
719 chrisfen 2325 (vterm1*riji*xhat - vterm2*rcuti2*d(1))) + sw*pref * ( &
720 chrisfen 2296 (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(1)) - &
721     qxx_j*2.0_dp*(xhat - rcuti*d(1))) + &
722     (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(1)) - &
723     qyy_j*2.0_dp*(xhat - rcuti*d(1))) + &
724     (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(1)) - &
725     qzz_j*2.0_dp*(xhat - rcuti*d(1))) )
726     dudy = dudy - (5.0_dp * &
727 chrisfen 2325 (vterm1*riji*yhat - vterm2*rcuti2*d(2))) + sw*pref * ( &
728 chrisfen 2296 (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(2)) - &
729     qxx_j*2.0_dp*(yhat - rcuti*d(2))) + &
730     (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(2)) - &
731     qyy_j*2.0_dp*(yhat - rcuti*d(2))) + &
732     (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(2)) - &
733     qzz_j*2.0_dp*(yhat - rcuti*d(2))) )
734     dudz = dudz - (5.0_dp * &
735 chrisfen 2325 (vterm1*riji*zhat - vterm2*rcuti2*d(3))) + sw*pref * ( &
736 chrisfen 2296 (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(3)) - &
737     qxx_j*2.0_dp*(zhat - rcuti*d(3))) + &
738     (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(3)) - &
739     qyy_j*2.0_dp*(zhat - rcuti*d(3))) + &
740     (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(3)) - &
741     qzz_j*2.0_dp*(zhat - rcuti*d(3))) )
742    
743 chrisfen 2325 dudux_j(1) = dudux_j(1) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*xhat) -&
744 chrisfen 2296 rcuti4*(qxx_j*6.0_dp*cx_j*d(1)))
745 chrisfen 2325 dudux_j(2) = dudux_j(2) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*yhat) -&
746 chrisfen 2296 rcuti4*(qxx_j*6.0_dp*cx_j*d(2)))
747 chrisfen 2325 dudux_j(3) = dudux_j(3) + sw*pref*(ri3*(qxx_j*6.0_dp*cx_j*zhat) -&
748 chrisfen 2296 rcuti4*(qxx_j*6.0_dp*cx_j*d(3)))
749    
750 chrisfen 2325 duduy_j(1) = duduy_j(1) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*xhat) -&
751 chrisfen 2296 rcuti4*(qyy_j*6.0_dp*cx_j*d(1)))
752 chrisfen 2325 duduy_j(2) = duduy_j(2) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*yhat) -&
753 chrisfen 2296 rcuti4*(qyy_j*6.0_dp*cx_j*d(2)))
754 chrisfen 2325 duduy_j(3) = duduy_j(3) + sw*pref*(ri3*(qyy_j*6.0_dp*cy_j*zhat) -&
755 chrisfen 2296 rcuti4*(qyy_j*6.0_dp*cx_j*d(3)))
756    
757 chrisfen 2325 duduz_j(1) = duduz_j(1) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*xhat) -&
758 chrisfen 2296 rcuti4*(qzz_j*6.0_dp*cx_j*d(1)))
759 chrisfen 2325 duduz_j(2) = duduz_j(2) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*yhat) -&
760 chrisfen 2296 rcuti4*(qzz_j*6.0_dp*cx_j*d(2)))
761 chrisfen 2325 duduz_j(3) = duduz_j(3) + sw*pref*(ri3*(qzz_j*6.0_dp*cz_j*zhat) -&
762 chrisfen 2296 rcuti4*(qzz_j*6.0_dp*cx_j*d(3)))
763    
764     else
765 chrisfen 2325 pref = pre14 * q_i / 3.0_dp
766 chrisfen 2296 vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
767     qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
768     qzz_j * (3.0_dp*cz2 - 1.0_dp))
769 chrisfen 2325 vpair = vpair + vterm
770     epot = epot + sw*vterm
771 chrisfen 2296
772 chrisfen 2325 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + sw*pref * ri4 * ( &
773 chrisfen 2296 qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
774     qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
775     qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
776 chrisfen 2325 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + sw*pref * ri4 * ( &
777 chrisfen 2296 qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
778     qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
779     qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
780 chrisfen 2325 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + sw*pref * ri4 * ( &
781 chrisfen 2296 qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
782     qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
783     qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
784    
785 chrisfen 2325 dudux_j(1) = dudux_j(1) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*xhat)
786     dudux_j(2) = dudux_j(2) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*yhat)
787     dudux_j(3) = dudux_j(3) + sw*pref * ri3*(qxx_j*6.0_dp*cx_j*zhat)
788 chrisfen 2296
789 chrisfen 2325 duduy_j(1) = duduy_j(1) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*xhat)
790     duduy_j(2) = duduy_j(2) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*yhat)
791     duduy_j(3) = duduy_j(3) + sw*pref * ri3*(qyy_j*6.0_dp*cy_j*zhat)
792 chrisfen 2296
793 chrisfen 2325 duduz_j(1) = duduz_j(1) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*xhat)
794     duduz_j(2) = duduz_j(2) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*yhat)
795     duduz_j(3) = duduz_j(3) + sw*pref * ri3*(qzz_j*6.0_dp*cz_j*zhat)
796 chrisfen 2296
797     endif
798 gezelter 2123 endif
799 gezelter 2095 endif
800 gezelter 2204
801 gezelter 2095 if (i_is_Dipole) then
802 gezelter 2204
803 gezelter 2095 if (j_is_Charge) then
804 chrisfen 2325
805     pref = pre12 * q_j * mu_i
806    
807 chrisfen 2310 if (summationMethod .eq. UNDAMPED_WOLF) then
808 chrisfen 2296 ri2 = riji * riji
809     ri3 = ri2 * riji
810 gezelter 2204
811 chrisfen 2325 pref = pre12 * q_j * mu_i
812 chrisfen 2296 vterm = pref * ct_i * (ri2 - rcuti2)
813 chrisfen 2325 vpair = vpair + vterm
814     epot = epot + sw*vterm
815 chrisfen 2296
816     !! this has a + sign in the () because the rij vector is
817     !! r_j - r_i and the charge-dipole potential takes the origin
818     !! as the point dipole, which is atom j in this case.
819    
820 chrisfen 2325 dudx = dudx + sw*pref * ( ri3*( uz_i(1) - 3.0d0*ct_i*xhat) &
821 chrisfen 2296 - rcuti3*( uz_i(1) - 3.0d0*ct_i*d(1)*rcuti ) )
822 chrisfen 2325 dudy = dudy + sw*pref * ( ri3*( uz_i(2) - 3.0d0*ct_i*yhat) &
823 chrisfen 2296 - rcuti3*( uz_i(2) - 3.0d0*ct_i*d(2)*rcuti ) )
824 chrisfen 2325 dudz = dudz + sw*pref * ( ri3*( uz_i(3) - 3.0d0*ct_i*zhat) &
825 chrisfen 2296 - rcuti3*( uz_i(3) - 3.0d0*ct_i*d(3)*rcuti ) )
826    
827 chrisfen 2325 duduz_i(1) = duduz_i(1) - sw*pref*( ri2*xhat - d(1)*rcuti3 )
828     duduz_i(2) = duduz_i(2) - sw*pref*( ri2*yhat - d(2)*rcuti3 )
829     duduz_i(3) = duduz_i(3) - sw*pref*( ri2*zhat - d(3)*rcuti3 )
830 gezelter 2095
831 chrisfen 2296 else
832     if (i_is_SplitDipole) then
833 gezelter 2105 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
834     ri = 1.0_dp / BigR
835 chrisfen 2296 scale = rij * ri
836     else
837 gezelter 2105 ri = riji
838     scale = 1.0_dp
839     endif
840 chrisfen 2296
841     ri2 = ri * ri
842     ri3 = ri2 * ri
843     sc2 = scale * scale
844 chrisfen 2325
845     pref = pre12 * q_j * mu_i
846 chrisfen 2296 vterm = pref * ct_i * ri2 * scale
847 chrisfen 2325 vpair = vpair + vterm
848     epot = epot + sw*vterm
849 chrisfen 2296
850 chrisfen 2325 dudx = dudx + sw*pref * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
851     dudy = dudy + sw*pref * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
852     dudz = dudz + sw*pref * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
853 chrisfen 2296
854 chrisfen 2325 duduz_i(1) = duduz_i(1) + sw*pref * ri2 * xhat * scale
855     duduz_i(2) = duduz_i(2) + sw*pref * ri2 * yhat * scale
856     duduz_i(3) = duduz_i(3) + sw*pref * ri2 * zhat * scale
857 gezelter 2105 endif
858 chrisfen 2296 endif
859 chrisfen 2325
860 chrisfen 2296 if (j_is_Dipole) then
861 gezelter 2105
862 chrisfen 2310 if (summationMethod .eq. UNDAMPED_WOLF) then
863 chrisfen 2296 ri2 = riji * riji
864     ri3 = ri2 * riji
865     ri4 = ri2 * ri2
866 gezelter 2204
867 chrisfen 2325 pref = pre22 * mu_i * mu_j
868 chrisfen 2296 vterm = pref * (ri3 - rcuti3) * (ct_ij - 3.0d0 * ct_i * ct_j)
869 chrisfen 2325 vpair = vpair + vterm
870     epot = epot + sw*vterm
871 chrisfen 2296
872     a1 = 5.0d0 * ct_i * ct_j - ct_ij
873    
874 chrisfen 2325 dudx = dudx + sw*pref*3.0d0*ri4 &
875     * (a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) &
876     - sw*pref*3.0d0*rcuti4 &
877     * (a1*rcuti*d(1)-ct_i*uz_j(1)-ct_j*uz_i(1))
878     dudy = dudy + sw*pref*3.0d0*ri4 &
879     * (a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) &
880     - sw*pref*3.0d0*rcuti4 &
881     * (a1*rcuti*d(2)-ct_i*uz_j(2)-ct_j*uz_i(2))
882     dudz = dudz + sw*pref*3.0d0*ri4 &
883     * (a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) &
884     - sw*pref*3.0d0*rcuti4 &
885     * (a1*rcuti*d(3)-ct_i*uz_j(3)-ct_j*uz_i(3))
886 chrisfen 2296
887 chrisfen 2325 duduz_i(1) = duduz_i(1) + sw*pref*(ri3*(uz_j(1)-3.0d0*ct_j*xhat) &
888 chrisfen 2296 - rcuti3*(uz_j(1) - 3.0d0*ct_j*d(1)*rcuti))
889 chrisfen 2325 duduz_i(2) = duduz_i(2) + sw*pref*(ri3*(uz_j(2)-3.0d0*ct_j*yhat) &
890 chrisfen 2296 - rcuti3*(uz_j(2) - 3.0d0*ct_j*d(2)*rcuti))
891 chrisfen 2325 duduz_i(3) = duduz_i(3) + sw*pref*(ri3*(uz_j(3)-3.0d0*ct_j*zhat) &
892 chrisfen 2296 - rcuti3*(uz_j(3) - 3.0d0*ct_j*d(3)*rcuti))
893 chrisfen 2325 duduz_j(1) = duduz_j(1) + sw*pref*(ri3*(uz_i(1)-3.0d0*ct_i*xhat) &
894 chrisfen 2296 - rcuti3*(uz_i(1) - 3.0d0*ct_i*d(1)*rcuti))
895 chrisfen 2325 duduz_j(2) = duduz_j(2) + sw*pref*(ri3*(uz_i(2)-3.0d0*ct_i*yhat) &
896 chrisfen 2296 - rcuti3*(uz_i(2) - 3.0d0*ct_i*d(2)*rcuti))
897 chrisfen 2325 duduz_j(3) = duduz_j(3) + sw*pref*(ri3*(uz_i(3)-3.0d0*ct_i*zhat) &
898 chrisfen 2296 - rcuti3*(uz_i(3) - 3.0d0*ct_i*d(3)*rcuti))
899 chrisfen 2325
900 chrisfen 2296 else
901     if (i_is_SplitDipole) then
902     if (j_is_SplitDipole) then
903     BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
904     else
905     BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
906     endif
907     ri = 1.0_dp / BigR
908     scale = rij * ri
909     else
910     if (j_is_SplitDipole) then
911     BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
912     ri = 1.0_dp / BigR
913     scale = rij * ri
914     else
915     ri = riji
916     scale = 1.0_dp
917     endif
918     endif
919    
920     ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
921    
922     ri2 = ri * ri
923     ri3 = ri2 * ri
924     ri4 = ri2 * ri2
925     sc2 = scale * scale
926    
927 chrisfen 2325 pref = pre22 * mu_i * mu_j
928 chrisfen 2296 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
929 chrisfen 2325 vpair = vpair + vterm
930     epot = epot + sw*vterm
931 chrisfen 2296
932     a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
933    
934 chrisfen 2325 dudx = dudx + sw*pref*3.0d0*ri4*scale &
935     *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
936     dudy = dudy + sw*pref*3.0d0*ri4*scale &
937     *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
938     dudz = dudz + sw*pref*3.0d0*ri4*scale &
939     *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
940 chrisfen 2296
941 chrisfen 2325 duduz_i(1) = duduz_i(1) + sw*pref*ri3 &
942     *(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
943     duduz_i(2) = duduz_i(2) + sw*pref*ri3 &
944     *(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
945     duduz_i(3) = duduz_i(3) + sw*pref*ri3 &
946     *(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
947 chrisfen 2296
948 chrisfen 2325 duduz_j(1) = duduz_j(1) + sw*pref*ri3 &
949     *(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
950     duduz_j(2) = duduz_j(2) + sw*pref*ri3 &
951     *(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
952     duduz_j(3) = duduz_j(3) + sw*pref*ri3 &
953     *(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
954 chrisfen 2296 endif
955 gezelter 2095 endif
956     endif
957 gezelter 2123
958     if (i_is_Quadrupole) then
959     if (j_is_Charge) then
960 gezelter 2204
961 gezelter 2123 ri2 = riji * riji
962     ri3 = ri2 * riji
963 gezelter 2124 ri4 = ri2 * ri2
964 gezelter 2123 cx2 = cx_i * cx_i
965     cy2 = cy_i * cy_i
966     cz2 = cz_i * cz_i
967 gezelter 2204
968 chrisfen 2310 if (summationMethod .eq. UNDAMPED_WOLF) then
969 chrisfen 2325 pref = pre14 * q_j / 3.0_dp
970 chrisfen 2296 vterm1 = pref * ri3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
971     qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
972     qzz_i * (3.0_dp*cz2 - 1.0_dp) )
973     vterm2 = pref * rcuti3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
974     qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
975     qzz_i * (3.0_dp*cz2 - 1.0_dp) )
976 chrisfen 2325 vpair = vpair + ( vterm1 - vterm2 )
977     epot = epot + sw*( vterm1 - vterm2 )
978 chrisfen 2296
979 chrisfen 2325 dudx = dudx - sw*(5.0_dp*(vterm1*riji*xhat-vterm2*rcuti2*d(1))) +&
980     sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(1)) - &
981 chrisfen 2296 qxx_i*2.0_dp*(xhat - rcuti*d(1))) + &
982     (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(1)) - &
983     qyy_i*2.0_dp*(xhat - rcuti*d(1))) + &
984     (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(1)) - &
985     qzz_i*2.0_dp*(xhat - rcuti*d(1))) )
986 chrisfen 2325 dudy = dudy - sw*(5.0_dp*(vterm1*riji*yhat-vterm2*rcuti2*d(2))) +&
987     sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(2)) - &
988 chrisfen 2296 qxx_i*2.0_dp*(yhat - rcuti*d(2))) + &
989     (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(2)) - &
990     qyy_i*2.0_dp*(yhat - rcuti*d(2))) + &
991     (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(2)) - &
992     qzz_i*2.0_dp*(yhat - rcuti*d(2))) )
993 chrisfen 2325 dudz = dudz - sw*(5.0_dp*(vterm1*riji*zhat-vterm2*rcuti2*d(3))) +&
994     sw*pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(3)) - &
995 chrisfen 2296 qxx_i*2.0_dp*(zhat - rcuti*d(3))) + &
996     (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(3)) - &
997     qyy_i*2.0_dp*(zhat - rcuti*d(3))) + &
998     (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(3)) - &
999     qzz_i*2.0_dp*(zhat - rcuti*d(3))) )
1000    
1001 chrisfen 2325 dudux_i(1) = dudux_i(1) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*xhat) -&
1002 chrisfen 2296 rcuti4*(qxx_i*6.0_dp*cx_i*d(1)))
1003 chrisfen 2325 dudux_i(2) = dudux_i(2) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*yhat) -&
1004 chrisfen 2296 rcuti4*(qxx_i*6.0_dp*cx_i*d(2)))
1005 chrisfen 2325 dudux_i(3) = dudux_i(3) + sw*pref*(ri3*(qxx_i*6.0_dp*cx_i*zhat) -&
1006 chrisfen 2296 rcuti4*(qxx_i*6.0_dp*cx_i*d(3)))
1007    
1008 chrisfen 2325 duduy_i(1) = duduy_i(1) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*xhat) -&
1009 chrisfen 2296 rcuti4*(qyy_i*6.0_dp*cx_i*d(1)))
1010 chrisfen 2325 duduy_i(2) = duduy_i(2) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*yhat) -&
1011 chrisfen 2296 rcuti4*(qyy_i*6.0_dp*cx_i*d(2)))
1012 chrisfen 2325 duduy_i(3) = duduy_i(3) + sw*pref*(ri3*(qyy_i*6.0_dp*cy_i*zhat) -&
1013 chrisfen 2296 rcuti4*(qyy_i*6.0_dp*cx_i*d(3)))
1014    
1015 chrisfen 2325 duduz_i(1) = duduz_i(1) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*xhat) -&
1016 chrisfen 2296 rcuti4*(qzz_i*6.0_dp*cx_i*d(1)))
1017 chrisfen 2325 duduz_i(2) = duduz_i(2) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*yhat) -&
1018 chrisfen 2296 rcuti4*(qzz_i*6.0_dp*cx_i*d(2)))
1019 chrisfen 2325 duduz_i(3) = duduz_i(3) + sw*pref*(ri3*(qzz_i*6.0_dp*cz_i*zhat) -&
1020 chrisfen 2296 rcuti4*(qzz_i*6.0_dp*cx_i*d(3)))
1021 gezelter 2204
1022 chrisfen 2296 else
1023 chrisfen 2325 pref = pre14 * q_j / 3.0_dp
1024 chrisfen 2296 vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
1025     qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
1026     qzz_i * (3.0_dp*cz2 - 1.0_dp))
1027 chrisfen 2325 vpair = vpair + vterm
1028     epot = epot + sw*vterm
1029 chrisfen 2296
1030 chrisfen 2325 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + sw*pref*ri4 * ( &
1031 chrisfen 2296 qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
1032     qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
1033     qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
1034 chrisfen 2325 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + sw*pref*ri4 * ( &
1035 chrisfen 2296 qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
1036     qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
1037     qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
1038 chrisfen 2325 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + sw*pref*ri4 * ( &
1039 chrisfen 2296 qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
1040     qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
1041     qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
1042    
1043 chrisfen 2325 dudux_i(1) = dudux_i(1) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*xhat)
1044     dudux_i(2) = dudux_i(2) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*yhat)
1045     dudux_i(3) = dudux_i(3) + sw*pref*ri3*(qxx_i*6.0_dp*cx_i*zhat)
1046 chrisfen 2296
1047 chrisfen 2325 duduy_i(1) = duduy_i(1) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*xhat)
1048     duduy_i(2) = duduy_i(2) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*yhat)
1049     duduy_i(3) = duduy_i(3) + sw*pref*ri3*(qyy_i*6.0_dp*cy_i*zhat)
1050 chrisfen 2296
1051 chrisfen 2325 duduz_i(1) = duduz_i(1) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*xhat)
1052     duduz_i(2) = duduz_i(2) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*yhat)
1053     duduz_i(3) = duduz_i(3) + sw*pref*ri3*(qzz_i*6.0_dp*cz_i*zhat)
1054 chrisfen 2296 endif
1055 gezelter 2123 endif
1056     endif
1057 gezelter 2204
1058    
1059 gezelter 2095 if (do_pot) then
1060     #ifdef IS_MPI
1061     pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
1062     pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
1063     #else
1064     pot = pot + epot
1065     #endif
1066     endif
1067 gezelter 2204
1068 gezelter 2095 #ifdef IS_MPI
1069     f_Row(1,atom1) = f_Row(1,atom1) + dudx
1070     f_Row(2,atom1) = f_Row(2,atom1) + dudy
1071     f_Row(3,atom1) = f_Row(3,atom1) + dudz
1072 gezelter 2204
1073 gezelter 2095 f_Col(1,atom2) = f_Col(1,atom2) - dudx
1074     f_Col(2,atom2) = f_Col(2,atom2) - dudy
1075     f_Col(3,atom2) = f_Col(3,atom2) - dudz
1076 gezelter 2204
1077 gezelter 2095 if (i_is_Dipole .or. i_is_Quadrupole) then
1078 gezelter 2123 t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1079     t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1080     t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1081 gezelter 2095 endif
1082 gezelter 2123 if (i_is_Quadrupole) then
1083     t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
1084     t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
1085     t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
1086 gezelter 2095
1087 gezelter 2123 t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
1088     t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
1089     t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
1090     endif
1091    
1092 gezelter 2095 if (j_is_Dipole .or. j_is_Quadrupole) then
1093 gezelter 2123 t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1094     t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1095     t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1096 gezelter 2095 endif
1097 gezelter 2123 if (j_is_Quadrupole) then
1098     t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
1099     t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
1100     t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
1101 gezelter 2095
1102 gezelter 2123 t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
1103     t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
1104     t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
1105     endif
1106    
1107 gezelter 2095 #else
1108     f(1,atom1) = f(1,atom1) + dudx
1109     f(2,atom1) = f(2,atom1) + dudy
1110     f(3,atom1) = f(3,atom1) + dudz
1111 gezelter 2204
1112 gezelter 2095 f(1,atom2) = f(1,atom2) - dudx
1113     f(2,atom2) = f(2,atom2) - dudy
1114     f(3,atom2) = f(3,atom2) - dudz
1115 gezelter 2204
1116 gezelter 2095 if (i_is_Dipole .or. i_is_Quadrupole) then
1117 gezelter 2123 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1118     t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1119     t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1120 gezelter 2095 endif
1121 gezelter 2123 if (i_is_Quadrupole) then
1122     t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
1123     t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
1124     t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
1125    
1126     t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
1127     t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
1128     t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
1129     endif
1130    
1131 gezelter 2095 if (j_is_Dipole .or. j_is_Quadrupole) then
1132 gezelter 2123 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1133     t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1134     t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1135 gezelter 2095 endif
1136 gezelter 2123 if (j_is_Quadrupole) then
1137     t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
1138     t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
1139     t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
1140    
1141     t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
1142     t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
1143     t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
1144     endif
1145    
1146 gezelter 2095 #endif
1147 gezelter 2204
1148 gezelter 2095 #ifdef IS_MPI
1149     id1 = AtomRowToGlobal(atom1)
1150     id2 = AtomColToGlobal(atom2)
1151     #else
1152     id1 = atom1
1153     id2 = atom2
1154     #endif
1155    
1156     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
1157 gezelter 2204
1158 gezelter 2095 fpair(1) = fpair(1) + dudx
1159     fpair(2) = fpair(2) + dudy
1160     fpair(3) = fpair(3) + dudz
1161    
1162     endif
1163    
1164     return
1165     end subroutine doElectrostaticPair
1166 chuckv 2189
1167 chrisfen 2229 !! calculates the switching functions and their derivatives for a given
1168     subroutine calc_switch(r, mu, scale, dscale)
1169 gezelter 2204
1170 chrisfen 2229 real (kind=dp), intent(in) :: r, mu
1171     real (kind=dp), intent(inout) :: scale, dscale
1172     real (kind=dp) :: rl, ru, mulow, minRatio, temp, scaleVal
1173    
1174     ! distances must be in angstroms
1175     rl = 2.75d0
1176 chrisfen 2231 ru = 3.75d0
1177     mulow = 0.0d0 !3.3856d0 ! 1.84 * 1.84
1178 chrisfen 2229 minRatio = mulow / (mu*mu)
1179     scaleVal = 1.0d0 - minRatio
1180    
1181     if (r.lt.rl) then
1182     scale = minRatio
1183     dscale = 0.0d0
1184     elseif (r.gt.ru) then
1185     scale = 1.0d0
1186     dscale = 0.0d0
1187     else
1188     scale = 1.0d0 - scaleVal*((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) &
1189     / ((ru - rl)**3)
1190     dscale = -scaleVal * 6.0d0 * (r-ru)*(r-rl)/((ru - rl)**3)
1191     endif
1192    
1193     return
1194     end subroutine calc_switch
1195    
1196 chuckv 2189 subroutine destroyElectrostaticTypes()
1197    
1198 gezelter 2204 if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap)
1199    
1200 chuckv 2189 end subroutine destroyElectrostaticTypes
1201    
1202 gezelter 2095 end module electrostatic_module