ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2279
Committed: Tue Aug 30 18:23:50 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 30885 byte(s)
Log Message:
made some changes for implementing the wolf potential

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