ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 3122
Committed: Wed Feb 28 00:53:14 2007 UTC (17 years, 4 months ago) by chuckv
File size: 53007 byte(s)
Log Message:
Fixed conversion bug in self_self, when adding self energy to mypot, energy needs to be converted by pre11.

File Contents

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