ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2189
Committed: Wed Apr 13 20:36:45 2005 UTC (19 years, 4 months ago) by chuckv
File size: 29189 byte(s)
Log Message:
Added destroy methods for Fortran modules.

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