ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 1467
Committed: Sat Jul 17 15:33:03 2010 UTC (14 years, 11 months ago) by gezelter
File size: 48954 byte(s)
Log Message:
well, it compiles, but still segfaults

File Contents

# User Rev Content
1 gezelter 411 !!
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 gezelter 1390 !! 1. Redistributions of source code must retain the above copyright
10 gezelter 411 !! notice, this list of conditions and the following disclaimer.
11     !!
12 gezelter 1390 !! 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 411 !! notice, this list of conditions and the following disclaimer in the
14     !! documentation and/or other materials provided with the
15     !! distribution.
16     !!
17     !! This software is provided "AS IS," without a warranty of any
18     !! kind. All express or implied conditions, representations and
19     !! warranties, including any implied warranty of merchantability,
20     !! fitness for a particular purpose or non-infringement, are hereby
21     !! excluded. The University of Notre Dame and its licensors shall not
22     !! be liable for any damages suffered by licensee as a result of
23     !! using, modifying or distributing the software or its
24     !! derivatives. In no event will the University of Notre Dame or its
25     !! licensors be liable for any lost revenue, profit or data, or for
26     !! direct, indirect, special, consequential, incidental or punitive
27     !! damages, however caused and regardless of the theory of liability,
28     !! arising out of the use of or inability to use software, even if the
29     !! University of Notre Dame has been advised of the possibility of
30     !! such damages.
31     !!
32 gezelter 1390 !! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     !! research, please cite the appropriate papers when you publish your
34     !! work. Good starting points are:
35     !!
36     !! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     !! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     !! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     !! [4] Vardeman & Gezelter, in progress (2009).
40     !!
41 gezelter 411
42     module electrostatic_module
43 gezelter 507
44 gezelter 411 use force_globals
45     use definitions
46     use atype_module
47     use vector_class
48     use simulation
49     use status
50 chrisfen 937 use interpolation
51 gezelter 411 implicit none
52    
53     PRIVATE
54    
55 chuckv 656
56 gezelter 602 #define __FORTRAN90
57 chuckv 656 #include "UseTheForce/DarkSide/fInteractionMap.h"
58 gezelter 602 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59 chrisfen 716 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
60 gezelter 602
61 chuckv 656
62 gezelter 434 !! these prefactors convert the multipole interactions into kcal / mol
63     !! all were computed assuming distances are measured in angstroms
64     !! Charge-Charge, assuming charges are measured in electrons
65 chrisfen 959 real(kind=dp), parameter :: pre11 = 332.0637778_dp
66 gezelter 434 !! Charge-Dipole, assuming charges are measured in electrons, and
67     !! dipoles are measured in debyes
68 chrisfen 959 real(kind=dp), parameter :: pre12 = 69.13373_dp
69 gezelter 434 !! Dipole-Dipole, assuming dipoles are measured in debyes
70 chrisfen 959 real(kind=dp), parameter :: pre22 = 14.39325_dp
71 gezelter 434 !! Charge-Quadrupole, assuming charges are measured in electrons, and
72     !! quadrupoles are measured in 10^-26 esu cm^2
73     !! This unit is also known affectionately as an esu centi-barn.
74 chrisfen 959 real(kind=dp), parameter :: pre14 = 69.13373_dp
75 gezelter 411
76 chrisfen 959 real(kind=dp), parameter :: zero = 0.0_dp
77 chrisfen 941
78 chrisfen 998 !! conversions for the simulation box dipole moment
79     real(kind=dp), parameter :: chargeToC = 1.60217733e-19_dp
80     real(kind=dp), parameter :: angstromToM = 1.0e-10_dp
81     real(kind=dp), parameter :: debyeToCm = 3.33564095198e-30_dp
82    
83 chrisfen 941 !! number of points for electrostatic splines
84     integer, parameter :: np = 100
85 gezelter 939
86 chrisfen 712 !! variables to handle different summation methods for long-range
87     !! electrostatics:
88 gezelter 602 integer, save :: summationMethod = NONE
89 chrisfen 710 integer, save :: screeningMethod = UNDAMPED
90 chrisfen 603 logical, save :: summationMethodChecked = .false.
91 gezelter 602 real(kind=DP), save :: defaultCutoff = 0.0_DP
92 chrisfen 682 real(kind=DP), save :: defaultCutoff2 = 0.0_DP
93 gezelter 602 logical, save :: haveDefaultCutoff = .false.
94     real(kind=DP), save :: dampingAlpha = 0.0_DP
95 chrisfen 987 real(kind=DP), save :: alpha2 = 0.0_DP
96     real(kind=DP), save :: alpha4 = 0.0_DP
97     real(kind=DP), save :: alpha6 = 0.0_DP
98     real(kind=DP), save :: alpha8 = 0.0_DP
99 gezelter 602 logical, save :: haveDampingAlpha = .false.
100 chrisfen 682 real(kind=DP), save :: dielectric = 1.0_DP
101 gezelter 602 logical, save :: haveDielectric = .false.
102     real(kind=DP), save :: constEXP = 0.0_DP
103 chrisfen 682 real(kind=dp), save :: rcuti = 0.0_DP
104     real(kind=dp), save :: rcuti2 = 0.0_DP
105     real(kind=dp), save :: rcuti3 = 0.0_DP
106     real(kind=dp), save :: rcuti4 = 0.0_DP
107     real(kind=dp), save :: alphaPi = 0.0_DP
108     real(kind=dp), save :: invRootPi = 0.0_DP
109     real(kind=dp), save :: rrf = 1.0_DP
110     real(kind=dp), save :: rt = 1.0_DP
111     real(kind=dp), save :: rrfsq = 1.0_DP
112     real(kind=dp), save :: preRF = 0.0_DP
113 chrisfen 695 real(kind=dp), save :: preRF2 = 0.0_DP
114 chrisfen 987 real(kind=dp), save :: erfcVal = 1.0_DP
115     real(kind=dp), save :: derfcVal = 0.0_DP
116     type(cubicSpline), save :: erfcSpline
117 chrisfen 941 logical, save :: haveElectroSpline = .false.
118 chrisfen 987 real(kind=dp), save :: c1 = 1.0_DP
119     real(kind=dp), save :: c2 = 1.0_DP
120     real(kind=dp), save :: c3 = 0.0_DP
121     real(kind=dp), save :: c4 = 0.0_DP
122     real(kind=dp), save :: c5 = 0.0_DP
123     real(kind=dp), save :: c6 = 0.0_DP
124     real(kind=dp), save :: c1c = 1.0_DP
125     real(kind=dp), save :: c2c = 1.0_DP
126     real(kind=dp), save :: c3c = 0.0_DP
127     real(kind=dp), save :: c4c = 0.0_DP
128     real(kind=dp), save :: c5c = 0.0_DP
129     real(kind=dp), save :: c6c = 0.0_DP
130 chrisfen 959 real(kind=dp), save :: one_third = 1.0_DP / 3.0_DP
131 chrisfen 716
132 gezelter 809 #if defined(__IFC) || defined(__PGI)
133 chuckv 632 ! error function for ifc version > 7.
134 gezelter 960 real(kind=dp), external :: erfc
135 chuckv 632 #endif
136    
137 chrisfen 853 public :: setElectrostaticSummationMethod
138 chrisfen 710 public :: setScreeningMethod
139 gezelter 602 public :: setElectrostaticCutoffRadius
140 chrisfen 710 public :: setDampingAlpha
141 gezelter 602 public :: setReactionFieldDielectric
142 chrisfen 941 public :: buildElectroSpline
143 gezelter 411 public :: newElectrostaticType
144     public :: setCharge
145     public :: setDipoleMoment
146     public :: setSplitDipoleDistance
147     public :: setQuadrupoleMoments
148     public :: doElectrostaticPair
149     public :: getCharge
150     public :: getDipoleMoment
151 chuckv 492 public :: destroyElectrostaticTypes
152 chrisfen 703 public :: self_self
153 chrisfen 700 public :: rf_self_excludes
154 chrisfen 998 public :: accumulate_box_dipole
155 gezelter 411
156     type :: Electrostatic
157     integer :: c_ident
158     logical :: is_Charge = .false.
159     logical :: is_Dipole = .false.
160     logical :: is_SplitDipole = .false.
161     logical :: is_Quadrupole = .false.
162 chrisfen 532 logical :: is_Tap = .false.
163 gezelter 411 real(kind=DP) :: charge = 0.0_DP
164     real(kind=DP) :: dipole_moment = 0.0_DP
165     real(kind=DP) :: split_dipole_distance = 0.0_DP
166     real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP
167     end type Electrostatic
168    
169     type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
170    
171 gezelter 938 logical, save :: hasElectrostaticMap
172    
173 gezelter 411 contains
174    
175 chrisfen 853 subroutine setElectrostaticSummationMethod(the_ESM)
176 gezelter 602 integer, intent(in) :: the_ESM
177    
178     if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then
179 chrisfen 853 call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method")
180 gezelter 602 endif
181    
182 chrisfen 610 summationMethod = the_ESM
183 chrisfen 626
184 chrisfen 853 end subroutine setElectrostaticSummationMethod
185 gezelter 602
186 chrisfen 710 subroutine setScreeningMethod(the_SM)
187     integer, intent(in) :: the_SM
188     screeningMethod = the_SM
189     end subroutine setScreeningMethod
190    
191 chrisfen 682 subroutine setElectrostaticCutoffRadius(thisRcut, thisRsw)
192 gezelter 602 real(kind=dp), intent(in) :: thisRcut
193 chrisfen 682 real(kind=dp), intent(in) :: thisRsw
194 gezelter 602 defaultCutoff = thisRcut
195 chrisfen 849 defaultCutoff2 = defaultCutoff*defaultCutoff
196 chrisfen 682 rrf = defaultCutoff
197     rt = thisRsw
198 gezelter 602 haveDefaultCutoff = .true.
199     end subroutine setElectrostaticCutoffRadius
200    
201 chrisfen 710 subroutine setDampingAlpha(thisAlpha)
202 gezelter 602 real(kind=dp), intent(in) :: thisAlpha
203     dampingAlpha = thisAlpha
204 chrisfen 716 alpha2 = dampingAlpha*dampingAlpha
205 chrisfen 987 alpha4 = alpha2*alpha2
206     alpha6 = alpha4*alpha2
207     alpha8 = alpha4*alpha4
208 gezelter 602 haveDampingAlpha = .true.
209 chrisfen 710 end subroutine setDampingAlpha
210 gezelter 602
211     subroutine setReactionFieldDielectric(thisDielectric)
212     real(kind=dp), intent(in) :: thisDielectric
213     dielectric = thisDielectric
214     haveDielectric = .true.
215     end subroutine setReactionFieldDielectric
216    
217 chrisfen 941 subroutine buildElectroSpline()
218     real( kind = dp ), dimension(np) :: xvals, yvals
219     real( kind = dp ) :: dx, rmin, rval
220     integer :: i
221 chrisfen 937
222 chrisfen 959 rmin = 0.0_dp
223 chrisfen 941
224     dx = (defaultCutoff-rmin) / dble(np-1)
225    
226     do i = 1, np
227     rval = rmin + dble(i-1)*dx
228     xvals(i) = rval
229 gezelter 960 yvals(i) = erfc(dampingAlpha*rval)
230 chrisfen 941 enddo
231    
232 chrisfen 987 call newSpline(erfcSpline, xvals, yvals, .true.)
233 chrisfen 941
234     haveElectroSpline = .true.
235     end subroutine buildElectroSpline
236    
237 gezelter 411 subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
238 chrisfen 532 is_SplitDipole, is_Quadrupole, is_Tap, status)
239 gezelter 507
240 gezelter 411 integer, intent(in) :: c_ident
241     logical, intent(in) :: is_Charge
242     logical, intent(in) :: is_Dipole
243     logical, intent(in) :: is_SplitDipole
244     logical, intent(in) :: is_Quadrupole
245 chrisfen 532 logical, intent(in) :: is_Tap
246 gezelter 411 integer, intent(out) :: status
247     integer :: nAtypes, myATID, i, j
248    
249     status = 0
250     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
251 gezelter 507
252 gezelter 411 !! Be simple-minded and assume that we need an ElectrostaticMap that
253     !! is the same size as the total number of atom types
254    
255     if (.not.allocated(ElectrostaticMap)) then
256 gezelter 507
257 gezelter 411 nAtypes = getSize(atypes)
258 gezelter 507
259 gezelter 411 if (nAtypes == 0) then
260     status = -1
261     return
262     end if
263 gezelter 507
264 gezelter 938 allocate(ElectrostaticMap(nAtypes))
265 gezelter 507
266 gezelter 411 end if
267    
268     if (myATID .gt. size(ElectrostaticMap)) then
269     status = -1
270     return
271     endif
272 gezelter 507
273 gezelter 411 ! set the values for ElectrostaticMap for this atom type:
274    
275     ElectrostaticMap(myATID)%c_ident = c_ident
276     ElectrostaticMap(myATID)%is_Charge = is_Charge
277     ElectrostaticMap(myATID)%is_Dipole = is_Dipole
278     ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole
279     ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole
280 chrisfen 532 ElectrostaticMap(myATID)%is_Tap = is_Tap
281 gezelter 507
282 gezelter 938 hasElectrostaticMap = .true.
283    
284 gezelter 411 end subroutine newElectrostaticType
285    
286     subroutine setCharge(c_ident, charge, status)
287     integer, intent(in) :: c_ident
288     real(kind=dp), intent(in) :: charge
289     integer, intent(out) :: status
290     integer :: myATID
291    
292     status = 0
293     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
294    
295 gezelter 938 if (.not.hasElectrostaticMap) then
296 gezelter 411 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!")
297     status = -1
298     return
299     end if
300    
301     if (myATID .gt. size(ElectrostaticMap)) then
302     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!")
303     status = -1
304     return
305     endif
306    
307     if (.not.ElectrostaticMap(myATID)%is_Charge) then
308     call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!")
309     status = -1
310     return
311 gezelter 507 endif
312 gezelter 411
313     ElectrostaticMap(myATID)%charge = charge
314     end subroutine setCharge
315    
316     subroutine setDipoleMoment(c_ident, dipole_moment, status)
317     integer, intent(in) :: c_ident
318     real(kind=dp), intent(in) :: dipole_moment
319     integer, intent(out) :: status
320     integer :: myATID
321    
322     status = 0
323     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
324    
325 gezelter 938 if (.not.hasElectrostaticMap) then
326 gezelter 411 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!")
327     status = -1
328     return
329     end if
330    
331     if (myATID .gt. size(ElectrostaticMap)) then
332     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!")
333     status = -1
334     return
335     endif
336    
337     if (.not.ElectrostaticMap(myATID)%is_Dipole) then
338     call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!")
339     status = -1
340     return
341     endif
342    
343     ElectrostaticMap(myATID)%dipole_moment = dipole_moment
344     end subroutine setDipoleMoment
345    
346     subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status)
347     integer, intent(in) :: c_ident
348     real(kind=dp), intent(in) :: split_dipole_distance
349     integer, intent(out) :: status
350     integer :: myATID
351    
352     status = 0
353     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
354    
355 gezelter 938 if (.not.hasElectrostaticMap) then
356 gezelter 411 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!")
357     status = -1
358     return
359     end if
360    
361     if (myATID .gt. size(ElectrostaticMap)) then
362     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!")
363     status = -1
364     return
365     endif
366    
367     if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then
368     call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!")
369     status = -1
370     return
371     endif
372    
373     ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance
374     end subroutine setSplitDipoleDistance
375    
376     subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status)
377     integer, intent(in) :: c_ident
378     real(kind=dp), intent(in), dimension(3) :: quadrupole_moments
379     integer, intent(out) :: status
380     integer :: myATID, i, j
381    
382     status = 0
383     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
384    
385 gezelter 938 if (.not.hasElectrostaticMap) then
386 gezelter 411 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!")
387     status = -1
388     return
389     end if
390    
391     if (myATID .gt. size(ElectrostaticMap)) then
392     call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!")
393     status = -1
394     return
395     endif
396    
397     if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then
398     call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!")
399     status = -1
400     return
401     endif
402 gezelter 507
403 gezelter 411 do i = 1, 3
404 gezelter 507 ElectrostaticMap(myATID)%quadrupole_moments(i) = &
405     quadrupole_moments(i)
406     enddo
407 gezelter 411
408     end subroutine setQuadrupoleMoments
409    
410 gezelter 507
411 gezelter 411 function getCharge(atid) result (c)
412     integer, intent(in) :: atid
413     integer :: localError
414     real(kind=dp) :: c
415 gezelter 507
416 gezelter 938 if (.not.hasElectrostaticMap) then
417 gezelter 411 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
418     return
419     end if
420 gezelter 507
421 gezelter 411 if (.not.ElectrostaticMap(atid)%is_Charge) then
422     call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!")
423     return
424     endif
425 gezelter 507
426 gezelter 411 c = ElectrostaticMap(atid)%charge
427     end function getCharge
428    
429     function getDipoleMoment(atid) result (dm)
430     integer, intent(in) :: atid
431     integer :: localError
432     real(kind=dp) :: dm
433 gezelter 507
434 gezelter 938 if (.not.hasElectrostaticMap) then
435 gezelter 411 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
436     return
437     end if
438 gezelter 507
439 gezelter 411 if (.not.ElectrostaticMap(atid)%is_Dipole) then
440     call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!")
441     return
442     endif
443 gezelter 507
444 gezelter 411 dm = ElectrostaticMap(atid)%dipole_moment
445     end function getDipoleMoment
446    
447 gezelter 602 subroutine checkSummationMethod()
448    
449 chrisfen 607 if (.not.haveDefaultCutoff) then
450     call handleError("checkSummationMethod", "no Default Cutoff set!")
451     endif
452    
453 chrisfen 959 rcuti = 1.0_dp / defaultCutoff
454 chrisfen 607 rcuti2 = rcuti*rcuti
455     rcuti3 = rcuti2*rcuti
456     rcuti4 = rcuti2*rcuti2
457    
458 chrisfen 710 if (screeningMethod .eq. DAMPED) then
459 chrisfen 703 if (.not.haveDampingAlpha) then
460     call handleError("checkSummationMethod", "no Damping Alpha set!")
461     endif
462    
463     if (.not.haveDefaultCutoff) then
464     call handleError("checkSummationMethod", "no Default Cutoff set!")
465     endif
466 chrisfen 603
467 chrisfen 849 constEXP = exp(-alpha2*defaultCutoff2)
468 chrisfen 959 invRootPi = 0.56418958354775628695_dp
469     alphaPi = 2.0_dp*dampingAlpha*invRootPi
470 chrisfen 987
471     c1c = erfc(dampingAlpha*defaultCutoff) * rcuti
472     c2c = alphaPi*constEXP*rcuti + c1c*rcuti
473     c3c = 2.0_dp*alphaPi*alpha2 + 3.0_dp*c2c*rcuti
474     c4c = 4.0_dp*alphaPi*alpha4 + 5.0_dp*c3c*rcuti2
475     c5c = 8.0_dp*alphaPi*alpha6 + 7.0_dp*c4c*rcuti2
476     c6c = 16.0_dp*alphaPi*alpha8 + 9.0_dp*c5c*rcuti2
477     else
478     c1c = rcuti
479     c2c = c1c*rcuti
480     c3c = 3.0_dp*c2c*rcuti
481     c4c = 5.0_dp*c3c*rcuti2
482     c5c = 7.0_dp*c4c*rcuti2
483     c6c = 9.0_dp*c5c*rcuti2
484 gezelter 602 endif
485    
486 chrisfen 603 if (summationMethod .eq. REACTION_FIELD) then
487 chrisfen 703 if (haveDielectric) then
488     defaultCutoff2 = defaultCutoff*defaultCutoff
489 chrisfen 959 preRF = (dielectric-1.0_dp) / &
490     ((2.0_dp*dielectric+1.0_dp)*defaultCutoff2*defaultCutoff)
491     preRF2 = 2.0_dp*preRF
492 chrisfen 703 else
493     call handleError("checkSummationMethod", "Dielectric not set")
494 chrisfen 603 endif
495 chrisfen 703
496 chrisfen 603 endif
497    
498 chrisfen 941 if (.not.haveElectroSpline) then
499     call buildElectroSpline()
500     end if
501    
502 chrisfen 603 summationMethodChecked = .true.
503 gezelter 602 end subroutine checkSummationMethod
504    
505 chrisfen 712
506 gezelter 1464 subroutine doElectrostaticPair(me1, me2, d, rij, r2, rcut, sw, &
507 gezelter 1467 electroMult, vpair, pot, eF1, eF2, f1, t1, t2)
508 gezelter 507
509 gezelter 1464 integer, intent(in) :: me1, me2
510 gezelter 411 integer :: localError
511    
512 gezelter 1286 real(kind=dp), intent(in) :: rij, r2, sw, rcut, electroMult
513 gezelter 411 real(kind=dp), intent(in), dimension(3) :: d
514     real(kind=dp), intent(inout) :: vpair
515    
516 chrisfen 626 real( kind = dp ) :: pot
517 gezelter 1386 real( kind = dp ), dimension(9) :: eF1, eF2 ! eFrame = electroFrame
518     real( kind = dp ), dimension(3) :: f1
519 chrisfen 710 real( kind = dp ), dimension(3,nLocal) :: felec
520 gezelter 1386 real( kind = dp ), dimension(3) :: t1, t2
521 gezelter 507
522 gezelter 439 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
523     real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
524     real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
525     real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
526 gezelter 411
527     logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
528     logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
529 chrisfen 532 logical :: i_is_Tap, j_is_Tap
530 gezelter 1386 integer :: id1, id2
531 gezelter 411 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
532 gezelter 439 real (kind=dp) :: qxx_i, qyy_i, qzz_i
533     real (kind=dp) :: qxx_j, qyy_j, qzz_j
534     real (kind=dp) :: cx_i, cy_i, cz_i
535     real (kind=dp) :: cx_j, cy_j, cz_j
536     real (kind=dp) :: cx2, cy2, cz2
537 chrisfen 719 real (kind=dp) :: ct_i, ct_j, ct_ij, a0, a1
538 gezelter 421 real (kind=dp) :: riji, ri, ri2, ri3, ri4
539 chrisfen 597 real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2
540 gezelter 421 real (kind=dp) :: xhat, yhat, zhat
541 gezelter 411 real (kind=dp) :: dudx, dudy, dudz
542 chrisfen 626 real (kind=dp) :: scale, sc2, bigR
543 chrisfen 716 real (kind=dp) :: varEXP
544 chrisfen 719 real (kind=dp) :: pot_term
545 chrisfen 695 real (kind=dp) :: preVal, rfVal
546 chrisfen 987 real (kind=dp) :: c2ri, c3ri, c4rij
547 chrisfen 959 real (kind=dp) :: cti3, ctj3, ctidotj
548 chrisfen 987 real (kind=dp) :: preSw, preSwSc
549 chrisfen 959 real (kind=dp) :: xhatdot2, yhatdot2, zhatdot2
550 chrisfen 987 real (kind=dp) :: xhatc4, yhatc4, zhatc4
551 gezelter 411
552 gezelter 602 if (.not.summationMethodChecked) then
553     call checkSummationMethod()
554     endif
555    
556 gezelter 411 !! some variables we'll need independent of electrostatic type:
557    
558 chrisfen 959 riji = 1.0_dp / rij
559 chrisfen 644
560 gezelter 421 xhat = d(1) * riji
561     yhat = d(2) * riji
562     zhat = d(3) * riji
563 gezelter 411
564     !! logicals
565     i_is_Charge = ElectrostaticMap(me1)%is_Charge
566     i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
567     i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
568     i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
569 chrisfen 532 i_is_Tap = ElectrostaticMap(me1)%is_Tap
570 gezelter 411
571     j_is_Charge = ElectrostaticMap(me2)%is_Charge
572     j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
573     j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
574     j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
575 chrisfen 532 j_is_Tap = ElectrostaticMap(me2)%is_Tap
576 gezelter 411
577     if (i_is_Charge) then
578     q_i = ElectrostaticMap(me1)%charge
579     endif
580 gezelter 507
581 gezelter 411 if (i_is_Dipole) then
582     mu_i = ElectrostaticMap(me1)%dipole_moment
583 gezelter 1386
584     uz_i(1) = eF1(3)
585     uz_i(2) = eF1(6)
586     uz_i(3) = eF1(9)
587    
588 gezelter 439 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
589 gezelter 411
590     if (i_is_SplitDipole) then
591     d_i = ElectrostaticMap(me1)%split_dipole_distance
592     endif
593 gezelter 939 duduz_i = zero
594 gezelter 411 endif
595    
596 gezelter 439 if (i_is_Quadrupole) then
597     qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
598     qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
599     qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
600 gezelter 1386
601     ux_i(1) = eF1(1)
602     ux_i(2) = eF1(4)
603     ux_i(3) = eF1(7)
604     uy_i(1) = eF1(2)
605     uy_i(2) = eF1(5)
606     uy_i(3) = eF1(8)
607     uz_i(1) = eF1(3)
608     uz_i(2) = eF1(6)
609     uz_i(3) = eF1(9)
610    
611 gezelter 439 cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
612     cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
613     cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
614 gezelter 939 dudux_i = zero
615     duduy_i = zero
616     duduz_i = zero
617 gezelter 439 endif
618    
619 gezelter 411 if (j_is_Charge) then
620     q_j = ElectrostaticMap(me2)%charge
621     endif
622 gezelter 507
623 gezelter 411 if (j_is_Dipole) then
624     mu_j = ElectrostaticMap(me2)%dipole_moment
625 gezelter 1386
626     uz_j(1) = eF2(3)
627     uz_j(2) = eF2(6)
628     uz_j(3) = eF2(9)
629    
630 chrisfen 465 ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
631 gezelter 411
632     if (j_is_SplitDipole) then
633     d_j = ElectrostaticMap(me2)%split_dipole_distance
634     endif
635 gezelter 939 duduz_j = zero
636 gezelter 411 endif
637    
638 gezelter 439 if (j_is_Quadrupole) then
639     qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
640     qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
641     qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
642 gezelter 1386
643     ux_j(1) = eF2(1)
644     ux_j(2) = eF2(4)
645     ux_j(3) = eF2(7)
646     uy_j(1) = eF2(2)
647     uy_j(2) = eF2(5)
648     uy_j(3) = eF2(8)
649     uz_j(1) = eF2(3)
650     uz_j(2) = eF2(6)
651     uz_j(3) = eF2(9)
652    
653 gezelter 439 cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
654     cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
655     cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
656 gezelter 939 dudux_j = zero
657     duduy_j = zero
658     duduz_j = zero
659 gezelter 439 endif
660 chrisfen 554
661 gezelter 939 epot = zero
662     dudx = zero
663     dudy = zero
664     dudz = zero
665 gezelter 411
666     if (i_is_Charge) then
667    
668     if (j_is_Charge) then
669 chrisfen 739 if (screeningMethod .eq. DAMPED) then
670 chrisfen 959 ! assemble the damping variables
671 chrisfen 987 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
672     c1 = erfcVal*riji
673     c2 = (-derfcVal + c1)*riji
674     else
675     c1 = riji
676     c2 = c1*riji
677 chrisfen 739 endif
678 gezelter 507
679 gezelter 1286 preVal = electroMult * pre11 * q_i * q_j
680 chrisfen 739
681 chrisfen 710 if (summationMethod .eq. SHIFTED_POTENTIAL) then
682 chrisfen 987 vterm = preVal * (c1 - c1c)
683 chrisfen 597
684 chrisfen 987 dudr = -sw * preVal * c2
685 chrisfen 739
686 chrisfen 710 elseif (summationMethod .eq. SHIFTED_FORCE) then
687 chrisfen 987 vterm = preVal * ( c1 - c1c + c2c*(rij - defaultCutoff) )
688 chrisfen 716
689 chrisfen 987 dudr = sw * preVal * (c2c - c2)
690 chrisfen 739
691 chrisfen 695 elseif (summationMethod .eq. REACTION_FIELD) then
692 gezelter 1286 rfVal = electroMult * preRF*rij*rij
693 chrisfen 695 vterm = preVal * ( riji + rfVal )
694 chrisfen 700
695 chrisfen 959 dudr = sw * preVal * ( 2.0_dp*rfVal - riji )*riji
696 chrisfen 739
697 chrisfen 597 else
698 chrisfen 987 vterm = preVal * riji*erfcVal
699 chrisfen 597
700 chrisfen 987 dudr = - sw * preVal * c2
701 chrisfen 739
702 chrisfen 597 endif
703    
704 chrisfen 739 vpair = vpair + vterm
705     epot = epot + sw*vterm
706    
707     dudx = dudx + dudr * xhat
708     dudy = dudy + dudr * yhat
709     dudz = dudz + dudr * zhat
710    
711 gezelter 411 endif
712    
713     if (j_is_Dipole) then
714 chrisfen 987 ! pref is used by all the possible methods
715 gezelter 1286 pref = electroMult * pre12 * q_i * mu_j
716 chrisfen 987 preSw = sw*pref
717 gezelter 411
718 chrisfen 710 if (summationMethod .eq. REACTION_FIELD) then
719 chrisfen 700 ri2 = riji * riji
720     ri3 = ri2 * riji
721 chrisfen 696
722     vterm = - pref * ct_j * ( ri2 - preRF2*rij )
723     vpair = vpair + vterm
724     epot = epot + sw*vterm
725    
726 chrisfen 987 dudx = dudx - preSw*( ri3*(uz_j(1) - 3.0_dp*ct_j*xhat) - &
727     preRF2*uz_j(1) )
728     dudy = dudy - preSw*( ri3*(uz_j(2) - 3.0_dp*ct_j*yhat) - &
729     preRF2*uz_j(2) )
730     dudz = dudz - preSw*( ri3*(uz_j(3) - 3.0_dp*ct_j*zhat) - &
731     preRF2*uz_j(3) )
732     duduz_j(1) = duduz_j(1) - preSw * xhat * ( ri2 - preRF2*rij )
733     duduz_j(2) = duduz_j(2) - preSw * yhat * ( ri2 - preRF2*rij )
734     duduz_j(3) = duduz_j(3) - preSw * zhat * ( ri2 - preRF2*rij )
735 chrisfen 696
736 chrisfen 597 else
737 chrisfen 987 ! determine the inverse r used if we have split dipoles
738 chrisfen 597 if (j_is_SplitDipole) then
739 chrisfen 959 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
740     ri = 1.0_dp / BigR
741 chrisfen 597 scale = rij * ri
742     else
743     ri = riji
744 chrisfen 959 scale = 1.0_dp
745 chrisfen 597 endif
746 chrisfen 987
747 chrisfen 597 sc2 = scale * scale
748 chrisfen 626
749 chrisfen 987 if (screeningMethod .eq. DAMPED) then
750     ! assemble the damping variables
751     call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
752     c1 = erfcVal*ri
753     c2 = (-derfcVal + c1)*ri
754     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
755     else
756     c1 = ri
757     c2 = c1*ri
758     c3 = 3.0_dp*c2*ri
759     endif
760    
761     c2ri = c2*ri
762    
763     ! calculate the potential
764     pot_term = scale * c2
765 chrisfen 959 vterm = -pref * ct_j * pot_term
766 chrisfen 626 vpair = vpair + vterm
767     epot = epot + sw*vterm
768 chrisfen 597
769 chrisfen 987 ! calculate derivatives for forces and torques
770     dudx = dudx - preSw*( uz_j(1)*c2ri - ct_j*xhat*sc2*c3 )
771     dudy = dudy - preSw*( uz_j(2)*c2ri - ct_j*yhat*sc2*c3 )
772     dudz = dudz - preSw*( uz_j(3)*c2ri - ct_j*zhat*sc2*c3 )
773 chrisfen 849
774 chrisfen 987 duduz_j(1) = duduz_j(1) - preSw * pot_term * xhat
775     duduz_j(2) = duduz_j(2) - preSw * pot_term * yhat
776     duduz_j(3) = duduz_j(3) - preSw * pot_term * zhat
777 gezelter 411
778 chrisfen 597 endif
779 gezelter 411 endif
780 gezelter 421
781 gezelter 439 if (j_is_Quadrupole) then
782 chrisfen 987 ! first precalculate some necessary variables
783     cx2 = cx_j * cx_j
784     cy2 = cy_j * cy_j
785     cz2 = cz_j * cz_j
786 gezelter 1286 pref = electroMult * pre14 * q_i * one_third
787 chrisfen 987
788 chrisfen 849 if (screeningMethod .eq. DAMPED) then
789 chrisfen 959 ! assemble the damping variables
790 chrisfen 987 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
791     c1 = erfcVal*riji
792     c2 = (-derfcVal + c1)*riji
793     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*riji
794     c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*riji*riji
795     else
796     c1 = riji
797     c2 = c1*riji
798     c3 = 3.0_dp*c2*riji
799     c4 = 5.0_dp*c3*riji*riji
800 chrisfen 849 endif
801    
802 chrisfen 987 ! precompute variables for convenience
803     preSw = sw*pref
804     c2ri = c2*riji
805     c3ri = c3*riji
806     c4rij = c4*rij
807     xhatdot2 = 2.0_dp*xhat*c3
808     yhatdot2 = 2.0_dp*yhat*c3
809     zhatdot2 = 2.0_dp*zhat*c3
810     xhatc4 = xhat*c4rij
811     yhatc4 = yhat*c4rij
812     zhatc4 = zhat*c4rij
813 gezelter 439
814 chrisfen 987 ! calculate the potential
815     pot_term = ( qxx_j*(cx2*c3 - c2ri) + qyy_j*(cy2*c3 - c2ri) + &
816     qzz_j*(cz2*c3 - c2ri) )
817 chrisfen 959 vterm = pref * pot_term
818 chrisfen 740 vpair = vpair + vterm
819     epot = epot + sw*vterm
820 chrisfen 959
821 chrisfen 987 ! calculate derivatives for the forces and torques
822     dudx = dudx - preSw * ( &
823     qxx_j*(cx2*xhatc4 - (2.0_dp*cx_j*ux_j(1) + xhat)*c3ri) + &
824     qyy_j*(cy2*xhatc4 - (2.0_dp*cy_j*uy_j(1) + xhat)*c3ri) + &
825     qzz_j*(cz2*xhatc4 - (2.0_dp*cz_j*uz_j(1) + xhat)*c3ri) )
826     dudy = dudy - preSw * ( &
827     qxx_j*(cx2*yhatc4 - (2.0_dp*cx_j*ux_j(2) + yhat)*c3ri) + &
828     qyy_j*(cy2*yhatc4 - (2.0_dp*cy_j*uy_j(2) + yhat)*c3ri) + &
829     qzz_j*(cz2*yhatc4 - (2.0_dp*cz_j*uz_j(2) + yhat)*c3ri) )
830     dudz = dudz - preSw * ( &
831     qxx_j*(cx2*zhatc4 - (2.0_dp*cx_j*ux_j(3) + zhat)*c3ri) + &
832     qyy_j*(cy2*zhatc4 - (2.0_dp*cy_j*uy_j(3) + zhat)*c3ri) + &
833     qzz_j*(cz2*zhatc4 - (2.0_dp*cz_j*uz_j(3) + zhat)*c3ri) )
834 chrisfen 597
835 chrisfen 987 dudux_j(1) = dudux_j(1) + preSw*(qxx_j*cx_j*xhatdot2)
836     dudux_j(2) = dudux_j(2) + preSw*(qxx_j*cx_j*yhatdot2)
837     dudux_j(3) = dudux_j(3) + preSw*(qxx_j*cx_j*zhatdot2)
838 chrisfen 740
839 chrisfen 987 duduy_j(1) = duduy_j(1) + preSw*(qyy_j*cy_j*xhatdot2)
840     duduy_j(2) = duduy_j(2) + preSw*(qyy_j*cy_j*yhatdot2)
841     duduy_j(3) = duduy_j(3) + preSw*(qyy_j*cy_j*zhatdot2)
842 chrisfen 740
843 chrisfen 987 duduz_j(1) = duduz_j(1) + preSw*(qzz_j*cz_j*xhatdot2)
844     duduz_j(2) = duduz_j(2) + preSw*(qzz_j*cz_j*yhatdot2)
845     duduz_j(3) = duduz_j(3) + preSw*(qzz_j*cz_j*zhatdot2)
846 chrisfen 959
847 chrisfen 849
848 gezelter 439 endif
849 gezelter 411 endif
850 chrisfen 740
851 gezelter 411 if (i_is_Dipole) then
852 gezelter 507
853 gezelter 411 if (j_is_Charge) then
854 chrisfen 987 ! variables used by all the methods
855 gezelter 1286 pref = electroMult * pre12 * q_j * mu_i
856 chrisfen 987 preSw = sw*pref
857    
858 chrisfen 959 if (summationMethod .eq. REACTION_FIELD) then
859 gezelter 507
860 chrisfen 719 ri2 = riji * riji
861     ri3 = ri2 * riji
862    
863 chrisfen 700 vterm = pref * ct_i * ( ri2 - preRF2*rij )
864 chrisfen 696 vpair = vpair + vterm
865     epot = epot + sw*vterm
866    
867 chrisfen 987 dudx = dudx + preSw * ( ri3*(uz_i(1) - 3.0_dp*ct_i*xhat) - &
868 chrisfen 700 preRF2*uz_i(1) )
869 chrisfen 987 dudy = dudy + preSw * ( ri3*(uz_i(2) - 3.0_dp*ct_i*yhat) - &
870 chrisfen 700 preRF2*uz_i(2) )
871 chrisfen 987 dudz = dudz + preSw * ( ri3*(uz_i(3) - 3.0_dp*ct_i*zhat) - &
872 chrisfen 700 preRF2*uz_i(3) )
873 chrisfen 696
874 chrisfen 987 duduz_i(1) = duduz_i(1) + preSw * xhat * ( ri2 - preRF2*rij )
875     duduz_i(2) = duduz_i(2) + preSw * yhat * ( ri2 - preRF2*rij )
876     duduz_i(3) = duduz_i(3) + preSw * zhat * ( ri2 - preRF2*rij )
877 chrisfen 696
878 chrisfen 597 else
879 chrisfen 987 ! determine inverse r if we are using split dipoles
880 chrisfen 597 if (i_is_SplitDipole) then
881 chrisfen 959 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
882     ri = 1.0_dp / BigR
883 chrisfen 597 scale = rij * ri
884     else
885 gezelter 421 ri = riji
886 chrisfen 959 scale = 1.0_dp
887 gezelter 421 endif
888 chrisfen 987
889 chrisfen 597 sc2 = scale * scale
890 chrisfen 987
891     if (screeningMethod .eq. DAMPED) then
892     ! assemble the damping variables
893     call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
894     c1 = erfcVal*ri
895     c2 = (-derfcVal + c1)*ri
896     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
897     else
898     c1 = ri
899     c2 = c1*ri
900     c3 = 3.0_dp*c2*ri
901     endif
902    
903     c2ri = c2*ri
904 chrisfen 626
905 chrisfen 987 ! calculate the potential
906     pot_term = c2 * scale
907 chrisfen 849 vterm = pref * ct_i * pot_term
908 chrisfen 626 vpair = vpair + vterm
909     epot = epot + sw*vterm
910 chrisfen 959
911 chrisfen 987 ! calculate derivatives for the forces and torques
912     dudx = dudx + preSw * ( uz_i(1)*c2ri - ct_i*xhat*sc2*c3 )
913     dudy = dudy + preSw * ( uz_i(2)*c2ri - ct_i*yhat*sc2*c3 )
914     dudz = dudz + preSw * ( uz_i(3)*c2ri - ct_i*zhat*sc2*c3 )
915    
916     duduz_i(1) = duduz_i(1) + preSw * pot_term * xhat
917     duduz_i(2) = duduz_i(2) + preSw * pot_term * yhat
918     duduz_i(3) = duduz_i(3) + preSw * pot_term * zhat
919 chrisfen 959
920 gezelter 421 endif
921 chrisfen 597 endif
922 chrisfen 626
923 chrisfen 597 if (j_is_Dipole) then
924 chrisfen 987 ! variables used by all methods
925 chrisfen 719 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
926 gezelter 1286 pref = electroMult * pre22 * mu_i * mu_j
927 chrisfen 987 preSw = sw*pref
928 gezelter 421
929 chrisfen 710 if (summationMethod .eq. REACTION_FIELD) then
930 chrisfen 987 ri2 = riji * riji
931     ri3 = ri2 * riji
932     ri4 = ri2 * ri2
933    
934 chrisfen 959 vterm = pref*( ri3*(ct_ij - 3.0_dp * ct_i * ct_j) - &
935 chrisfen 695 preRF2*ct_ij )
936     vpair = vpair + vterm
937     epot = epot + sw*vterm
938    
939 chrisfen 959 a1 = 5.0_dp * ct_i * ct_j - ct_ij
940 chrisfen 695
941 chrisfen 987 dudx = dudx + preSw*3.0_dp*ri4*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
942     dudy = dudy + preSw*3.0_dp*ri4*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
943     dudz = dudz + preSw*3.0_dp*ri4*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
944 chrisfen 695
945 chrisfen 987 duduz_i(1) = duduz_i(1) + preSw*(ri3*(uz_j(1)-3.0_dp*ct_j*xhat) &
946 chrisfen 695 - preRF2*uz_j(1))
947 chrisfen 987 duduz_i(2) = duduz_i(2) + preSw*(ri3*(uz_j(2)-3.0_dp*ct_j*yhat) &
948 chrisfen 695 - preRF2*uz_j(2))
949 chrisfen 987 duduz_i(3) = duduz_i(3) + preSw*(ri3*(uz_j(3)-3.0_dp*ct_j*zhat) &
950 chrisfen 695 - preRF2*uz_j(3))
951 chrisfen 987 duduz_j(1) = duduz_j(1) + preSw*(ri3*(uz_i(1)-3.0_dp*ct_i*xhat) &
952 chrisfen 695 - preRF2*uz_i(1))
953 chrisfen 987 duduz_j(2) = duduz_j(2) + preSw*(ri3*(uz_i(2)-3.0_dp*ct_i*yhat) &
954 chrisfen 695 - preRF2*uz_i(2))
955 chrisfen 987 duduz_j(3) = duduz_j(3) + preSw*(ri3*(uz_i(3)-3.0_dp*ct_i*zhat) &
956 chrisfen 695 - preRF2*uz_i(3))
957    
958 chrisfen 597 else
959     if (i_is_SplitDipole) then
960     if (j_is_SplitDipole) then
961 chrisfen 959 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
962 chrisfen 597 else
963 chrisfen 959 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
964 chrisfen 597 endif
965 chrisfen 959 ri = 1.0_dp / BigR
966 chrisfen 597 scale = rij * ri
967     else
968     if (j_is_SplitDipole) then
969 chrisfen 959 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
970     ri = 1.0_dp / BigR
971 chrisfen 597 scale = rij * ri
972     else
973     ri = riji
974 chrisfen 959 scale = 1.0_dp
975 chrisfen 597 endif
976     endif
977 chrisfen 719
978 chrisfen 987 if (screeningMethod .eq. DAMPED) then
979     ! assemble the damping variables
980     call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
981     c1 = erfcVal*ri
982     c2 = (-derfcVal + c1)*ri
983     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
984     c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*ri*ri
985     else
986     c1 = ri
987     c2 = c1*ri
988     c3 = 3.0_dp*c2*ri
989     c4 = 5.0_dp*c3*ri*ri
990     endif
991    
992     ! precompute variables for convenience
993 chrisfen 986 sc2 = scale * scale
994 chrisfen 987 cti3 = ct_i*sc2*c3
995     ctj3 = ct_j*sc2*c3
996 chrisfen 986 ctidotj = ct_i * ct_j * sc2
997 chrisfen 987 preSwSc = preSw*scale
998     c2ri = c2*ri
999     c3ri = c3*ri
1000     c4rij = c4*rij
1001 chrisfen 959
1002 chrisfen 987
1003 chrisfen 986 ! calculate the potential
1004 chrisfen 987 pot_term = (ct_ij*c2ri - ctidotj*c3)
1005     vterm = pref * pot_term
1006 chrisfen 986 vpair = vpair + vterm
1007     epot = epot + sw*vterm
1008 chrisfen 959
1009 chrisfen 986 ! calculate derivatives for the forces and torques
1010 chrisfen 987 dudx = dudx + preSwSc * ( ctidotj*xhat*c4rij - &
1011     (ct_i*uz_j(1) + ct_j*uz_i(1) + ct_ij*xhat)*c3ri )
1012     dudy = dudy + preSwSc * ( ctidotj*yhat*c4rij - &
1013     (ct_i*uz_j(2) + ct_j*uz_i(2) + ct_ij*yhat)*c3ri )
1014     dudz = dudz + preSwSc * ( ctidotj*zhat*c4rij - &
1015     (ct_i*uz_j(3) + ct_j*uz_i(3) + ct_ij*zhat)*c3ri )
1016 chrisfen 986
1017 chrisfen 987 duduz_i(1) = duduz_i(1) + preSw * ( uz_j(1)*c2ri - ctj3*xhat )
1018     duduz_i(2) = duduz_i(2) + preSw * ( uz_j(2)*c2ri - ctj3*yhat )
1019     duduz_i(3) = duduz_i(3) + preSw * ( uz_j(3)*c2ri - ctj3*zhat )
1020 chrisfen 597
1021 chrisfen 987 duduz_j(1) = duduz_j(1) + preSw * ( uz_i(1)*c2ri - cti3*xhat )
1022     duduz_j(2) = duduz_j(2) + preSw * ( uz_i(2)*c2ri - cti3*yhat )
1023     duduz_j(3) = duduz_j(3) + preSw * ( uz_i(3)*c2ri - cti3*zhat )
1024 chrisfen 849
1025 chrisfen 597 endif
1026 gezelter 411 endif
1027     endif
1028 gezelter 439
1029     if (i_is_Quadrupole) then
1030     if (j_is_Charge) then
1031 chrisfen 1022 ! precompute some necessary variables
1032     cx2 = cx_i * cx_i
1033     cy2 = cy_i * cy_i
1034     cz2 = cz_i * cz_i
1035 gezelter 1286 pref = electroMult * pre14 * q_j * one_third
1036 chrisfen 1022
1037 chrisfen 849 if (screeningMethod .eq. DAMPED) then
1038 chrisfen 959 ! assemble the damping variables
1039 chrisfen 987 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
1040     c1 = erfcVal*riji
1041     c2 = (-derfcVal + c1)*riji
1042     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*riji
1043     c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*riji*riji
1044     else
1045     c1 = riji
1046     c2 = c1*riji
1047     c3 = 3.0_dp*c2*riji
1048     c4 = 5.0_dp*c3*riji*riji
1049 chrisfen 849 endif
1050 chrisfen 987
1051 chrisfen 1022 ! precompute some variables for convenience
1052 chrisfen 987 preSw = sw*pref
1053     c2ri = c2*riji
1054     c3ri = c3*riji
1055     c4rij = c4*rij
1056     xhatdot2 = 2.0_dp*xhat*c3
1057     yhatdot2 = 2.0_dp*yhat*c3
1058     zhatdot2 = 2.0_dp*zhat*c3
1059     xhatc4 = xhat*c4rij
1060     yhatc4 = yhat*c4rij
1061     zhatc4 = zhat*c4rij
1062 chrisfen 959
1063 chrisfen 1022 ! calculate the potential
1064     pot_term = ( qxx_i * (cx2*c3 - c2ri) + qyy_i * (cy2*c3 - c2ri) + &
1065     qzz_i * (cz2*c3 - c2ri) )
1066    
1067     vterm = pref * pot_term
1068     vpair = vpair + vterm
1069     epot = epot + sw*vterm
1070    
1071 chrisfen 987 ! calculate the derivatives for the forces and torques
1072     dudx = dudx - preSw * ( &
1073     qxx_i*(cx2*xhatc4 - (2.0_dp*cx_i*ux_i(1) + xhat)*c3ri) + &
1074     qyy_i*(cy2*xhatc4 - (2.0_dp*cy_i*uy_i(1) + xhat)*c3ri) + &
1075     qzz_i*(cz2*xhatc4 - (2.0_dp*cz_i*uz_i(1) + xhat)*c3ri) )
1076     dudy = dudy - preSw * ( &
1077     qxx_i*(cx2*yhatc4 - (2.0_dp*cx_i*ux_i(2) + yhat)*c3ri) + &
1078     qyy_i*(cy2*yhatc4 - (2.0_dp*cy_i*uy_i(2) + yhat)*c3ri) + &
1079     qzz_i*(cz2*yhatc4 - (2.0_dp*cz_i*uz_i(2) + yhat)*c3ri) )
1080     dudz = dudz - preSw * ( &
1081     qxx_i*(cx2*zhatc4 - (2.0_dp*cx_i*ux_i(3) + zhat)*c3ri) + &
1082     qyy_i*(cy2*zhatc4 - (2.0_dp*cy_i*uy_i(3) + zhat)*c3ri) + &
1083     qzz_i*(cz2*zhatc4 - (2.0_dp*cz_i*uz_i(3) + zhat)*c3ri) )
1084 chrisfen 740
1085 chrisfen 987 dudux_i(1) = dudux_i(1) + preSw*(qxx_i*cx_i*xhatdot2)
1086     dudux_i(2) = dudux_i(2) + preSw*(qxx_i*cx_i*yhatdot2)
1087     dudux_i(3) = dudux_i(3) + preSw*(qxx_i*cx_i*zhatdot2)
1088 chrisfen 740
1089 chrisfen 987 duduy_i(1) = duduy_i(1) + preSw*(qyy_i*cy_i*xhatdot2)
1090     duduy_i(2) = duduy_i(2) + preSw*(qyy_i*cy_i*yhatdot2)
1091     duduy_i(3) = duduy_i(3) + preSw*(qyy_i*cy_i*zhatdot2)
1092 chrisfen 740
1093 chrisfen 987 duduz_i(1) = duduz_i(1) + preSw*(qzz_i*cz_i*xhatdot2)
1094     duduz_i(2) = duduz_i(2) + preSw*(qzz_i*cz_i*yhatdot2)
1095     duduz_i(3) = duduz_i(3) + preSw*(qzz_i*cz_i*zhatdot2)
1096 gezelter 439 endif
1097     endif
1098 gezelter 507
1099 gezelter 1386 pot = pot + epot
1100 gezelter 507
1101 gezelter 1386 f1(1) = f1(1) + dudx
1102     f1(2) = f1(2) + dudy
1103     f1(3) = f1(3) + dudz
1104 gezelter 507
1105 gezelter 411 if (i_is_Dipole .or. i_is_Quadrupole) then
1106 gezelter 1386 t1(1) = t1(1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1107     t1(2) = t1(2) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1108     t1(3) = t1(3) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1109 gezelter 411 endif
1110 gezelter 439 if (i_is_Quadrupole) then
1111 gezelter 1386 t1(1) = t1(1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
1112     t1(2) = t1(2) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
1113     t1(3) = t1(3) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
1114 gezelter 411
1115 gezelter 1386 t1(1) = t1(1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
1116     t1(2) = t1(2) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
1117     t1(3) = t1(3) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
1118 gezelter 439 endif
1119    
1120 gezelter 411 if (j_is_Dipole .or. j_is_Quadrupole) then
1121 gezelter 1386 t2(1) = t2(1) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1122     t2(2) = t2(2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1123     t2(3) = t2(3) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1124 gezelter 411 endif
1125 gezelter 439 if (j_is_Quadrupole) then
1126 gezelter 1386 t2(1) = t2(1) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
1127     t2(2) = t2(2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
1128     t2(3) = t2(3) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
1129 gezelter 411
1130 gezelter 1386 t2(1) = t2(1) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
1131     t2(2) = t2(2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
1132     t2(3) = t2(3) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
1133 gezelter 439 endif
1134    
1135 gezelter 411 return
1136     end subroutine doElectrostaticPair
1137 chuckv 492
1138     subroutine destroyElectrostaticTypes()
1139    
1140 gezelter 507 if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap)
1141    
1142 chuckv 492 end subroutine destroyElectrostaticTypes
1143 gezelter 1340
1144 gezelter 1464 subroutine self_self(atom1, eFrame, skch, mypot, t)
1145 chrisfen 682 integer, intent(in) :: atom1
1146 chrisfen 695 integer :: atid1
1147 chrisfen 682 real(kind=dp), dimension(9,nLocal) :: eFrame
1148 chrisfen 695 real(kind=dp), dimension(3,nLocal) :: t
1149 gezelter 1340 real(kind=dp) :: mu1, chg1, c1e
1150     real(kind=dp) :: preVal, epot, mypot, skch
1151     real(kind=dp) :: eix, eiy, eiz, self
1152 chrisfen 682
1153 chrisfen 695 ! this is a local only array, so we use the local atom type id's:
1154     atid1 = atid(atom1)
1155 chrisfen 703
1156     if (.not.summationMethodChecked) then
1157     call checkSummationMethod()
1158     endif
1159 chrisfen 695
1160 chrisfen 703 if (summationMethod .eq. REACTION_FIELD) then
1161     if (ElectrostaticMap(atid1)%is_Dipole) then
1162     mu1 = getDipoleMoment(atid1)
1163    
1164     preVal = pre22 * preRF2 * mu1*mu1
1165 chrisfen 959 mypot = mypot - 0.5_dp*preVal
1166 chrisfen 703
1167     ! The self-correction term adds into the reaction field vector
1168    
1169     eix = preVal * eFrame(3,atom1)
1170     eiy = preVal * eFrame(6,atom1)
1171     eiz = preVal * eFrame(9,atom1)
1172    
1173     ! once again, this is self-self, so only the local arrays are needed
1174     ! even for MPI jobs:
1175    
1176     t(1,atom1)=t(1,atom1) - eFrame(6,atom1)*eiz + &
1177     eFrame(9,atom1)*eiy
1178     t(2,atom1)=t(2,atom1) - eFrame(9,atom1)*eix + &
1179     eFrame(3,atom1)*eiz
1180     t(3,atom1)=t(3,atom1) - eFrame(3,atom1)*eiy + &
1181     eFrame(6,atom1)*eix
1182    
1183     endif
1184    
1185 chrisfen 743 elseif ( (summationMethod .eq. SHIFTED_FORCE) .or. &
1186     (summationMethod .eq. SHIFTED_POTENTIAL) ) then
1187 chrisfen 717 if (ElectrostaticMap(atid1)%is_Charge) then
1188 gezelter 1340 chg1 = getCharge(atid1)
1189 chrisfen 717 if (screeningMethod .eq. DAMPED) then
1190 gezelter 1340 self = - 0.5_dp * (c1c+alphaPi) * chg1 * (chg1+skch) * pre11
1191 chrisfen 717 else
1192 gezelter 1340 self = - 0.5_dp * rcuti * chg1 * (chg1+skch) * pre11
1193 chrisfen 717 endif
1194 gezelter 1340
1195     mypot = mypot + self
1196 chrisfen 717 endif
1197     endif
1198 chrisfen 695
1199 gezelter 1337
1200    
1201 chrisfen 682 return
1202 chrisfen 703 end subroutine self_self
1203 chrisfen 682
1204 gezelter 1286 subroutine rf_self_excludes(atom1, atom2, sw, electroMult, eFrame, d, &
1205 gezelter 1464 rij, vpair, myPot, f, t)
1206 chrisfen 700 integer, intent(in) :: atom1
1207     integer, intent(in) :: atom2
1208     logical :: i_is_Charge, j_is_Charge
1209     logical :: i_is_Dipole, j_is_Dipole
1210     integer :: atid1
1211     integer :: atid2
1212     real(kind=dp), intent(in) :: rij
1213 gezelter 1286 real(kind=dp), intent(in) :: sw, electroMult
1214 chrisfen 700 real(kind=dp), intent(in), dimension(3) :: d
1215     real(kind=dp), intent(inout) :: vpair
1216     real(kind=dp), dimension(9,nLocal) :: eFrame
1217     real(kind=dp), dimension(3,nLocal) :: f
1218     real(kind=dp), dimension(3,nLocal) :: t
1219     real (kind = dp), dimension(3) :: duduz_i
1220     real (kind = dp), dimension(3) :: duduz_j
1221     real (kind = dp), dimension(3) :: uz_i
1222     real (kind = dp), dimension(3) :: uz_j
1223     real(kind=dp) :: q_i, q_j, mu_i, mu_j
1224     real(kind=dp) :: xhat, yhat, zhat
1225     real(kind=dp) :: ct_i, ct_j
1226     real(kind=dp) :: ri2, ri3, riji, vterm
1227 chrisfen 703 real(kind=dp) :: pref, preVal, rfVal, myPot
1228 chrisfen 700 real(kind=dp) :: dudx, dudy, dudz, dudr
1229    
1230 chrisfen 703 if (.not.summationMethodChecked) then
1231     call checkSummationMethod()
1232 chrisfen 700 endif
1233    
1234 gezelter 939 dudx = zero
1235     dudy = zero
1236     dudz = zero
1237 chrisfen 700
1238 chrisfen 959 riji = 1.0_dp/rij
1239 chrisfen 700
1240     xhat = d(1) * riji
1241     yhat = d(2) * riji
1242     zhat = d(3) * riji
1243    
1244     ! this is a local only array, so we use the local atom type id's:
1245     atid1 = atid(atom1)
1246     atid2 = atid(atom2)
1247     i_is_Charge = ElectrostaticMap(atid1)%is_Charge
1248     j_is_Charge = ElectrostaticMap(atid2)%is_Charge
1249     i_is_Dipole = ElectrostaticMap(atid1)%is_Dipole
1250     j_is_Dipole = ElectrostaticMap(atid2)%is_Dipole
1251    
1252     if (i_is_Charge.and.j_is_Charge) then
1253     q_i = ElectrostaticMap(atid1)%charge
1254     q_j = ElectrostaticMap(atid2)%charge
1255    
1256 gezelter 1286 preVal = electroMult * pre11 * q_i * q_j
1257 chrisfen 700 rfVal = preRF*rij*rij
1258     vterm = preVal * rfVal
1259    
1260 chrisfen 703 myPot = myPot + sw*vterm
1261    
1262 chrisfen 959 dudr = sw*preVal * 2.0_dp*rfVal*riji
1263 chrisfen 703
1264 chrisfen 700 dudx = dudx + dudr * xhat
1265     dudy = dudy + dudr * yhat
1266     dudz = dudz + dudr * zhat
1267 chrisfen 703
1268 chrisfen 700 elseif (i_is_Charge.and.j_is_Dipole) then
1269     q_i = ElectrostaticMap(atid1)%charge
1270     mu_j = ElectrostaticMap(atid2)%dipole_moment
1271     uz_j(1) = eFrame(3,atom2)
1272     uz_j(2) = eFrame(6,atom2)
1273     uz_j(3) = eFrame(9,atom2)
1274     ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
1275 chrisfen 703
1276 chrisfen 700 ri2 = riji * riji
1277     ri3 = ri2 * riji
1278    
1279 gezelter 1286 pref = electroMult * pre12 * q_i * mu_j
1280 chrisfen 700 vterm = - pref * ct_j * ( ri2 - preRF2*rij )
1281 chrisfen 703 myPot = myPot + sw*vterm
1282    
1283 chrisfen 959 dudx = dudx - sw*pref*( ri3*(uz_j(1)-3.0_dp*ct_j*xhat) &
1284 chrisfen 703 - preRF2*uz_j(1) )
1285 chrisfen 959 dudy = dudy - sw*pref*( ri3*(uz_j(2)-3.0_dp*ct_j*yhat) &
1286 chrisfen 703 - preRF2*uz_j(2) )
1287 chrisfen 959 dudz = dudz - sw*pref*( ri3*(uz_j(3)-3.0_dp*ct_j*zhat) &
1288 chrisfen 703 - preRF2*uz_j(3) )
1289    
1290 chrisfen 700 duduz_j(1) = duduz_j(1) - sw * pref * xhat * ( ri2 - preRF2*rij )
1291     duduz_j(2) = duduz_j(2) - sw * pref * yhat * ( ri2 - preRF2*rij )
1292     duduz_j(3) = duduz_j(3) - sw * pref * zhat * ( ri2 - preRF2*rij )
1293 chrisfen 703
1294 chrisfen 700 elseif (i_is_Dipole.and.j_is_Charge) then
1295     mu_i = ElectrostaticMap(atid1)%dipole_moment
1296     q_j = ElectrostaticMap(atid2)%charge
1297     uz_i(1) = eFrame(3,atom1)
1298     uz_i(2) = eFrame(6,atom1)
1299     uz_i(3) = eFrame(9,atom1)
1300     ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
1301 chrisfen 703
1302 chrisfen 700 ri2 = riji * riji
1303     ri3 = ri2 * riji
1304    
1305 gezelter 1286 pref = electroMult * pre12 * q_j * mu_i
1306 chrisfen 700 vterm = pref * ct_i * ( ri2 - preRF2*rij )
1307 chrisfen 703 myPot = myPot + sw*vterm
1308 chrisfen 700
1309 chrisfen 959 dudx = dudx + sw*pref*( ri3*(uz_i(1)-3.0_dp*ct_i*xhat) &
1310 chrisfen 703 - preRF2*uz_i(1) )
1311 chrisfen 959 dudy = dudy + sw*pref*( ri3*(uz_i(2)-3.0_dp*ct_i*yhat) &
1312 chrisfen 703 - preRF2*uz_i(2) )
1313 chrisfen 959 dudz = dudz + sw*pref*( ri3*(uz_i(3)-3.0_dp*ct_i*zhat) &
1314 chrisfen 703 - preRF2*uz_i(3) )
1315 chrisfen 700
1316     duduz_i(1) = duduz_i(1) + sw * pref * xhat * ( ri2 - preRF2*rij )
1317     duduz_i(2) = duduz_i(2) + sw * pref * yhat * ( ri2 - preRF2*rij )
1318     duduz_i(3) = duduz_i(3) + sw * pref * zhat * ( ri2 - preRF2*rij )
1319    
1320     endif
1321 chrisfen 703
1322    
1323     ! accumulate the forces and torques resulting from the self term
1324 chrisfen 700 f(1,atom1) = f(1,atom1) + dudx
1325     f(2,atom1) = f(2,atom1) + dudy
1326     f(3,atom1) = f(3,atom1) + dudz
1327    
1328     f(1,atom2) = f(1,atom2) - dudx
1329     f(2,atom2) = f(2,atom2) - dudy
1330     f(3,atom2) = f(3,atom2) - dudz
1331    
1332     if (i_is_Dipole) then
1333     t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1334     t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1335     t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1336     elseif (j_is_Dipole) then
1337     t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1338     t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1339     t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1340     endif
1341    
1342     return
1343     end subroutine rf_self_excludes
1344    
1345 chrisfen 998 subroutine accumulate_box_dipole(atom1, eFrame, d, pChg, nChg, pChgPos, &
1346     nChgPos, dipVec, pChgCount, nChgCount)
1347     integer, intent(in) :: atom1
1348     logical :: i_is_Charge
1349     logical :: i_is_Dipole
1350     integer :: atid1
1351     integer :: pChgCount
1352     integer :: nChgCount
1353     real(kind=dp), intent(in), dimension(3) :: d
1354     real(kind=dp), dimension(9,nLocal) :: eFrame
1355     real(kind=dp) :: pChg
1356     real(kind=dp) :: nChg
1357     real(kind=dp), dimension(3) :: pChgPos
1358     real(kind=dp), dimension(3) :: nChgPos
1359     real(kind=dp), dimension(3) :: dipVec
1360     real(kind=dp), dimension(3) :: uz_i
1361     real(kind=dp), dimension(3) :: pos
1362     real(kind=dp) :: q_i, mu_i
1363     real(kind=dp) :: pref, preVal
1364    
1365     if (.not.summationMethodChecked) then
1366     call checkSummationMethod()
1367     endif
1368    
1369     ! this is a local only array, so we use the local atom type id's:
1370     atid1 = atid(atom1)
1371     i_is_Charge = ElectrostaticMap(atid1)%is_Charge
1372     i_is_Dipole = ElectrostaticMap(atid1)%is_Dipole
1373    
1374     if (i_is_Charge) then
1375     q_i = ElectrostaticMap(atid1)%charge
1376     ! convert to the proper units
1377     q_i = q_i * chargeToC
1378     pos = d * angstromToM
1379    
1380     if (q_i.le.0.0_dp) then
1381     nChg = nChg - q_i
1382     nChgPos(1) = nChgPos(1) + pos(1)
1383     nChgPos(2) = nChgPos(2) + pos(2)
1384     nChgPos(3) = nChgPos(3) + pos(3)
1385     nChgCount = nChgCount + 1
1386    
1387     else
1388     pChg = pChg + q_i
1389     pChgPos(1) = pChgPos(1) + pos(1)
1390     pChgPos(2) = pChgPos(2) + pos(2)
1391     pChgPos(3) = pChgPos(3) + pos(3)
1392     pChgCount = pChgCount + 1
1393    
1394     endif
1395    
1396     endif
1397    
1398     if (i_is_Dipole) then
1399     mu_i = ElectrostaticMap(atid1)%dipole_moment
1400     uz_i(1) = eFrame(3,atom1)
1401     uz_i(2) = eFrame(6,atom1)
1402     uz_i(3) = eFrame(9,atom1)
1403     ! convert to the proper units
1404     mu_i = mu_i * debyeToCm
1405    
1406     dipVec(1) = dipVec(1) + uz_i(1)*mu_i
1407     dipVec(2) = dipVec(2) + uz_i(2)*mu_i
1408     dipVec(3) = dipVec(3) + uz_i(3)*mu_i
1409    
1410     endif
1411    
1412     return
1413     end subroutine accumulate_box_dipole
1414    
1415 gezelter 411 end module electrostatic_module

Properties

Name Value
svn:keywords Author Id Revision Date