ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2755
Committed: Wed May 17 14:03:15 2006 UTC (18 years, 3 months ago) by chrisfen
File size: 50347 byte(s)
Log Message:
electrostatic refinement

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