ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2127
Committed: Mon Mar 21 19:23:27 2005 UTC (19 years, 3 months ago) by gezelter
File size: 29075 byte(s)
Log Message:
constant back to correct value

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    
44     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 2118 !! these prefactors convert the multipole interactions into kcal / mol
58     !! all were computed assuming distances are measured in angstroms
59     !! Charge-Charge, assuming charges are measured in electrons
60 gezelter 2095 real(kind=dp), parameter :: pre11 = 332.0637778_dp
61 gezelter 2118 !! Charge-Dipole, assuming charges are measured in electrons, and
62     !! dipoles are measured in debyes
63     real(kind=dp), parameter :: pre12 = 69.13373_dp
64     !! Dipole-Dipole, assuming dipoles are measured in debyes
65     real(kind=dp), parameter :: pre22 = 14.39325_dp
66     !! Charge-Quadrupole, assuming charges are measured in electrons, and
67     !! quadrupoles are measured in 10^-26 esu cm^2
68     !! This unit is also known affectionately as an esu centi-barn.
69     real(kind=dp), parameter :: pre14 = 69.13373_dp
70 gezelter 2095
71     public :: newElectrostaticType
72     public :: setCharge
73     public :: setDipoleMoment
74     public :: setSplitDipoleDistance
75     public :: setQuadrupoleMoments
76     public :: doElectrostaticPair
77     public :: getCharge
78     public :: getDipoleMoment
79    
80     type :: Electrostatic
81     integer :: c_ident
82     logical :: is_Charge = .false.
83     logical :: is_Dipole = .false.
84     logical :: is_SplitDipole = .false.
85     logical :: is_Quadrupole = .false.
86     real(kind=DP) :: charge = 0.0_DP
87     real(kind=DP) :: dipole_moment = 0.0_DP
88     real(kind=DP) :: split_dipole_distance = 0.0_DP
89     real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP
90     end type Electrostatic
91    
92     type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
93    
94     contains
95    
96     subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
97     is_SplitDipole, is_Quadrupole, status)
98    
99     integer, intent(in) :: c_ident
100     logical, intent(in) :: is_Charge
101     logical, intent(in) :: is_Dipole
102     logical, intent(in) :: is_SplitDipole
103     logical, intent(in) :: is_Quadrupole
104     integer, intent(out) :: status
105     integer :: nAtypes, myATID, i, j
106    
107     status = 0
108     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
109    
110     !! Be simple-minded and assume that we need an ElectrostaticMap that
111     !! is the same size as the total number of atom types
112    
113     if (.not.allocated(ElectrostaticMap)) then
114    
115     nAtypes = getSize(atypes)
116    
117     if (nAtypes == 0) then
118     status = -1
119     return
120     end if
121    
122     if (.not. allocated(ElectrostaticMap)) then
123     allocate(ElectrostaticMap(nAtypes))
124     endif
125    
126     end if
127    
128     if (myATID .gt. size(ElectrostaticMap)) then
129     status = -1
130     return
131     endif
132    
133     ! set the values for ElectrostaticMap for this atom type:
134    
135     ElectrostaticMap(myATID)%c_ident = c_ident
136     ElectrostaticMap(myATID)%is_Charge = is_Charge
137     ElectrostaticMap(myATID)%is_Dipole = is_Dipole
138     ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole
139     ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole
140    
141     end subroutine newElectrostaticType
142    
143     subroutine setCharge(c_ident, charge, status)
144     integer, intent(in) :: c_ident
145     real(kind=dp), intent(in) :: charge
146     integer, intent(out) :: status
147     integer :: myATID
148    
149     status = 0
150     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
151    
152     if (.not.allocated(ElectrostaticMap)) then
153     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!")
154     status = -1
155     return
156     end if
157    
158     if (myATID .gt. size(ElectrostaticMap)) then
159     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!")
160     status = -1
161     return
162     endif
163    
164     if (.not.ElectrostaticMap(myATID)%is_Charge) then
165     call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!")
166     status = -1
167     return
168     endif
169    
170     ElectrostaticMap(myATID)%charge = charge
171     end subroutine setCharge
172    
173     subroutine setDipoleMoment(c_ident, dipole_moment, status)
174     integer, intent(in) :: c_ident
175     real(kind=dp), intent(in) :: dipole_moment
176     integer, intent(out) :: status
177     integer :: myATID
178    
179     status = 0
180     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
181    
182     if (.not.allocated(ElectrostaticMap)) then
183     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!")
184     status = -1
185     return
186     end if
187    
188     if (myATID .gt. size(ElectrostaticMap)) then
189     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!")
190     status = -1
191     return
192     endif
193    
194     if (.not.ElectrostaticMap(myATID)%is_Dipole) then
195     call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!")
196     status = -1
197     return
198     endif
199    
200     ElectrostaticMap(myATID)%dipole_moment = dipole_moment
201     end subroutine setDipoleMoment
202    
203     subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status)
204     integer, intent(in) :: c_ident
205     real(kind=dp), intent(in) :: split_dipole_distance
206     integer, intent(out) :: status
207     integer :: myATID
208    
209     status = 0
210     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
211    
212     if (.not.allocated(ElectrostaticMap)) then
213     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!")
214     status = -1
215     return
216     end if
217    
218     if (myATID .gt. size(ElectrostaticMap)) then
219     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!")
220     status = -1
221     return
222     endif
223    
224     if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then
225     call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!")
226     status = -1
227     return
228     endif
229    
230     ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance
231     end subroutine setSplitDipoleDistance
232    
233     subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status)
234     integer, intent(in) :: c_ident
235     real(kind=dp), intent(in), dimension(3) :: quadrupole_moments
236     integer, intent(out) :: status
237     integer :: myATID, i, j
238    
239     status = 0
240     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
241    
242     if (.not.allocated(ElectrostaticMap)) then
243     call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!")
244     status = -1
245     return
246     end if
247    
248     if (myATID .gt. size(ElectrostaticMap)) then
249     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!")
250     status = -1
251     return
252     endif
253    
254     if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then
255     call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!")
256     status = -1
257     return
258     endif
259    
260     do i = 1, 3
261     ElectrostaticMap(myATID)%quadrupole_moments(i) = &
262     quadrupole_moments(i)
263     enddo
264    
265     end subroutine setQuadrupoleMoments
266    
267    
268     function getCharge(atid) result (c)
269     integer, intent(in) :: atid
270     integer :: localError
271     real(kind=dp) :: c
272    
273     if (.not.allocated(ElectrostaticMap)) then
274     call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
275     return
276     end if
277    
278     if (.not.ElectrostaticMap(atid)%is_Charge) then
279     call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!")
280     return
281     endif
282    
283     c = ElectrostaticMap(atid)%charge
284     end function getCharge
285    
286     function getDipoleMoment(atid) result (dm)
287     integer, intent(in) :: atid
288     integer :: localError
289     real(kind=dp) :: dm
290    
291     if (.not.allocated(ElectrostaticMap)) then
292     call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
293     return
294     end if
295    
296     if (.not.ElectrostaticMap(atid)%is_Dipole) then
297     call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!")
298     return
299     endif
300    
301     dm = ElectrostaticMap(atid)%dipole_moment
302     end function getDipoleMoment
303    
304     subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, &
305     vpair, fpair, pot, eFrame, f, t, do_pot)
306    
307     logical, intent(in) :: do_pot
308    
309     integer, intent(in) :: atom1, atom2
310     integer :: localError
311    
312     real(kind=dp), intent(in) :: rij, r2, sw
313     real(kind=dp), intent(in), dimension(3) :: d
314     real(kind=dp), intent(inout) :: vpair
315     real(kind=dp), intent(inout), dimension(3) :: fpair
316    
317     real( kind = dp ) :: pot
318     real( kind = dp ), dimension(9,nLocal) :: eFrame
319     real( kind = dp ), dimension(3,nLocal) :: f
320     real( kind = dp ), dimension(3,nLocal) :: t
321    
322 gezelter 2123 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
323     real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
324     real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
325     real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
326 gezelter 2095
327     logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
328     logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
329     integer :: me1, me2, id1, id2
330     real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
331 gezelter 2123 real (kind=dp) :: qxx_i, qyy_i, qzz_i
332     real (kind=dp) :: qxx_j, qyy_j, qzz_j
333     real (kind=dp) :: cx_i, cy_i, cz_i
334     real (kind=dp) :: cx_j, cy_j, cz_j
335     real (kind=dp) :: cx2, cy2, cz2
336 gezelter 2095 real (kind=dp) :: ct_i, ct_j, ct_ij, a1
337 gezelter 2105 real (kind=dp) :: riji, ri, ri2, ri3, ri4
338 gezelter 2095 real (kind=dp) :: pref, vterm, epot, dudr
339 gezelter 2105 real (kind=dp) :: xhat, yhat, zhat
340 gezelter 2095 real (kind=dp) :: dudx, dudy, dudz
341     real (kind=dp) :: drdxj, drdyj, drdzj
342 gezelter 2105 real (kind=dp) :: scale, sc2, bigR
343 gezelter 2095
344     if (.not.allocated(ElectrostaticMap)) then
345     call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
346     return
347     end if
348    
349     #ifdef IS_MPI
350     me1 = atid_Row(atom1)
351     me2 = atid_Col(atom2)
352     #else
353     me1 = atid(atom1)
354     me2 = atid(atom2)
355     #endif
356    
357     !! some variables we'll need independent of electrostatic type:
358    
359     riji = 1.0d0 / rij
360    
361 gezelter 2105 xhat = d(1) * riji
362     yhat = d(2) * riji
363     zhat = d(3) * riji
364 gezelter 2095
365 gezelter 2105 drdxj = xhat
366     drdyj = yhat
367     drdzj = zhat
368 gezelter 2095
369     !! logicals
370    
371     i_is_Charge = ElectrostaticMap(me1)%is_Charge
372     i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
373     i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
374     i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
375    
376     j_is_Charge = ElectrostaticMap(me2)%is_Charge
377     j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
378     j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
379     j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
380    
381     if (i_is_Charge) then
382     q_i = ElectrostaticMap(me1)%charge
383     endif
384    
385     if (i_is_Dipole) then
386     mu_i = ElectrostaticMap(me1)%dipole_moment
387     #ifdef IS_MPI
388 gezelter 2123 uz_i(1) = eFrame_Row(3,atom1)
389     uz_i(2) = eFrame_Row(6,atom1)
390     uz_i(3) = eFrame_Row(9,atom1)
391 gezelter 2095 #else
392 gezelter 2123 uz_i(1) = eFrame(3,atom1)
393     uz_i(2) = eFrame(6,atom1)
394     uz_i(3) = eFrame(9,atom1)
395 gezelter 2095 #endif
396 gezelter 2123 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
397 gezelter 2095
398     if (i_is_SplitDipole) then
399     d_i = ElectrostaticMap(me1)%split_dipole_distance
400     endif
401    
402     endif
403    
404 gezelter 2123 if (i_is_Quadrupole) then
405     qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
406     qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
407     qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
408     #ifdef IS_MPI
409     ux_i(1) = eFrame_Row(1,atom1)
410     ux_i(2) = eFrame_Row(4,atom1)
411     ux_i(3) = eFrame_Row(7,atom1)
412     uy_i(1) = eFrame_Row(2,atom1)
413     uy_i(2) = eFrame_Row(5,atom1)
414     uy_i(3) = eFrame_Row(8,atom1)
415     uz_i(1) = eFrame_Row(3,atom1)
416     uz_i(2) = eFrame_Row(6,atom1)
417     uz_i(3) = eFrame_Row(9,atom1)
418     #else
419     ux_i(1) = eFrame(1,atom1)
420     ux_i(2) = eFrame(4,atom1)
421     ux_i(3) = eFrame(7,atom1)
422     uy_i(1) = eFrame(2,atom1)
423     uy_i(2) = eFrame(5,atom1)
424     uy_i(3) = eFrame(8,atom1)
425     uz_i(1) = eFrame(3,atom1)
426     uz_i(2) = eFrame(6,atom1)
427     uz_i(3) = eFrame(9,atom1)
428     #endif
429     cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
430     cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
431     cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
432     endif
433    
434    
435 gezelter 2095 if (j_is_Charge) then
436     q_j = ElectrostaticMap(me2)%charge
437     endif
438    
439     if (j_is_Dipole) then
440     mu_j = ElectrostaticMap(me2)%dipole_moment
441     #ifdef IS_MPI
442 gezelter 2123 uz_j(1) = eFrame_Col(3,atom2)
443     uz_j(2) = eFrame_Col(6,atom2)
444     uz_j(3) = eFrame_Col(9,atom2)
445 gezelter 2095 #else
446 gezelter 2123 uz_j(1) = eFrame(3,atom2)
447     uz_j(2) = eFrame(6,atom2)
448     uz_j(3) = eFrame(9,atom2)
449 gezelter 2095 #endif
450 gezelter 2123 ct_j = uz_j(1)*drdxj + uz_j(2)*drdyj + uz_j(3)*drdzj
451 gezelter 2095
452     if (j_is_SplitDipole) then
453     d_j = ElectrostaticMap(me2)%split_dipole_distance
454     endif
455     endif
456    
457 gezelter 2123 if (j_is_Quadrupole) then
458     qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
459     qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
460     qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
461     #ifdef IS_MPI
462     ux_j(1) = eFrame_Col(1,atom2)
463     ux_j(2) = eFrame_Col(4,atom2)
464     ux_j(3) = eFrame_Col(7,atom2)
465     uy_j(1) = eFrame_Col(2,atom2)
466     uy_j(2) = eFrame_Col(5,atom2)
467     uy_j(3) = eFrame_Col(8,atom2)
468     uz_j(1) = eFrame_Col(3,atom2)
469     uz_j(2) = eFrame_Col(6,atom2)
470     uz_j(3) = eFrame_Col(9,atom2)
471     #else
472     ux_j(1) = eFrame(1,atom2)
473     ux_j(2) = eFrame(4,atom2)
474     ux_j(3) = eFrame(7,atom2)
475     uy_j(1) = eFrame(2,atom2)
476     uy_j(2) = eFrame(5,atom2)
477     uy_j(3) = eFrame(8,atom2)
478     uz_j(1) = eFrame(3,atom2)
479     uz_j(2) = eFrame(6,atom2)
480     uz_j(3) = eFrame(9,atom2)
481     #endif
482     cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
483     cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
484     cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
485     endif
486    
487 gezelter 2095 epot = 0.0_dp
488     dudx = 0.0_dp
489     dudy = 0.0_dp
490     dudz = 0.0_dp
491    
492 gezelter 2123 dudux_i = 0.0_dp
493     duduy_i = 0.0_dp
494     duduz_i = 0.0_dp
495 gezelter 2095
496 gezelter 2123 dudux_j = 0.0_dp
497     duduy_j = 0.0_dp
498     duduz_j = 0.0_dp
499 gezelter 2095
500     if (i_is_Charge) then
501    
502     if (j_is_Charge) then
503    
504     vterm = pre11 * q_i * q_j * riji
505     vpair = vpair + vterm
506     epot = epot + sw*vterm
507    
508     dudr = - sw * vterm * riji
509    
510     dudx = dudx + dudr * drdxj
511     dudy = dudy + dudr * drdyj
512     dudz = dudz + dudr * drdzj
513    
514     endif
515    
516     if (j_is_Dipole) then
517    
518 gezelter 2105 if (j_is_SplitDipole) then
519     BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
520     ri = 1.0_dp / BigR
521     scale = rij * ri
522     else
523     ri = riji
524     scale = 1.0_dp
525     endif
526 gezelter 2095
527 gezelter 2105 ri2 = ri * ri
528     ri3 = ri2 * ri
529     sc2 = scale * scale
530    
531 gezelter 2095 pref = pre12 * q_i * mu_j
532 gezelter 2105 vterm = pref * ct_j * ri2 * scale
533 gezelter 2095 vpair = vpair + vterm
534     epot = epot + sw * vterm
535    
536 gezelter 2105 !! this has a + sign in the () because the rij vector is
537     !! r_j - r_i and the charge-dipole potential takes the origin
538     !! as the point dipole, which is atom j in this case.
539 gezelter 2095
540 gezelter 2123 dudx = dudx + pref * sw * ri3 * ( uz_j(1) + 3.0d0*ct_j*xhat*sc2)
541     dudy = dudy + pref * sw * ri3 * ( uz_j(2) + 3.0d0*ct_j*yhat*sc2)
542     dudz = dudz + pref * sw * ri3 * ( uz_j(3) + 3.0d0*ct_j*zhat*sc2)
543 gezelter 2105
544 gezelter 2123 duduz_j(1) = duduz_j(1) - pref * sw * ri2 * xhat * scale
545     duduz_j(2) = duduz_j(2) - pref * sw * ri2 * yhat * scale
546     duduz_j(3) = duduz_j(3) - pref * sw * ri2 * zhat * scale
547 gezelter 2095
548     endif
549 gezelter 2105
550 gezelter 2123 if (j_is_Quadrupole) then
551     ri2 = riji * riji
552     ri3 = ri2 * riji
553 gezelter 2124 ri4 = ri2 * ri2
554 gezelter 2123 cx2 = cx_j * cx_j
555     cy2 = cy_j * cy_j
556     cz2 = cz_j * cz_j
557    
558 gezelter 2127
559     pref = pre14 * q_i / 6.0_dp
560 gezelter 2123 vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
561     qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
562     qzz_j * (3.0_dp*cz2 - 1.0_dp))
563     vpair = vpair + vterm
564     epot = epot + sw * vterm
565    
566     dudx = dudx - 5.0_dp*sw*vterm*riji*xhat - pref * sw * ri4 * ( &
567     qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
568     qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
569     qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
570     dudy = dudy - 5.0_dp*sw*vterm*riji*yhat - pref * sw * ri4 * ( &
571     qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
572     qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
573     qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
574     dudz = dudz - 5.0_dp*sw*vterm*riji*zhat - pref * sw * ri4 * ( &
575     qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
576     qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
577     qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
578    
579     dudux_j(1) = dudux_j(1) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*xhat)
580     dudux_j(2) = dudux_j(2) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*yhat)
581     dudux_j(3) = dudux_j(3) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*zhat)
582    
583     duduy_j(1) = duduy_j(1) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*xhat)
584     duduy_j(2) = duduy_j(2) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*yhat)
585     duduy_j(3) = duduy_j(3) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*zhat)
586    
587     duduz_j(1) = duduz_j(1) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*xhat)
588     duduz_j(2) = duduz_j(2) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*yhat)
589     duduz_j(3) = duduz_j(3) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*zhat)
590     endif
591    
592 gezelter 2095 endif
593    
594     if (i_is_Dipole) then
595    
596     if (j_is_Charge) then
597    
598 gezelter 2105 if (i_is_SplitDipole) then
599     BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
600     ri = 1.0_dp / BigR
601     scale = rij * ri
602     else
603     ri = riji
604     scale = 1.0_dp
605     endif
606 gezelter 2095
607 gezelter 2105 ri2 = ri * ri
608     ri3 = ri2 * ri
609     sc2 = scale * scale
610    
611 gezelter 2095 pref = pre12 * q_j * mu_i
612 gezelter 2105 vterm = pref * ct_i * ri2 * scale
613 gezelter 2095 vpair = vpair + vterm
614     epot = epot + sw * vterm
615    
616 gezelter 2123 dudx = dudx + pref * sw * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
617     dudy = dudy + pref * sw * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
618     dudz = dudz + pref * sw * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
619 gezelter 2095
620 gezelter 2123 duduz_i(1) = duduz_i(1) + pref * sw * ri2 * xhat * scale
621     duduz_i(2) = duduz_i(2) + pref * sw * ri2 * yhat * scale
622     duduz_i(3) = duduz_i(3) + pref * sw * ri2 * zhat * scale
623 gezelter 2095 endif
624    
625     if (j_is_Dipole) then
626    
627 gezelter 2105 if (i_is_SplitDipole) then
628     if (j_is_SplitDipole) then
629     BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
630     else
631     BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
632     endif
633     ri = 1.0_dp / BigR
634     scale = rij * ri
635     else
636     if (j_is_SplitDipole) then
637     BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
638     ri = 1.0_dp / BigR
639     scale = rij * ri
640     else
641     ri = riji
642     scale = 1.0_dp
643     endif
644     endif
645    
646 gezelter 2123 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
647 gezelter 2105
648     ri2 = ri * ri
649     ri3 = ri2 * ri
650 gezelter 2095 ri4 = ri2 * ri2
651 gezelter 2105 sc2 = scale * scale
652 gezelter 2095
653     pref = pre22 * mu_i * mu_j
654 gezelter 2105 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
655 gezelter 2095 vpair = vpair + vterm
656     epot = epot + sw * vterm
657    
658 gezelter 2105 a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
659 gezelter 2095
660 gezelter 2123 dudx=dudx+pref*sw*3.0d0*ri4*scale*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
661     dudy=dudy+pref*sw*3.0d0*ri4*scale*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
662     dudz=dudz+pref*sw*3.0d0*ri4*scale*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
663 gezelter 2095
664 gezelter 2123 duduz_i(1) = duduz_i(1) + pref*sw*ri3*(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
665     duduz_i(2) = duduz_i(2) + pref*sw*ri3*(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
666     duduz_i(3) = duduz_i(3) + pref*sw*ri3*(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
667 gezelter 2095
668 gezelter 2123 duduz_j(1) = duduz_j(1) + pref*sw*ri3*(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
669     duduz_j(2) = duduz_j(2) + pref*sw*ri3*(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
670     duduz_j(3) = duduz_j(3) + pref*sw*ri3*(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
671 gezelter 2095 endif
672    
673     endif
674 gezelter 2123
675     if (i_is_Quadrupole) then
676     if (j_is_Charge) then
677    
678     ri2 = riji * riji
679     ri3 = ri2 * riji
680 gezelter 2124 ri4 = ri2 * ri2
681 gezelter 2123 cx2 = cx_i * cx_i
682     cy2 = cy_i * cy_i
683     cz2 = cz_i * cz_i
684    
685 gezelter 2124 pref = pre14 * q_j / 6.0_dp
686 gezelter 2123 vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
687     qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
688     qzz_i * (3.0_dp*cz2 - 1.0_dp))
689     vpair = vpair + vterm
690     epot = epot + sw * vterm
691    
692     dudx = dudx - 5.0_dp*sw*vterm*riji*xhat - pref * sw * ri4 * ( &
693     qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
694     qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
695     qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
696     dudy = dudy - 5.0_dp*sw*vterm*riji*yhat - pref * sw * ri4 * ( &
697     qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
698     qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
699     qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
700     dudz = dudz - 5.0_dp*sw*vterm*riji*zhat - pref * sw * ri4 * ( &
701     qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
702     qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
703     qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
704    
705     dudux_i(1) = dudux_i(1) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*xhat)
706     dudux_i(2) = dudux_i(2) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*yhat)
707     dudux_i(3) = dudux_i(3) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*zhat)
708    
709     duduy_i(1) = duduy_i(1) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*xhat)
710     duduy_i(2) = duduy_i(2) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*yhat)
711     duduy_i(3) = duduy_i(3) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*zhat)
712    
713     duduz_i(1) = duduz_i(1) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*xhat)
714     duduz_i(2) = duduz_i(2) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*yhat)
715     duduz_i(3) = duduz_i(3) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*zhat)
716     endif
717     endif
718    
719 gezelter 2095
720     if (do_pot) then
721     #ifdef IS_MPI
722     pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
723     pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
724     #else
725     pot = pot + epot
726     #endif
727     endif
728    
729     #ifdef IS_MPI
730     f_Row(1,atom1) = f_Row(1,atom1) + dudx
731     f_Row(2,atom1) = f_Row(2,atom1) + dudy
732     f_Row(3,atom1) = f_Row(3,atom1) + dudz
733    
734     f_Col(1,atom2) = f_Col(1,atom2) - dudx
735     f_Col(2,atom2) = f_Col(2,atom2) - dudy
736     f_Col(3,atom2) = f_Col(3,atom2) - dudz
737    
738     if (i_is_Dipole .or. i_is_Quadrupole) then
739 gezelter 2123 t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
740     t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
741     t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
742 gezelter 2095 endif
743 gezelter 2123 if (i_is_Quadrupole) then
744     t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
745     t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
746     t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
747 gezelter 2095
748 gezelter 2123 t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
749     t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
750     t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
751     endif
752    
753 gezelter 2095 if (j_is_Dipole .or. j_is_Quadrupole) then
754 gezelter 2123 t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
755     t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
756     t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
757 gezelter 2095 endif
758 gezelter 2123 if (j_is_Quadrupole) then
759     t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
760     t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
761     t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
762 gezelter 2095
763 gezelter 2123 t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
764     t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
765     t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
766     endif
767    
768 gezelter 2095 #else
769     f(1,atom1) = f(1,atom1) + dudx
770     f(2,atom1) = f(2,atom1) + dudy
771     f(3,atom1) = f(3,atom1) + dudz
772    
773     f(1,atom2) = f(1,atom2) - dudx
774     f(2,atom2) = f(2,atom2) - dudy
775     f(3,atom2) = f(3,atom2) - dudz
776    
777     if (i_is_Dipole .or. i_is_Quadrupole) then
778 gezelter 2123 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
779     t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
780     t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
781 gezelter 2095 endif
782 gezelter 2123 if (i_is_Quadrupole) then
783     t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
784     t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
785     t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
786    
787     t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
788     t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
789     t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
790     endif
791    
792 gezelter 2095 if (j_is_Dipole .or. j_is_Quadrupole) then
793 gezelter 2123 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
794     t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
795     t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
796 gezelter 2095 endif
797 gezelter 2123 if (j_is_Quadrupole) then
798     t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
799     t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
800     t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
801    
802     t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
803     t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
804     t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
805     endif
806    
807 gezelter 2095 #endif
808    
809     #ifdef IS_MPI
810     id1 = AtomRowToGlobal(atom1)
811     id2 = AtomColToGlobal(atom2)
812     #else
813     id1 = atom1
814     id2 = atom2
815     #endif
816    
817     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
818    
819     fpair(1) = fpair(1) + dudx
820     fpair(2) = fpair(2) + dudy
821     fpair(3) = fpair(3) + dudz
822    
823     endif
824    
825     return
826     end subroutine doElectrostaticPair
827    
828     end module electrostatic_module