ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2330
Committed: Mon Sep 26 17:07:54 2005 UTC (18 years, 11 months ago) by chuckv
File size: 45810 byte(s)
Log Message:
Changed erfc to derfc and declared it to be external to fix issure with ifc7. Hopefully this will not cause a problem with other compilers where derfc is an intrinsic function.

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