ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2339
Committed: Wed Sep 28 18:47:17 2005 UTC (18 years, 11 months ago) by chrisfen
File size: 46573 byte(s)
Log Message:
working on wolf

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