ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2843
Committed: Fri Jun 9 18:26:18 2006 UTC (18 years, 1 month ago) by chrisfen
File size: 50686 byte(s)
Log Message:
reformulated some of the electrostatics

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