ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 3559
Committed: Fri Oct 23 18:41:09 2009 UTC (15 years, 9 months ago) by gezelter
File size: 49159 byte(s)
Log Message:
removing MPI responsibilities from the lowest level force routines.  This is
in preparation for migrating these routines (LJ, electrostatic, eam, suttonchen,
gay-berne, sticky) to C++

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 chrisfen 2715 use interpolation
51 gezelter 2095 implicit none
52    
53     PRIVATE
54    
55 chuckv 2355
56 gezelter 2301 #define __FORTRAN90
57 chuckv 2355 #include "UseTheForce/DarkSide/fInteractionMap.h"
58 gezelter 2301 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59 chrisfen 2415 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
60 gezelter 2301
61 chuckv 2355
62 gezelter 2118 !! 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 2755 real(kind=dp), parameter :: pre11 = 332.0637778_dp
66 gezelter 2118 !! Charge-Dipole, assuming charges are measured in electrons, and
67     !! dipoles are measured in debyes
68 chrisfen 2755 real(kind=dp), parameter :: pre12 = 69.13373_dp
69 gezelter 2118 !! Dipole-Dipole, assuming dipoles are measured in debyes
70 chrisfen 2755 real(kind=dp), parameter :: pre22 = 14.39325_dp
71 gezelter 2118 !! 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 2755 real(kind=dp), parameter :: pre14 = 69.13373_dp
75 gezelter 2095
76 chrisfen 2755 real(kind=dp), parameter :: zero = 0.0_dp
77 chrisfen 2724
78 chrisfen 2917 !! 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 2724 !! number of points for electrostatic splines
84     integer, parameter :: np = 100
85 gezelter 2722
86 chrisfen 2411 !! variables to handle different summation methods for long-range
87     !! electrostatics:
88 gezelter 2301 integer, save :: summationMethod = NONE
89 chrisfen 2409 integer, save :: screeningMethod = UNDAMPED
90 chrisfen 2302 logical, save :: summationMethodChecked = .false.
91 gezelter 2301 real(kind=DP), save :: defaultCutoff = 0.0_DP
92 chrisfen 2381 real(kind=DP), save :: defaultCutoff2 = 0.0_DP
93 gezelter 2301 logical, save :: haveDefaultCutoff = .false.
94     real(kind=DP), save :: dampingAlpha = 0.0_DP
95 chrisfen 2843 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 2301 logical, save :: haveDampingAlpha = .false.
100 chrisfen 2381 real(kind=DP), save :: dielectric = 1.0_DP
101 gezelter 2301 logical, save :: haveDielectric = .false.
102     real(kind=DP), save :: constEXP = 0.0_DP
103 chrisfen 2381 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 2394 real(kind=dp), save :: preRF2 = 0.0_DP
114 chrisfen 2843 real(kind=dp), save :: erfcVal = 1.0_DP
115     real(kind=dp), save :: derfcVal = 0.0_DP
116     type(cubicSpline), save :: erfcSpline
117 chrisfen 2724 logical, save :: haveElectroSpline = .false.
118 chrisfen 2843 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 2755 real(kind=dp), save :: one_third = 1.0_DP / 3.0_DP
131 chrisfen 2415
132 gezelter 2508 #if defined(__IFC) || defined(__PGI)
133 chuckv 2331 ! error function for ifc version > 7.
134 gezelter 2756 real(kind=dp), external :: erfc
135 chuckv 2331 #endif
136    
137 chrisfen 2552 public :: setElectrostaticSummationMethod
138 chrisfen 2409 public :: setScreeningMethod
139 gezelter 2301 public :: setElectrostaticCutoffRadius
140 chrisfen 2409 public :: setDampingAlpha
141 gezelter 2301 public :: setReactionFieldDielectric
142 chrisfen 2724 public :: buildElectroSpline
143 gezelter 2095 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 2189 public :: destroyElectrostaticTypes
152 chrisfen 2402 public :: self_self
153 chrisfen 2399 public :: rf_self_excludes
154 chrisfen 2917 public :: accumulate_box_dipole
155 gezelter 2095
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 2229 logical :: is_Tap = .false.
163 gezelter 2095 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 2717 logical, save :: hasElectrostaticMap
172    
173 gezelter 2095 contains
174    
175 chrisfen 2552 subroutine setElectrostaticSummationMethod(the_ESM)
176 gezelter 2301 integer, intent(in) :: the_ESM
177    
178     if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then
179 chrisfen 2552 call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method")
180 gezelter 2301 endif
181    
182 chrisfen 2309 summationMethod = the_ESM
183 chrisfen 2325
184 chrisfen 2552 end subroutine setElectrostaticSummationMethod
185 gezelter 2301
186 chrisfen 2409 subroutine setScreeningMethod(the_SM)
187     integer, intent(in) :: the_SM
188     screeningMethod = the_SM
189     end subroutine setScreeningMethod
190    
191 chrisfen 2381 subroutine setElectrostaticCutoffRadius(thisRcut, thisRsw)
192 gezelter 2301 real(kind=dp), intent(in) :: thisRcut
193 chrisfen 2381 real(kind=dp), intent(in) :: thisRsw
194 gezelter 2301 defaultCutoff = thisRcut
195 chrisfen 2548 defaultCutoff2 = defaultCutoff*defaultCutoff
196 chrisfen 2381 rrf = defaultCutoff
197     rt = thisRsw
198 gezelter 2301 haveDefaultCutoff = .true.
199     end subroutine setElectrostaticCutoffRadius
200    
201 chrisfen 2409 subroutine setDampingAlpha(thisAlpha)
202 gezelter 2301 real(kind=dp), intent(in) :: thisAlpha
203     dampingAlpha = thisAlpha
204 chrisfen 2415 alpha2 = dampingAlpha*dampingAlpha
205 chrisfen 2843 alpha4 = alpha2*alpha2
206     alpha6 = alpha4*alpha2
207     alpha8 = alpha4*alpha4
208 gezelter 2301 haveDampingAlpha = .true.
209 chrisfen 2409 end subroutine setDampingAlpha
210 gezelter 2301
211     subroutine setReactionFieldDielectric(thisDielectric)
212     real(kind=dp), intent(in) :: thisDielectric
213     dielectric = thisDielectric
214     haveDielectric = .true.
215     end subroutine setReactionFieldDielectric
216    
217 chrisfen 2724 subroutine buildElectroSpline()
218     real( kind = dp ), dimension(np) :: xvals, yvals
219     real( kind = dp ) :: dx, rmin, rval
220     integer :: i
221 chrisfen 2715
222 chrisfen 2755 rmin = 0.0_dp
223 chrisfen 2724
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 2756 yvals(i) = erfc(dampingAlpha*rval)
230 chrisfen 2724 enddo
231    
232 chrisfen 2843 call newSpline(erfcSpline, xvals, yvals, .true.)
233 chrisfen 2724
234     haveElectroSpline = .true.
235     end subroutine buildElectroSpline
236    
237 gezelter 2095 subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
238 chrisfen 2229 is_SplitDipole, is_Quadrupole, is_Tap, status)
239 gezelter 2204
240 gezelter 2095 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 2229 logical, intent(in) :: is_Tap
246 gezelter 2095 integer, intent(out) :: status
247     integer :: nAtypes, myATID, i, j
248    
249     status = 0
250     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
251 gezelter 2204
252 gezelter 2095 !! 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 2204
257 gezelter 2095 nAtypes = getSize(atypes)
258 gezelter 2204
259 gezelter 2095 if (nAtypes == 0) then
260     status = -1
261     return
262     end if
263 gezelter 2204
264 gezelter 2717 allocate(ElectrostaticMap(nAtypes))
265 gezelter 2204
266 gezelter 2095 end if
267    
268     if (myATID .gt. size(ElectrostaticMap)) then
269     status = -1
270     return
271     endif
272 gezelter 2204
273 gezelter 2095 ! 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 2229 ElectrostaticMap(myATID)%is_Tap = is_Tap
281 gezelter 2204
282 gezelter 2717 hasElectrostaticMap = .true.
283    
284 gezelter 2095 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 2717 if (.not.hasElectrostaticMap) then
296 gezelter 2095 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 2204 endif
312 gezelter 2095
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 2717 if (.not.hasElectrostaticMap) then
326 gezelter 2095 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 2717 if (.not.hasElectrostaticMap) then
356 gezelter 2095 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 2717 if (.not.hasElectrostaticMap) then
386 gezelter 2095 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 2204
403 gezelter 2095 do i = 1, 3
404 gezelter 2204 ElectrostaticMap(myATID)%quadrupole_moments(i) = &
405     quadrupole_moments(i)
406     enddo
407 gezelter 2095
408     end subroutine setQuadrupoleMoments
409    
410 gezelter 2204
411 gezelter 2095 function getCharge(atid) result (c)
412     integer, intent(in) :: atid
413     integer :: localError
414     real(kind=dp) :: c
415 gezelter 2204
416 gezelter 2717 if (.not.hasElectrostaticMap) then
417 gezelter 2095 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
418     return
419     end if
420 gezelter 2204
421 gezelter 2095 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 2204
426 gezelter 2095 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 2204
434 gezelter 2717 if (.not.hasElectrostaticMap) then
435 gezelter 2095 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
436     return
437     end if
438 gezelter 2204
439 gezelter 2095 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 2204
444 gezelter 2095 dm = ElectrostaticMap(atid)%dipole_moment
445     end function getDipoleMoment
446    
447 gezelter 2301 subroutine checkSummationMethod()
448    
449 chrisfen 2306 if (.not.haveDefaultCutoff) then
450     call handleError("checkSummationMethod", "no Default Cutoff set!")
451     endif
452    
453 chrisfen 2755 rcuti = 1.0_dp / defaultCutoff
454 chrisfen 2306 rcuti2 = rcuti*rcuti
455     rcuti3 = rcuti2*rcuti
456     rcuti4 = rcuti2*rcuti2
457    
458 chrisfen 2409 if (screeningMethod .eq. DAMPED) then
459 chrisfen 2402 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 2302
467 chrisfen 2548 constEXP = exp(-alpha2*defaultCutoff2)
468 chrisfen 2755 invRootPi = 0.56418958354775628695_dp
469     alphaPi = 2.0_dp*dampingAlpha*invRootPi
470 chrisfen 2843
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 2301 endif
485    
486 chrisfen 2302 if (summationMethod .eq. REACTION_FIELD) then
487 chrisfen 2402 if (haveDielectric) then
488     defaultCutoff2 = defaultCutoff*defaultCutoff
489 chrisfen 2755 preRF = (dielectric-1.0_dp) / &
490     ((2.0_dp*dielectric+1.0_dp)*defaultCutoff2*defaultCutoff)
491     preRF2 = 2.0_dp*preRF
492 chrisfen 2402 else
493     call handleError("checkSummationMethod", "Dielectric not set")
494 chrisfen 2302 endif
495 chrisfen 2402
496 chrisfen 2302 endif
497    
498 chrisfen 2724 if (.not.haveElectroSpline) then
499     call buildElectroSpline()
500     end if
501    
502 chrisfen 2302 summationMethodChecked = .true.
503 gezelter 2301 end subroutine checkSummationMethod
504    
505 chrisfen 2411
506 gezelter 3559 subroutine doElectrostaticPair(atom1, atom2, me1, me2, d, rij, r2, rcut, sw, &
507     electroMult, vpair, fpair, pot, eF1, eF2, f1, t1, t2, do_pot)
508 gezelter 2204
509 chrisfen 2402 logical, intent(in) :: do_pot
510 gezelter 2204
511 gezelter 3559 integer, intent(in) :: atom1, atom2, me1, me2
512 gezelter 2095 integer :: localError
513    
514 gezelter 3441 real(kind=dp), intent(in) :: rij, r2, sw, rcut, electroMult
515 gezelter 2095 real(kind=dp), intent(in), dimension(3) :: d
516     real(kind=dp), intent(inout) :: vpair
517 chrisfen 2402 real(kind=dp), intent(inout), dimension(3) :: fpair
518 gezelter 2095
519 chrisfen 2325 real( kind = dp ) :: pot
520 gezelter 3559 real( kind = dp ), dimension(9) :: eF1, eF2 ! eFrame = electroFrame
521     real( kind = dp ), dimension(3) :: f1
522 chrisfen 2409 real( kind = dp ), dimension(3,nLocal) :: felec
523 gezelter 3559 real( kind = dp ), dimension(3) :: t1, t2
524 gezelter 2204
525 gezelter 2123 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
526     real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
527     real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
528     real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
529 gezelter 2095
530     logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
531     logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
532 chrisfen 2229 logical :: i_is_Tap, j_is_Tap
533 gezelter 3559 integer :: id1, id2
534 gezelter 2095 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
535 gezelter 2123 real (kind=dp) :: qxx_i, qyy_i, qzz_i
536     real (kind=dp) :: qxx_j, qyy_j, qzz_j
537     real (kind=dp) :: cx_i, cy_i, cz_i
538     real (kind=dp) :: cx_j, cy_j, cz_j
539     real (kind=dp) :: cx2, cy2, cz2
540 chrisfen 2418 real (kind=dp) :: ct_i, ct_j, ct_ij, a0, a1
541 gezelter 2105 real (kind=dp) :: riji, ri, ri2, ri3, ri4
542 chrisfen 2296 real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2
543 gezelter 2105 real (kind=dp) :: xhat, yhat, zhat
544 gezelter 2095 real (kind=dp) :: dudx, dudy, dudz
545 chrisfen 2325 real (kind=dp) :: scale, sc2, bigR
546 chrisfen 2415 real (kind=dp) :: varEXP
547 chrisfen 2418 real (kind=dp) :: pot_term
548 chrisfen 2394 real (kind=dp) :: preVal, rfVal
549 chrisfen 2843 real (kind=dp) :: c2ri, c3ri, c4rij
550 chrisfen 2755 real (kind=dp) :: cti3, ctj3, ctidotj
551 chrisfen 2843 real (kind=dp) :: preSw, preSwSc
552 chrisfen 2755 real (kind=dp) :: xhatdot2, yhatdot2, zhatdot2
553 chrisfen 2843 real (kind=dp) :: xhatc4, yhatc4, zhatc4
554 gezelter 2095
555 gezelter 2301 if (.not.summationMethodChecked) then
556     call checkSummationMethod()
557     endif
558    
559 gezelter 2095 !! some variables we'll need independent of electrostatic type:
560    
561 chrisfen 2755 riji = 1.0_dp / rij
562 chrisfen 2343
563 gezelter 2105 xhat = d(1) * riji
564     yhat = d(2) * riji
565     zhat = d(3) * riji
566 gezelter 2095
567     !! logicals
568     i_is_Charge = ElectrostaticMap(me1)%is_Charge
569     i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
570     i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
571     i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
572 chrisfen 2229 i_is_Tap = ElectrostaticMap(me1)%is_Tap
573 gezelter 2095
574     j_is_Charge = ElectrostaticMap(me2)%is_Charge
575     j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
576     j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
577     j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
578 chrisfen 2229 j_is_Tap = ElectrostaticMap(me2)%is_Tap
579 gezelter 2095
580     if (i_is_Charge) then
581     q_i = ElectrostaticMap(me1)%charge
582     endif
583 gezelter 2204
584 gezelter 2095 if (i_is_Dipole) then
585     mu_i = ElectrostaticMap(me1)%dipole_moment
586 gezelter 3559
587     uz_i(1) = eF1(3)
588     uz_i(2) = eF1(6)
589     uz_i(3) = eF1(9)
590    
591 gezelter 2123 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
592 gezelter 2095
593     if (i_is_SplitDipole) then
594     d_i = ElectrostaticMap(me1)%split_dipole_distance
595     endif
596 gezelter 2722 duduz_i = zero
597 gezelter 2095 endif
598    
599 gezelter 2123 if (i_is_Quadrupole) then
600     qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
601     qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
602     qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
603 gezelter 3559
604     ux_i(1) = eF1(1)
605     ux_i(2) = eF1(4)
606     ux_i(3) = eF1(7)
607     uy_i(1) = eF1(2)
608     uy_i(2) = eF1(5)
609     uy_i(3) = eF1(8)
610     uz_i(1) = eF1(3)
611     uz_i(2) = eF1(6)
612     uz_i(3) = eF1(9)
613    
614 gezelter 2123 cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
615     cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
616     cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
617 gezelter 2722 dudux_i = zero
618     duduy_i = zero
619     duduz_i = zero
620 gezelter 2123 endif
621    
622 gezelter 2095 if (j_is_Charge) then
623     q_j = ElectrostaticMap(me2)%charge
624     endif
625 gezelter 2204
626 gezelter 2095 if (j_is_Dipole) then
627     mu_j = ElectrostaticMap(me2)%dipole_moment
628 gezelter 3559
629     uz_j(1) = eF2(3)
630     uz_j(2) = eF2(6)
631     uz_j(3) = eF2(9)
632    
633 chrisfen 2162 ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
634 gezelter 2095
635     if (j_is_SplitDipole) then
636     d_j = ElectrostaticMap(me2)%split_dipole_distance
637     endif
638 gezelter 2722 duduz_j = zero
639 gezelter 2095 endif
640    
641 gezelter 2123 if (j_is_Quadrupole) then
642     qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
643     qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
644     qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
645 gezelter 3559
646     ux_j(1) = eF2(1)
647     ux_j(2) = eF2(4)
648     ux_j(3) = eF2(7)
649     uy_j(1) = eF2(2)
650     uy_j(2) = eF2(5)
651     uy_j(3) = eF2(8)
652     uz_j(1) = eF2(3)
653     uz_j(2) = eF2(6)
654     uz_j(3) = eF2(9)
655    
656 gezelter 2123 cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
657     cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
658     cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
659 gezelter 2722 dudux_j = zero
660     duduy_j = zero
661     duduz_j = zero
662 gezelter 2123 endif
663 chrisfen 2251
664 gezelter 2722 epot = zero
665     dudx = zero
666     dudy = zero
667     dudz = zero
668 gezelter 2095
669     if (i_is_Charge) then
670    
671     if (j_is_Charge) then
672 chrisfen 2438 if (screeningMethod .eq. DAMPED) then
673 chrisfen 2755 ! assemble the damping variables
674 chrisfen 2843 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
675     c1 = erfcVal*riji
676     c2 = (-derfcVal + c1)*riji
677     else
678     c1 = riji
679     c2 = c1*riji
680 chrisfen 2438 endif
681 gezelter 2204
682 gezelter 3441 preVal = electroMult * pre11 * q_i * q_j
683 chrisfen 2438
684 chrisfen 2409 if (summationMethod .eq. SHIFTED_POTENTIAL) then
685 chrisfen 2843 vterm = preVal * (c1 - c1c)
686 chrisfen 2296
687 chrisfen 2843 dudr = -sw * preVal * c2
688 chrisfen 2438
689 chrisfen 2409 elseif (summationMethod .eq. SHIFTED_FORCE) then
690 chrisfen 2843 vterm = preVal * ( c1 - c1c + c2c*(rij - defaultCutoff) )
691 chrisfen 2415
692 chrisfen 2843 dudr = sw * preVal * (c2c - c2)
693 chrisfen 2438
694 chrisfen 2394 elseif (summationMethod .eq. REACTION_FIELD) then
695 gezelter 3441 rfVal = electroMult * preRF*rij*rij
696 chrisfen 2394 vterm = preVal * ( riji + rfVal )
697 chrisfen 2399
698 chrisfen 2755 dudr = sw * preVal * ( 2.0_dp*rfVal - riji )*riji
699 chrisfen 2438
700 chrisfen 2296 else
701 chrisfen 2843 vterm = preVal * riji*erfcVal
702 chrisfen 2296
703 chrisfen 2843 dudr = - sw * preVal * c2
704 chrisfen 2438
705 chrisfen 2296 endif
706    
707 chrisfen 2438 vpair = vpair + vterm
708     epot = epot + sw*vterm
709    
710     dudx = dudx + dudr * xhat
711     dudy = dudy + dudr * yhat
712     dudz = dudz + dudr * zhat
713    
714 gezelter 2095 endif
715    
716     if (j_is_Dipole) then
717 chrisfen 2843 ! pref is used by all the possible methods
718 gezelter 3441 pref = electroMult * pre12 * q_i * mu_j
719 chrisfen 2843 preSw = sw*pref
720 gezelter 2095
721 chrisfen 2409 if (summationMethod .eq. REACTION_FIELD) then
722 chrisfen 2399 ri2 = riji * riji
723     ri3 = ri2 * riji
724 chrisfen 2395
725     vterm = - pref * ct_j * ( ri2 - preRF2*rij )
726     vpair = vpair + vterm
727     epot = epot + sw*vterm
728    
729 chrisfen 2843 dudx = dudx - preSw*( ri3*(uz_j(1) - 3.0_dp*ct_j*xhat) - &
730     preRF2*uz_j(1) )
731     dudy = dudy - preSw*( ri3*(uz_j(2) - 3.0_dp*ct_j*yhat) - &
732     preRF2*uz_j(2) )
733     dudz = dudz - preSw*( ri3*(uz_j(3) - 3.0_dp*ct_j*zhat) - &
734     preRF2*uz_j(3) )
735     duduz_j(1) = duduz_j(1) - preSw * xhat * ( ri2 - preRF2*rij )
736     duduz_j(2) = duduz_j(2) - preSw * yhat * ( ri2 - preRF2*rij )
737     duduz_j(3) = duduz_j(3) - preSw * zhat * ( ri2 - preRF2*rij )
738 chrisfen 2395
739 chrisfen 2296 else
740 chrisfen 2843 ! determine the inverse r used if we have split dipoles
741 chrisfen 2296 if (j_is_SplitDipole) then
742 chrisfen 2755 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
743     ri = 1.0_dp / BigR
744 chrisfen 2296 scale = rij * ri
745     else
746     ri = riji
747 chrisfen 2755 scale = 1.0_dp
748 chrisfen 2296 endif
749 chrisfen 2843
750 chrisfen 2296 sc2 = scale * scale
751 chrisfen 2325
752 chrisfen 2843 if (screeningMethod .eq. DAMPED) then
753     ! assemble the damping variables
754     call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
755     c1 = erfcVal*ri
756     c2 = (-derfcVal + c1)*ri
757     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
758     else
759     c1 = ri
760     c2 = c1*ri
761     c3 = 3.0_dp*c2*ri
762     endif
763    
764     c2ri = c2*ri
765    
766     ! calculate the potential
767     pot_term = scale * c2
768 chrisfen 2755 vterm = -pref * ct_j * pot_term
769 chrisfen 2325 vpair = vpair + vterm
770     epot = epot + sw*vterm
771 chrisfen 2296
772 chrisfen 2843 ! calculate derivatives for forces and torques
773     dudx = dudx - preSw*( uz_j(1)*c2ri - ct_j*xhat*sc2*c3 )
774     dudy = dudy - preSw*( uz_j(2)*c2ri - ct_j*yhat*sc2*c3 )
775     dudz = dudz - preSw*( uz_j(3)*c2ri - ct_j*zhat*sc2*c3 )
776 chrisfen 2548
777 chrisfen 2843 duduz_j(1) = duduz_j(1) - preSw * pot_term * xhat
778     duduz_j(2) = duduz_j(2) - preSw * pot_term * yhat
779     duduz_j(3) = duduz_j(3) - preSw * pot_term * zhat
780 gezelter 2095
781 chrisfen 2296 endif
782 gezelter 2095 endif
783 gezelter 2105
784 gezelter 2123 if (j_is_Quadrupole) then
785 chrisfen 2843 ! first precalculate some necessary variables
786     cx2 = cx_j * cx_j
787     cy2 = cy_j * cy_j
788     cz2 = cz_j * cz_j
789 gezelter 3441 pref = electroMult * pre14 * q_i * one_third
790 chrisfen 2843
791 chrisfen 2548 if (screeningMethod .eq. DAMPED) then
792 chrisfen 2755 ! assemble the damping variables
793 chrisfen 2843 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
794     c1 = erfcVal*riji
795     c2 = (-derfcVal + c1)*riji
796     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*riji
797     c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*riji*riji
798     else
799     c1 = riji
800     c2 = c1*riji
801     c3 = 3.0_dp*c2*riji
802     c4 = 5.0_dp*c3*riji*riji
803 chrisfen 2548 endif
804    
805 chrisfen 2843 ! precompute variables for convenience
806     preSw = sw*pref
807     c2ri = c2*riji
808     c3ri = c3*riji
809     c4rij = c4*rij
810     xhatdot2 = 2.0_dp*xhat*c3
811     yhatdot2 = 2.0_dp*yhat*c3
812     zhatdot2 = 2.0_dp*zhat*c3
813     xhatc4 = xhat*c4rij
814     yhatc4 = yhat*c4rij
815     zhatc4 = zhat*c4rij
816 gezelter 2123
817 chrisfen 2843 ! calculate the potential
818     pot_term = ( qxx_j*(cx2*c3 - c2ri) + qyy_j*(cy2*c3 - c2ri) + &
819     qzz_j*(cz2*c3 - c2ri) )
820 chrisfen 2755 vterm = pref * pot_term
821 chrisfen 2439 vpair = vpair + vterm
822     epot = epot + sw*vterm
823 chrisfen 2755
824 chrisfen 2843 ! calculate derivatives for the forces and torques
825     dudx = dudx - preSw * ( &
826     qxx_j*(cx2*xhatc4 - (2.0_dp*cx_j*ux_j(1) + xhat)*c3ri) + &
827     qyy_j*(cy2*xhatc4 - (2.0_dp*cy_j*uy_j(1) + xhat)*c3ri) + &
828     qzz_j*(cz2*xhatc4 - (2.0_dp*cz_j*uz_j(1) + xhat)*c3ri) )
829     dudy = dudy - preSw * ( &
830     qxx_j*(cx2*yhatc4 - (2.0_dp*cx_j*ux_j(2) + yhat)*c3ri) + &
831     qyy_j*(cy2*yhatc4 - (2.0_dp*cy_j*uy_j(2) + yhat)*c3ri) + &
832     qzz_j*(cz2*yhatc4 - (2.0_dp*cz_j*uz_j(2) + yhat)*c3ri) )
833     dudz = dudz - preSw * ( &
834     qxx_j*(cx2*zhatc4 - (2.0_dp*cx_j*ux_j(3) + zhat)*c3ri) + &
835     qyy_j*(cy2*zhatc4 - (2.0_dp*cy_j*uy_j(3) + zhat)*c3ri) + &
836     qzz_j*(cz2*zhatc4 - (2.0_dp*cz_j*uz_j(3) + zhat)*c3ri) )
837 chrisfen 2296
838 chrisfen 2843 dudux_j(1) = dudux_j(1) + preSw*(qxx_j*cx_j*xhatdot2)
839     dudux_j(2) = dudux_j(2) + preSw*(qxx_j*cx_j*yhatdot2)
840     dudux_j(3) = dudux_j(3) + preSw*(qxx_j*cx_j*zhatdot2)
841 chrisfen 2439
842 chrisfen 2843 duduy_j(1) = duduy_j(1) + preSw*(qyy_j*cy_j*xhatdot2)
843     duduy_j(2) = duduy_j(2) + preSw*(qyy_j*cy_j*yhatdot2)
844     duduy_j(3) = duduy_j(3) + preSw*(qyy_j*cy_j*zhatdot2)
845 chrisfen 2439
846 chrisfen 2843 duduz_j(1) = duduz_j(1) + preSw*(qzz_j*cz_j*xhatdot2)
847     duduz_j(2) = duduz_j(2) + preSw*(qzz_j*cz_j*yhatdot2)
848     duduz_j(3) = duduz_j(3) + preSw*(qzz_j*cz_j*zhatdot2)
849 chrisfen 2755
850 chrisfen 2548
851 gezelter 2123 endif
852 gezelter 2095 endif
853 chrisfen 2439
854 gezelter 2095 if (i_is_Dipole) then
855 gezelter 2204
856 gezelter 2095 if (j_is_Charge) then
857 chrisfen 2843 ! variables used by all the methods
858 gezelter 3441 pref = electroMult * pre12 * q_j * mu_i
859 chrisfen 2843 preSw = sw*pref
860    
861 chrisfen 2755 if (summationMethod .eq. REACTION_FIELD) then
862 gezelter 2204
863 chrisfen 2418 ri2 = riji * riji
864     ri3 = ri2 * riji
865    
866 chrisfen 2399 vterm = pref * ct_i * ( ri2 - preRF2*rij )
867 chrisfen 2395 vpair = vpair + vterm
868     epot = epot + sw*vterm
869    
870 chrisfen 2843 dudx = dudx + preSw * ( ri3*(uz_i(1) - 3.0_dp*ct_i*xhat) - &
871 chrisfen 2399 preRF2*uz_i(1) )
872 chrisfen 2843 dudy = dudy + preSw * ( ri3*(uz_i(2) - 3.0_dp*ct_i*yhat) - &
873 chrisfen 2399 preRF2*uz_i(2) )
874 chrisfen 2843 dudz = dudz + preSw * ( ri3*(uz_i(3) - 3.0_dp*ct_i*zhat) - &
875 chrisfen 2399 preRF2*uz_i(3) )
876 chrisfen 2395
877 chrisfen 2843 duduz_i(1) = duduz_i(1) + preSw * xhat * ( ri2 - preRF2*rij )
878     duduz_i(2) = duduz_i(2) + preSw * yhat * ( ri2 - preRF2*rij )
879     duduz_i(3) = duduz_i(3) + preSw * zhat * ( ri2 - preRF2*rij )
880 chrisfen 2395
881 chrisfen 2296 else
882 chrisfen 2843 ! determine inverse r if we are using split dipoles
883 chrisfen 2296 if (i_is_SplitDipole) then
884 chrisfen 2755 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
885     ri = 1.0_dp / BigR
886 chrisfen 2296 scale = rij * ri
887     else
888 gezelter 2105 ri = riji
889 chrisfen 2755 scale = 1.0_dp
890 gezelter 2105 endif
891 chrisfen 2843
892 chrisfen 2296 sc2 = scale * scale
893 chrisfen 2843
894     if (screeningMethod .eq. DAMPED) then
895     ! assemble the damping variables
896     call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
897     c1 = erfcVal*ri
898     c2 = (-derfcVal + c1)*ri
899     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
900     else
901     c1 = ri
902     c2 = c1*ri
903     c3 = 3.0_dp*c2*ri
904     endif
905    
906     c2ri = c2*ri
907 chrisfen 2325
908 chrisfen 2843 ! calculate the potential
909     pot_term = c2 * scale
910 chrisfen 2548 vterm = pref * ct_i * pot_term
911 chrisfen 2325 vpair = vpair + vterm
912     epot = epot + sw*vterm
913 chrisfen 2755
914 chrisfen 2843 ! calculate derivatives for the forces and torques
915     dudx = dudx + preSw * ( uz_i(1)*c2ri - ct_i*xhat*sc2*c3 )
916     dudy = dudy + preSw * ( uz_i(2)*c2ri - ct_i*yhat*sc2*c3 )
917     dudz = dudz + preSw * ( uz_i(3)*c2ri - ct_i*zhat*sc2*c3 )
918    
919     duduz_i(1) = duduz_i(1) + preSw * pot_term * xhat
920     duduz_i(2) = duduz_i(2) + preSw * pot_term * yhat
921     duduz_i(3) = duduz_i(3) + preSw * pot_term * zhat
922 chrisfen 2755
923 gezelter 2105 endif
924 chrisfen 2296 endif
925 chrisfen 2325
926 chrisfen 2296 if (j_is_Dipole) then
927 chrisfen 2843 ! variables used by all methods
928 chrisfen 2418 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
929 gezelter 3441 pref = electroMult * pre22 * mu_i * mu_j
930 chrisfen 2843 preSw = sw*pref
931 gezelter 2105
932 chrisfen 2409 if (summationMethod .eq. REACTION_FIELD) then
933 chrisfen 2843 ri2 = riji * riji
934     ri3 = ri2 * riji
935     ri4 = ri2 * ri2
936    
937 chrisfen 2755 vterm = pref*( ri3*(ct_ij - 3.0_dp * ct_i * ct_j) - &
938 chrisfen 2394 preRF2*ct_ij )
939     vpair = vpair + vterm
940     epot = epot + sw*vterm
941    
942 chrisfen 2755 a1 = 5.0_dp * ct_i * ct_j - ct_ij
943 chrisfen 2394
944 chrisfen 2843 dudx = dudx + preSw*3.0_dp*ri4*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
945     dudy = dudy + preSw*3.0_dp*ri4*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
946     dudz = dudz + preSw*3.0_dp*ri4*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
947 chrisfen 2394
948 chrisfen 2843 duduz_i(1) = duduz_i(1) + preSw*(ri3*(uz_j(1)-3.0_dp*ct_j*xhat) &
949 chrisfen 2394 - preRF2*uz_j(1))
950 chrisfen 2843 duduz_i(2) = duduz_i(2) + preSw*(ri3*(uz_j(2)-3.0_dp*ct_j*yhat) &
951 chrisfen 2394 - preRF2*uz_j(2))
952 chrisfen 2843 duduz_i(3) = duduz_i(3) + preSw*(ri3*(uz_j(3)-3.0_dp*ct_j*zhat) &
953 chrisfen 2394 - preRF2*uz_j(3))
954 chrisfen 2843 duduz_j(1) = duduz_j(1) + preSw*(ri3*(uz_i(1)-3.0_dp*ct_i*xhat) &
955 chrisfen 2394 - preRF2*uz_i(1))
956 chrisfen 2843 duduz_j(2) = duduz_j(2) + preSw*(ri3*(uz_i(2)-3.0_dp*ct_i*yhat) &
957 chrisfen 2394 - preRF2*uz_i(2))
958 chrisfen 2843 duduz_j(3) = duduz_j(3) + preSw*(ri3*(uz_i(3)-3.0_dp*ct_i*zhat) &
959 chrisfen 2394 - preRF2*uz_i(3))
960    
961 chrisfen 2296 else
962     if (i_is_SplitDipole) then
963     if (j_is_SplitDipole) then
964 chrisfen 2755 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
965 chrisfen 2296 else
966 chrisfen 2755 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
967 chrisfen 2296 endif
968 chrisfen 2755 ri = 1.0_dp / BigR
969 chrisfen 2296 scale = rij * ri
970     else
971     if (j_is_SplitDipole) then
972 chrisfen 2755 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
973     ri = 1.0_dp / BigR
974 chrisfen 2296 scale = rij * ri
975     else
976     ri = riji
977 chrisfen 2755 scale = 1.0_dp
978 chrisfen 2296 endif
979     endif
980 chrisfen 2418
981 chrisfen 2843 if (screeningMethod .eq. DAMPED) then
982     ! assemble the damping variables
983     call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
984     c1 = erfcVal*ri
985     c2 = (-derfcVal + c1)*ri
986     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*ri
987     c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*ri*ri
988     else
989     c1 = ri
990     c2 = c1*ri
991     c3 = 3.0_dp*c2*ri
992     c4 = 5.0_dp*c3*ri*ri
993     endif
994    
995     ! precompute variables for convenience
996 chrisfen 2820 sc2 = scale * scale
997 chrisfen 2843 cti3 = ct_i*sc2*c3
998     ctj3 = ct_j*sc2*c3
999 chrisfen 2820 ctidotj = ct_i * ct_j * sc2
1000 chrisfen 2843 preSwSc = preSw*scale
1001     c2ri = c2*ri
1002     c3ri = c3*ri
1003     c4rij = c4*rij
1004 chrisfen 2755
1005 chrisfen 2843
1006 chrisfen 2820 ! calculate the potential
1007 chrisfen 2843 pot_term = (ct_ij*c2ri - ctidotj*c3)
1008     vterm = pref * pot_term
1009 chrisfen 2820 vpair = vpair + vterm
1010     epot = epot + sw*vterm
1011 chrisfen 2755
1012 chrisfen 2820 ! calculate derivatives for the forces and torques
1013 chrisfen 2843 dudx = dudx + preSwSc * ( ctidotj*xhat*c4rij - &
1014     (ct_i*uz_j(1) + ct_j*uz_i(1) + ct_ij*xhat)*c3ri )
1015     dudy = dudy + preSwSc * ( ctidotj*yhat*c4rij - &
1016     (ct_i*uz_j(2) + ct_j*uz_i(2) + ct_ij*yhat)*c3ri )
1017     dudz = dudz + preSwSc * ( ctidotj*zhat*c4rij - &
1018     (ct_i*uz_j(3) + ct_j*uz_i(3) + ct_ij*zhat)*c3ri )
1019 chrisfen 2820
1020 chrisfen 2843 duduz_i(1) = duduz_i(1) + preSw * ( uz_j(1)*c2ri - ctj3*xhat )
1021     duduz_i(2) = duduz_i(2) + preSw * ( uz_j(2)*c2ri - ctj3*yhat )
1022     duduz_i(3) = duduz_i(3) + preSw * ( uz_j(3)*c2ri - ctj3*zhat )
1023 chrisfen 2296
1024 chrisfen 2843 duduz_j(1) = duduz_j(1) + preSw * ( uz_i(1)*c2ri - cti3*xhat )
1025     duduz_j(2) = duduz_j(2) + preSw * ( uz_i(2)*c2ri - cti3*yhat )
1026     duduz_j(3) = duduz_j(3) + preSw * ( uz_i(3)*c2ri - cti3*zhat )
1027 chrisfen 2548
1028 chrisfen 2296 endif
1029 gezelter 2095 endif
1030     endif
1031 gezelter 2123
1032     if (i_is_Quadrupole) then
1033     if (j_is_Charge) then
1034 chrisfen 2972 ! precompute some necessary variables
1035     cx2 = cx_i * cx_i
1036     cy2 = cy_i * cy_i
1037     cz2 = cz_i * cz_i
1038 gezelter 3441 pref = electroMult * pre14 * q_j * one_third
1039 chrisfen 2972
1040 chrisfen 2548 if (screeningMethod .eq. DAMPED) then
1041 chrisfen 2755 ! assemble the damping variables
1042 chrisfen 2843 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
1043     c1 = erfcVal*riji
1044     c2 = (-derfcVal + c1)*riji
1045     c3 = -2.0_dp*derfcVal*alpha2 + 3.0_dp*c2*riji
1046     c4 = -4.0_dp*derfcVal*alpha4 + 5.0_dp*c3*riji*riji
1047     else
1048     c1 = riji
1049     c2 = c1*riji
1050     c3 = 3.0_dp*c2*riji
1051     c4 = 5.0_dp*c3*riji*riji
1052 chrisfen 2548 endif
1053 chrisfen 2843
1054 chrisfen 2972 ! precompute some variables for convenience
1055 chrisfen 2843 preSw = sw*pref
1056     c2ri = c2*riji
1057     c3ri = c3*riji
1058     c4rij = c4*rij
1059     xhatdot2 = 2.0_dp*xhat*c3
1060     yhatdot2 = 2.0_dp*yhat*c3
1061     zhatdot2 = 2.0_dp*zhat*c3
1062     xhatc4 = xhat*c4rij
1063     yhatc4 = yhat*c4rij
1064     zhatc4 = zhat*c4rij
1065 chrisfen 2755
1066 chrisfen 2972 ! calculate the potential
1067     pot_term = ( qxx_i * (cx2*c3 - c2ri) + qyy_i * (cy2*c3 - c2ri) + &
1068     qzz_i * (cz2*c3 - c2ri) )
1069    
1070     vterm = pref * pot_term
1071     vpair = vpair + vterm
1072     epot = epot + sw*vterm
1073    
1074 chrisfen 2843 ! calculate the derivatives for the forces and torques
1075     dudx = dudx - preSw * ( &
1076     qxx_i*(cx2*xhatc4 - (2.0_dp*cx_i*ux_i(1) + xhat)*c3ri) + &
1077     qyy_i*(cy2*xhatc4 - (2.0_dp*cy_i*uy_i(1) + xhat)*c3ri) + &
1078     qzz_i*(cz2*xhatc4 - (2.0_dp*cz_i*uz_i(1) + xhat)*c3ri) )
1079     dudy = dudy - preSw * ( &
1080     qxx_i*(cx2*yhatc4 - (2.0_dp*cx_i*ux_i(2) + yhat)*c3ri) + &
1081     qyy_i*(cy2*yhatc4 - (2.0_dp*cy_i*uy_i(2) + yhat)*c3ri) + &
1082     qzz_i*(cz2*yhatc4 - (2.0_dp*cz_i*uz_i(2) + yhat)*c3ri) )
1083     dudz = dudz - preSw * ( &
1084     qxx_i*(cx2*zhatc4 - (2.0_dp*cx_i*ux_i(3) + zhat)*c3ri) + &
1085     qyy_i*(cy2*zhatc4 - (2.0_dp*cy_i*uy_i(3) + zhat)*c3ri) + &
1086     qzz_i*(cz2*zhatc4 - (2.0_dp*cz_i*uz_i(3) + zhat)*c3ri) )
1087 chrisfen 2439
1088 chrisfen 2843 dudux_i(1) = dudux_i(1) + preSw*(qxx_i*cx_i*xhatdot2)
1089     dudux_i(2) = dudux_i(2) + preSw*(qxx_i*cx_i*yhatdot2)
1090     dudux_i(3) = dudux_i(3) + preSw*(qxx_i*cx_i*zhatdot2)
1091 chrisfen 2439
1092 chrisfen 2843 duduy_i(1) = duduy_i(1) + preSw*(qyy_i*cy_i*xhatdot2)
1093     duduy_i(2) = duduy_i(2) + preSw*(qyy_i*cy_i*yhatdot2)
1094     duduy_i(3) = duduy_i(3) + preSw*(qyy_i*cy_i*zhatdot2)
1095 chrisfen 2439
1096 chrisfen 2843 duduz_i(1) = duduz_i(1) + preSw*(qzz_i*cz_i*xhatdot2)
1097     duduz_i(2) = duduz_i(2) + preSw*(qzz_i*cz_i*yhatdot2)
1098     duduz_i(3) = duduz_i(3) + preSw*(qzz_i*cz_i*zhatdot2)
1099 gezelter 2123 endif
1100     endif
1101 gezelter 2204
1102 gezelter 3559 pot = pot + epot
1103 gezelter 2204
1104 gezelter 3559 f1(1) = f1(1) + dudx
1105     f1(2) = f1(2) + dudy
1106     f1(3) = f1(3) + dudz
1107 gezelter 2204
1108 gezelter 2095 if (i_is_Dipole .or. i_is_Quadrupole) then
1109 gezelter 3559 t1(1) = t1(1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1110     t1(2) = t1(2) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1111     t1(3) = t1(3) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1112 gezelter 2095 endif
1113 gezelter 2123 if (i_is_Quadrupole) then
1114 gezelter 3559 t1(1) = t1(1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
1115     t1(2) = t1(2) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
1116     t1(3) = t1(3) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
1117 gezelter 2095
1118 gezelter 3559 t1(1) = t1(1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
1119     t1(2) = t1(2) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
1120     t1(3) = t1(3) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
1121 gezelter 2123 endif
1122    
1123 gezelter 2095 if (j_is_Dipole .or. j_is_Quadrupole) then
1124 gezelter 3559 t2(1) = t2(1) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1125     t2(2) = t2(2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1126     t2(3) = t2(3) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1127 gezelter 2095 endif
1128 gezelter 2123 if (j_is_Quadrupole) then
1129 gezelter 3559 t2(1) = t2(1) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
1130     t2(2) = t2(2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
1131     t2(3) = t2(3) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
1132 gezelter 2095
1133 gezelter 3559 t2(1) = t2(1) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
1134     t2(2) = t2(2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
1135     t2(3) = t2(3) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
1136 gezelter 2123 endif
1137    
1138 gezelter 2095 return
1139     end subroutine doElectrostaticPair
1140 chuckv 2189
1141     subroutine destroyElectrostaticTypes()
1142    
1143 gezelter 2204 if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap)
1144    
1145 chuckv 2189 end subroutine destroyElectrostaticTypes
1146 gezelter 3500
1147     subroutine self_self(atom1, eFrame, skch, mypot, t, do_pot)
1148 chrisfen 2394 logical, intent(in) :: do_pot
1149 chrisfen 2381 integer, intent(in) :: atom1
1150 chrisfen 2394 integer :: atid1
1151 chrisfen 2381 real(kind=dp), dimension(9,nLocal) :: eFrame
1152 chrisfen 2394 real(kind=dp), dimension(3,nLocal) :: t
1153 gezelter 3500 real(kind=dp) :: mu1, chg1, c1e
1154     real(kind=dp) :: preVal, epot, mypot, skch
1155     real(kind=dp) :: eix, eiy, eiz, self
1156 chrisfen 2381
1157 chrisfen 2394 ! this is a local only array, so we use the local atom type id's:
1158     atid1 = atid(atom1)
1159 chrisfen 2402
1160     if (.not.summationMethodChecked) then
1161     call checkSummationMethod()
1162     endif
1163 chrisfen 2394
1164 chrisfen 2402 if (summationMethod .eq. REACTION_FIELD) then
1165     if (ElectrostaticMap(atid1)%is_Dipole) then
1166     mu1 = getDipoleMoment(atid1)
1167    
1168     preVal = pre22 * preRF2 * mu1*mu1
1169 chrisfen 2755 mypot = mypot - 0.5_dp*preVal
1170 chrisfen 2402
1171     ! The self-correction term adds into the reaction field vector
1172    
1173     eix = preVal * eFrame(3,atom1)
1174     eiy = preVal * eFrame(6,atom1)
1175     eiz = preVal * eFrame(9,atom1)
1176    
1177     ! once again, this is self-self, so only the local arrays are needed
1178     ! even for MPI jobs:
1179    
1180     t(1,atom1)=t(1,atom1) - eFrame(6,atom1)*eiz + &
1181     eFrame(9,atom1)*eiy
1182     t(2,atom1)=t(2,atom1) - eFrame(9,atom1)*eix + &
1183     eFrame(3,atom1)*eiz
1184     t(3,atom1)=t(3,atom1) - eFrame(3,atom1)*eiy + &
1185     eFrame(6,atom1)*eix
1186    
1187     endif
1188    
1189 chrisfen 2442 elseif ( (summationMethod .eq. SHIFTED_FORCE) .or. &
1190     (summationMethod .eq. SHIFTED_POTENTIAL) ) then
1191 chrisfen 2416 if (ElectrostaticMap(atid1)%is_Charge) then
1192 gezelter 3500 chg1 = getCharge(atid1)
1193 chrisfen 2416 if (screeningMethod .eq. DAMPED) then
1194 gezelter 3500 self = - 0.5_dp * (c1c+alphaPi) * chg1 * (chg1+skch) * pre11
1195 chrisfen 2416 else
1196 gezelter 3500 self = - 0.5_dp * rcuti * chg1 * (chg1+skch) * pre11
1197 chrisfen 2416 endif
1198 gezelter 3500
1199     mypot = mypot + self
1200 chrisfen 2416 endif
1201     endif
1202 chrisfen 2394
1203 gezelter 3497
1204    
1205 chrisfen 2381 return
1206 chrisfen 2402 end subroutine self_self
1207 chrisfen 2381
1208 gezelter 3441 subroutine rf_self_excludes(atom1, atom2, sw, electroMult, eFrame, d, &
1209     rij, vpair, myPot, f, t, do_pot)
1210 chrisfen 2399 logical, intent(in) :: do_pot
1211     integer, intent(in) :: atom1
1212     integer, intent(in) :: atom2
1213     logical :: i_is_Charge, j_is_Charge
1214     logical :: i_is_Dipole, j_is_Dipole
1215     integer :: atid1
1216     integer :: atid2
1217     real(kind=dp), intent(in) :: rij
1218 gezelter 3441 real(kind=dp), intent(in) :: sw, electroMult
1219 chrisfen 2399 real(kind=dp), intent(in), dimension(3) :: d
1220     real(kind=dp), intent(inout) :: vpair
1221     real(kind=dp), dimension(9,nLocal) :: eFrame
1222     real(kind=dp), dimension(3,nLocal) :: f
1223     real(kind=dp), dimension(3,nLocal) :: t
1224     real (kind = dp), dimension(3) :: duduz_i
1225     real (kind = dp), dimension(3) :: duduz_j
1226     real (kind = dp), dimension(3) :: uz_i
1227     real (kind = dp), dimension(3) :: uz_j
1228     real(kind=dp) :: q_i, q_j, mu_i, mu_j
1229     real(kind=dp) :: xhat, yhat, zhat
1230     real(kind=dp) :: ct_i, ct_j
1231     real(kind=dp) :: ri2, ri3, riji, vterm
1232 chrisfen 2402 real(kind=dp) :: pref, preVal, rfVal, myPot
1233 chrisfen 2399 real(kind=dp) :: dudx, dudy, dudz, dudr
1234    
1235 chrisfen 2402 if (.not.summationMethodChecked) then
1236     call checkSummationMethod()
1237 chrisfen 2399 endif
1238    
1239 gezelter 2722 dudx = zero
1240     dudy = zero
1241     dudz = zero
1242 chrisfen 2399
1243 chrisfen 2755 riji = 1.0_dp/rij
1244 chrisfen 2399
1245     xhat = d(1) * riji
1246     yhat = d(2) * riji
1247     zhat = d(3) * riji
1248    
1249     ! this is a local only array, so we use the local atom type id's:
1250     atid1 = atid(atom1)
1251     atid2 = atid(atom2)
1252     i_is_Charge = ElectrostaticMap(atid1)%is_Charge
1253     j_is_Charge = ElectrostaticMap(atid2)%is_Charge
1254     i_is_Dipole = ElectrostaticMap(atid1)%is_Dipole
1255     j_is_Dipole = ElectrostaticMap(atid2)%is_Dipole
1256    
1257     if (i_is_Charge.and.j_is_Charge) then
1258     q_i = ElectrostaticMap(atid1)%charge
1259     q_j = ElectrostaticMap(atid2)%charge
1260    
1261 gezelter 3441 preVal = electroMult * pre11 * q_i * q_j
1262 chrisfen 2399 rfVal = preRF*rij*rij
1263     vterm = preVal * rfVal
1264    
1265 chrisfen 2402 myPot = myPot + sw*vterm
1266    
1267 chrisfen 2755 dudr = sw*preVal * 2.0_dp*rfVal*riji
1268 chrisfen 2402
1269 chrisfen 2399 dudx = dudx + dudr * xhat
1270     dudy = dudy + dudr * yhat
1271     dudz = dudz + dudr * zhat
1272 chrisfen 2402
1273 chrisfen 2399 elseif (i_is_Charge.and.j_is_Dipole) then
1274     q_i = ElectrostaticMap(atid1)%charge
1275     mu_j = ElectrostaticMap(atid2)%dipole_moment
1276     uz_j(1) = eFrame(3,atom2)
1277     uz_j(2) = eFrame(6,atom2)
1278     uz_j(3) = eFrame(9,atom2)
1279     ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
1280 chrisfen 2402
1281 chrisfen 2399 ri2 = riji * riji
1282     ri3 = ri2 * riji
1283    
1284 gezelter 3441 pref = electroMult * pre12 * q_i * mu_j
1285 chrisfen 2399 vterm = - pref * ct_j * ( ri2 - preRF2*rij )
1286 chrisfen 2402 myPot = myPot + sw*vterm
1287    
1288 chrisfen 2755 dudx = dudx - sw*pref*( ri3*(uz_j(1)-3.0_dp*ct_j*xhat) &
1289 chrisfen 2402 - preRF2*uz_j(1) )
1290 chrisfen 2755 dudy = dudy - sw*pref*( ri3*(uz_j(2)-3.0_dp*ct_j*yhat) &
1291 chrisfen 2402 - preRF2*uz_j(2) )
1292 chrisfen 2755 dudz = dudz - sw*pref*( ri3*(uz_j(3)-3.0_dp*ct_j*zhat) &
1293 chrisfen 2402 - preRF2*uz_j(3) )
1294    
1295 chrisfen 2399 duduz_j(1) = duduz_j(1) - sw * pref * xhat * ( ri2 - preRF2*rij )
1296     duduz_j(2) = duduz_j(2) - sw * pref * yhat * ( ri2 - preRF2*rij )
1297     duduz_j(3) = duduz_j(3) - sw * pref * zhat * ( ri2 - preRF2*rij )
1298 chrisfen 2402
1299 chrisfen 2399 elseif (i_is_Dipole.and.j_is_Charge) then
1300     mu_i = ElectrostaticMap(atid1)%dipole_moment
1301     q_j = ElectrostaticMap(atid2)%charge
1302     uz_i(1) = eFrame(3,atom1)
1303     uz_i(2) = eFrame(6,atom1)
1304     uz_i(3) = eFrame(9,atom1)
1305     ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
1306 chrisfen 2402
1307 chrisfen 2399 ri2 = riji * riji
1308     ri3 = ri2 * riji
1309    
1310 gezelter 3441 pref = electroMult * pre12 * q_j * mu_i
1311 chrisfen 2399 vterm = pref * ct_i * ( ri2 - preRF2*rij )
1312 chrisfen 2402 myPot = myPot + sw*vterm
1313 chrisfen 2399
1314 chrisfen 2755 dudx = dudx + sw*pref*( ri3*(uz_i(1)-3.0_dp*ct_i*xhat) &
1315 chrisfen 2402 - preRF2*uz_i(1) )
1316 chrisfen 2755 dudy = dudy + sw*pref*( ri3*(uz_i(2)-3.0_dp*ct_i*yhat) &
1317 chrisfen 2402 - preRF2*uz_i(2) )
1318 chrisfen 2755 dudz = dudz + sw*pref*( ri3*(uz_i(3)-3.0_dp*ct_i*zhat) &
1319 chrisfen 2402 - preRF2*uz_i(3) )
1320 chrisfen 2399
1321     duduz_i(1) = duduz_i(1) + sw * pref * xhat * ( ri2 - preRF2*rij )
1322     duduz_i(2) = duduz_i(2) + sw * pref * yhat * ( ri2 - preRF2*rij )
1323     duduz_i(3) = duduz_i(3) + sw * pref * zhat * ( ri2 - preRF2*rij )
1324    
1325     endif
1326 chrisfen 2402
1327    
1328     ! accumulate the forces and torques resulting from the self term
1329 chrisfen 2399 f(1,atom1) = f(1,atom1) + dudx
1330     f(2,atom1) = f(2,atom1) + dudy
1331     f(3,atom1) = f(3,atom1) + dudz
1332    
1333     f(1,atom2) = f(1,atom2) - dudx
1334     f(2,atom2) = f(2,atom2) - dudy
1335     f(3,atom2) = f(3,atom2) - dudz
1336    
1337     if (i_is_Dipole) then
1338     t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1339     t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1340     t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1341     elseif (j_is_Dipole) then
1342     t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1343     t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1344     t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1345     endif
1346    
1347     return
1348     end subroutine rf_self_excludes
1349    
1350 chrisfen 2917 subroutine accumulate_box_dipole(atom1, eFrame, d, pChg, nChg, pChgPos, &
1351     nChgPos, dipVec, pChgCount, nChgCount)
1352     integer, intent(in) :: atom1
1353     logical :: i_is_Charge
1354     logical :: i_is_Dipole
1355     integer :: atid1
1356     integer :: pChgCount
1357     integer :: nChgCount
1358     real(kind=dp), intent(in), dimension(3) :: d
1359     real(kind=dp), dimension(9,nLocal) :: eFrame
1360     real(kind=dp) :: pChg
1361     real(kind=dp) :: nChg
1362     real(kind=dp), dimension(3) :: pChgPos
1363     real(kind=dp), dimension(3) :: nChgPos
1364     real(kind=dp), dimension(3) :: dipVec
1365     real(kind=dp), dimension(3) :: uz_i
1366     real(kind=dp), dimension(3) :: pos
1367     real(kind=dp) :: q_i, mu_i
1368     real(kind=dp) :: pref, preVal
1369    
1370     if (.not.summationMethodChecked) then
1371     call checkSummationMethod()
1372     endif
1373    
1374     ! this is a local only array, so we use the local atom type id's:
1375     atid1 = atid(atom1)
1376     i_is_Charge = ElectrostaticMap(atid1)%is_Charge
1377     i_is_Dipole = ElectrostaticMap(atid1)%is_Dipole
1378    
1379     if (i_is_Charge) then
1380     q_i = ElectrostaticMap(atid1)%charge
1381     ! convert to the proper units
1382     q_i = q_i * chargeToC
1383     pos = d * angstromToM
1384    
1385     if (q_i.le.0.0_dp) then
1386     nChg = nChg - q_i
1387     nChgPos(1) = nChgPos(1) + pos(1)
1388     nChgPos(2) = nChgPos(2) + pos(2)
1389     nChgPos(3) = nChgPos(3) + pos(3)
1390     nChgCount = nChgCount + 1
1391    
1392     else
1393     pChg = pChg + q_i
1394     pChgPos(1) = pChgPos(1) + pos(1)
1395     pChgPos(2) = pChgPos(2) + pos(2)
1396     pChgPos(3) = pChgPos(3) + pos(3)
1397     pChgCount = pChgCount + 1
1398    
1399     endif
1400    
1401     endif
1402    
1403     if (i_is_Dipole) then
1404     mu_i = ElectrostaticMap(atid1)%dipole_moment
1405     uz_i(1) = eFrame(3,atom1)
1406     uz_i(2) = eFrame(6,atom1)
1407     uz_i(3) = eFrame(9,atom1)
1408     ! convert to the proper units
1409     mu_i = mu_i * debyeToCm
1410    
1411     dipVec(1) = dipVec(1) + uz_i(1)*mu_i
1412     dipVec(2) = dipVec(2) + uz_i(2)*mu_i
1413     dipVec(3) = dipVec(3) + uz_i(3)*mu_i
1414    
1415     endif
1416    
1417     return
1418     end subroutine accumulate_box_dipole
1419    
1420 gezelter 2095 end module electrostatic_module