ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2301
Committed: Thu Sep 15 22:05:21 2005 UTC (18 years, 10 months ago) by gezelter
File size: 44239 byte(s)
Log Message:
Working on adding WOLF

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 #ifdef IS_MPI
51 use mpiSimulation
52 #endif
53 implicit none
54
55 PRIVATE
56
57 #define __FORTRAN90
58 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59
60 !! these prefactors convert the multipole interactions into kcal / mol
61 !! all were computed assuming distances are measured in angstroms
62 !! Charge-Charge, assuming charges are measured in electrons
63 real(kind=dp), parameter :: pre11 = 332.0637778_dp
64 !! Charge-Dipole, assuming charges are measured in electrons, and
65 !! dipoles are measured in debyes
66 real(kind=dp), parameter :: pre12 = 69.13373_dp
67 !! Dipole-Dipole, assuming dipoles are measured in debyes
68 real(kind=dp), parameter :: pre22 = 14.39325_dp
69 !! Charge-Quadrupole, assuming charges are measured in electrons, and
70 !! quadrupoles are measured in 10^-26 esu cm^2
71 !! This unit is also known affectionately as an esu centi-barn.
72 real(kind=dp), parameter :: pre14 = 69.13373_dp
73
74 !! variables to handle different summation methods for long-range electrostatics:
75 integer, save :: summationMethod = NONE
76 real(kind=DP), save :: defaultCutoff = 0.0_DP
77 logical, save :: haveDefaultCutoff = .false.
78 real(kind=DP), save :: dampingAlpha = 0.0_DP
79 logical, save :: haveDampingAlpha = .false.
80 real(kind=DP), save :: dielectric = 0.0_DP
81 logical, save :: haveDielectric = .false.
82 real(kind=DP), save :: constERFC = 0.0_DP
83 real(kind=DP), save :: constEXP = 0.0_DP
84 logical, save :: haveDWAconstants = .false.
85
86
87 public :: setElectrostaticSummationMethod
88 public :: setElectrostaticCutoffRadius
89 public :: setDampedWolfAlpha
90 public :: setReactionFieldDielectric
91 public :: newElectrostaticType
92 public :: setCharge
93 public :: setDipoleMoment
94 public :: setSplitDipoleDistance
95 public :: setQuadrupoleMoments
96 public :: doElectrostaticPair
97 public :: getCharge
98 public :: getDipoleMoment
99 public :: pre22
100 public :: destroyElectrostaticTypes
101
102 type :: Electrostatic
103 integer :: c_ident
104 logical :: is_Charge = .false.
105 logical :: is_Dipole = .false.
106 logical :: is_SplitDipole = .false.
107 logical :: is_Quadrupole = .false.
108 logical :: is_Tap = .false.
109 real(kind=DP) :: charge = 0.0_DP
110 real(kind=DP) :: dipole_moment = 0.0_DP
111 real(kind=DP) :: split_dipole_distance = 0.0_DP
112 real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP
113 end type Electrostatic
114
115 type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
116
117 contains
118
119 subroutine setElectrostaticSummationMethod(the_ESM)
120
121 integer, intent(in) :: the_ESM
122
123 if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then
124 call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method")
125 endif
126
127 end subroutine setElectrostaticSummationMethod
128
129 subroutine setElectrostaticCutoffRadius(thisRcut)
130 real(kind=dp), intent(in) :: thisRcut
131 defaultCutoff = thisRcut
132 haveDefaultCutoff = .true.
133 end subroutine setElectrostaticCutoffRadius
134
135 subroutine setDampedWolfAlpha(thisAlpha)
136 real(kind=dp), intent(in) :: thisAlpha
137 dampingAlpha = thisAlpha
138 haveDampingAlpha = .true.
139 end subroutine setDampedWolfAlpha
140
141 subroutine setReactionFieldDielectric(thisDielectric)
142 real(kind=dp), intent(in) :: thisDielectric
143 dielectric = thisDielectric
144 haveDielectric = .true.
145 end subroutine setReactionFieldDielectric
146
147 subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
148 is_SplitDipole, is_Quadrupole, is_Tap, status)
149
150 integer, intent(in) :: c_ident
151 logical, intent(in) :: is_Charge
152 logical, intent(in) :: is_Dipole
153 logical, intent(in) :: is_SplitDipole
154 logical, intent(in) :: is_Quadrupole
155 logical, intent(in) :: is_Tap
156 integer, intent(out) :: status
157 integer :: nAtypes, myATID, i, j
158
159 status = 0
160 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
161
162 !! Be simple-minded and assume that we need an ElectrostaticMap that
163 !! is the same size as the total number of atom types
164
165 if (.not.allocated(ElectrostaticMap)) then
166
167 nAtypes = getSize(atypes)
168
169 if (nAtypes == 0) then
170 status = -1
171 return
172 end if
173
174 if (.not. allocated(ElectrostaticMap)) then
175 allocate(ElectrostaticMap(nAtypes))
176 endif
177
178 end if
179
180 if (myATID .gt. size(ElectrostaticMap)) then
181 status = -1
182 return
183 endif
184
185 ! set the values for ElectrostaticMap for this atom type:
186
187 ElectrostaticMap(myATID)%c_ident = c_ident
188 ElectrostaticMap(myATID)%is_Charge = is_Charge
189 ElectrostaticMap(myATID)%is_Dipole = is_Dipole
190 ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole
191 ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole
192 ElectrostaticMap(myATID)%is_Tap = is_Tap
193
194 end subroutine newElectrostaticType
195
196 subroutine setCharge(c_ident, charge, status)
197 integer, intent(in) :: c_ident
198 real(kind=dp), intent(in) :: charge
199 integer, intent(out) :: status
200 integer :: myATID
201
202 status = 0
203 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
204
205 if (.not.allocated(ElectrostaticMap)) then
206 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!")
207 status = -1
208 return
209 end if
210
211 if (myATID .gt. size(ElectrostaticMap)) then
212 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!")
213 status = -1
214 return
215 endif
216
217 if (.not.ElectrostaticMap(myATID)%is_Charge) then
218 call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!")
219 status = -1
220 return
221 endif
222
223 ElectrostaticMap(myATID)%charge = charge
224 end subroutine setCharge
225
226 subroutine setDipoleMoment(c_ident, dipole_moment, status)
227 integer, intent(in) :: c_ident
228 real(kind=dp), intent(in) :: dipole_moment
229 integer, intent(out) :: status
230 integer :: myATID
231
232 status = 0
233 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
234
235 if (.not.allocated(ElectrostaticMap)) then
236 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!")
237 status = -1
238 return
239 end if
240
241 if (myATID .gt. size(ElectrostaticMap)) then
242 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!")
243 status = -1
244 return
245 endif
246
247 if (.not.ElectrostaticMap(myATID)%is_Dipole) then
248 call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!")
249 status = -1
250 return
251 endif
252
253 ElectrostaticMap(myATID)%dipole_moment = dipole_moment
254 end subroutine setDipoleMoment
255
256 subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status)
257 integer, intent(in) :: c_ident
258 real(kind=dp), intent(in) :: split_dipole_distance
259 integer, intent(out) :: status
260 integer :: myATID
261
262 status = 0
263 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
264
265 if (.not.allocated(ElectrostaticMap)) then
266 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!")
267 status = -1
268 return
269 end if
270
271 if (myATID .gt. size(ElectrostaticMap)) then
272 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!")
273 status = -1
274 return
275 endif
276
277 if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then
278 call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!")
279 status = -1
280 return
281 endif
282
283 ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance
284 end subroutine setSplitDipoleDistance
285
286 subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status)
287 integer, intent(in) :: c_ident
288 real(kind=dp), intent(in), dimension(3) :: quadrupole_moments
289 integer, intent(out) :: status
290 integer :: myATID, i, j
291
292 status = 0
293 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
294
295 if (.not.allocated(ElectrostaticMap)) then
296 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!")
297 status = -1
298 return
299 end if
300
301 if (myATID .gt. size(ElectrostaticMap)) then
302 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!")
303 status = -1
304 return
305 endif
306
307 if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then
308 call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!")
309 status = -1
310 return
311 endif
312
313 do i = 1, 3
314 ElectrostaticMap(myATID)%quadrupole_moments(i) = &
315 quadrupole_moments(i)
316 enddo
317
318 end subroutine setQuadrupoleMoments
319
320
321 function getCharge(atid) result (c)
322 integer, intent(in) :: atid
323 integer :: localError
324 real(kind=dp) :: c
325
326 if (.not.allocated(ElectrostaticMap)) then
327 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
328 return
329 end if
330
331 if (.not.ElectrostaticMap(atid)%is_Charge) then
332 call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!")
333 return
334 endif
335
336 c = ElectrostaticMap(atid)%charge
337 end function getCharge
338
339 function getDipoleMoment(atid) result (dm)
340 integer, intent(in) :: atid
341 integer :: localError
342 real(kind=dp) :: dm
343
344 if (.not.allocated(ElectrostaticMap)) then
345 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
346 return
347 end if
348
349 if (.not.ElectrostaticMap(atid)%is_Dipole) then
350 call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!")
351 return
352 endif
353
354 dm = ElectrostaticMap(atid)%dipole_moment
355 end function getDipoleMoment
356
357 subroutine checkSummationMethod()
358
359 if (summationMethod .eq. DAMPED_WOLF) then
360 if (.not.haveDWAconstants) then
361
362 if (.not.haveDampingAlpha) then
363 call handleError("checkSummationMethod", "no Damping Alpha set!")
364 endif
365
366 if (.not.have....)
367 constEXP =
368 constERFC =
369
370 haveDWAconstants = .true.
371 endif
372 endif
373
374 end subroutine checkSummationMethod
375
376
377
378 subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, &
379 vpair, fpair, pot, eFrame, f, t, do_pot, corrMethod, rcuti)
380
381 logical, intent(in) :: do_pot
382
383 integer, intent(in) :: atom1, atom2
384 integer :: localError
385 integer, intent(in) :: corrMethod
386
387 real(kind=dp), intent(in) :: rij, r2, sw, rcuti
388 real(kind=dp), intent(in), dimension(3) :: d
389 real(kind=dp), intent(inout) :: vpair
390 real(kind=dp), intent(inout), dimension(3) :: fpair
391
392 real( kind = dp ) :: pot, swi
393 real( kind = dp ), dimension(9,nLocal) :: eFrame
394 real( kind = dp ), dimension(3,nLocal) :: f
395 real( kind = dp ), dimension(3,nLocal) :: t
396
397 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
398 real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
399 real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
400 real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
401
402 logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
403 logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
404 logical :: i_is_Tap, j_is_Tap
405 integer :: me1, me2, id1, id2
406 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
407 real (kind=dp) :: qxx_i, qyy_i, qzz_i
408 real (kind=dp) :: qxx_j, qyy_j, qzz_j
409 real (kind=dp) :: cx_i, cy_i, cz_i
410 real (kind=dp) :: cx_j, cy_j, cz_j
411 real (kind=dp) :: cx2, cy2, cz2
412 real (kind=dp) :: ct_i, ct_j, ct_ij, a1
413 real (kind=dp) :: riji, ri, ri2, ri3, ri4
414 real (kind=dp) :: pref, vterm, epot, dudr, vterm1, vterm2
415 real (kind=dp) :: xhat, yhat, zhat
416 real (kind=dp) :: dudx, dudy, dudz
417 real (kind=dp) :: scale, sc2, bigR, switcher, dswitcher
418 real (kind=dp) :: rcuti2, rcuti3, rcuti4
419
420 if (.not.allocated(ElectrostaticMap)) then
421 call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
422 return
423 end if
424
425 if (.not.summationMethodChecked) then
426 call checkSummationMethod()
427 endif
428
429
430 #ifdef IS_MPI
431 me1 = atid_Row(atom1)
432 me2 = atid_Col(atom2)
433 #else
434 me1 = atid(atom1)
435 me2 = atid(atom2)
436 #endif
437
438 !! some variables we'll need independent of electrostatic type:
439
440 riji = 1.0d0 / rij
441
442 xhat = d(1) * riji
443 yhat = d(2) * riji
444 zhat = d(3) * riji
445
446 rcuti2 = rcuti*rcuti
447 rcuti3 = rcuti2*rcuti
448 rcuti4 = rcuti2*rcuti2
449
450 swi = 1.0d0 / sw
451
452 !! logicals
453 i_is_Charge = ElectrostaticMap(me1)%is_Charge
454 i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
455 i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
456 i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
457 i_is_Tap = ElectrostaticMap(me1)%is_Tap
458
459 j_is_Charge = ElectrostaticMap(me2)%is_Charge
460 j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
461 j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
462 j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
463 j_is_Tap = ElectrostaticMap(me2)%is_Tap
464
465 if (i_is_Charge) then
466 q_i = ElectrostaticMap(me1)%charge
467 endif
468
469 if (i_is_Dipole) then
470 mu_i = ElectrostaticMap(me1)%dipole_moment
471 #ifdef IS_MPI
472 uz_i(1) = eFrame_Row(3,atom1)
473 uz_i(2) = eFrame_Row(6,atom1)
474 uz_i(3) = eFrame_Row(9,atom1)
475 #else
476 uz_i(1) = eFrame(3,atom1)
477 uz_i(2) = eFrame(6,atom1)
478 uz_i(3) = eFrame(9,atom1)
479 #endif
480 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
481
482 if (i_is_SplitDipole) then
483 d_i = ElectrostaticMap(me1)%split_dipole_distance
484 endif
485
486 endif
487
488 if (i_is_Quadrupole) then
489 qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
490 qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
491 qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
492 #ifdef IS_MPI
493 ux_i(1) = eFrame_Row(1,atom1)
494 ux_i(2) = eFrame_Row(4,atom1)
495 ux_i(3) = eFrame_Row(7,atom1)
496 uy_i(1) = eFrame_Row(2,atom1)
497 uy_i(2) = eFrame_Row(5,atom1)
498 uy_i(3) = eFrame_Row(8,atom1)
499 uz_i(1) = eFrame_Row(3,atom1)
500 uz_i(2) = eFrame_Row(6,atom1)
501 uz_i(3) = eFrame_Row(9,atom1)
502 #else
503 ux_i(1) = eFrame(1,atom1)
504 ux_i(2) = eFrame(4,atom1)
505 ux_i(3) = eFrame(7,atom1)
506 uy_i(1) = eFrame(2,atom1)
507 uy_i(2) = eFrame(5,atom1)
508 uy_i(3) = eFrame(8,atom1)
509 uz_i(1) = eFrame(3,atom1)
510 uz_i(2) = eFrame(6,atom1)
511 uz_i(3) = eFrame(9,atom1)
512 #endif
513 cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
514 cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
515 cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
516 endif
517
518 if (j_is_Charge) then
519 q_j = ElectrostaticMap(me2)%charge
520 endif
521
522 if (j_is_Dipole) then
523 mu_j = ElectrostaticMap(me2)%dipole_moment
524 #ifdef IS_MPI
525 uz_j(1) = eFrame_Col(3,atom2)
526 uz_j(2) = eFrame_Col(6,atom2)
527 uz_j(3) = eFrame_Col(9,atom2)
528 #else
529 uz_j(1) = eFrame(3,atom2)
530 uz_j(2) = eFrame(6,atom2)
531 uz_j(3) = eFrame(9,atom2)
532 #endif
533 ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
534
535 if (j_is_SplitDipole) then
536 d_j = ElectrostaticMap(me2)%split_dipole_distance
537 endif
538 endif
539
540 if (j_is_Quadrupole) then
541 qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
542 qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
543 qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
544 #ifdef IS_MPI
545 ux_j(1) = eFrame_Col(1,atom2)
546 ux_j(2) = eFrame_Col(4,atom2)
547 ux_j(3) = eFrame_Col(7,atom2)
548 uy_j(1) = eFrame_Col(2,atom2)
549 uy_j(2) = eFrame_Col(5,atom2)
550 uy_j(3) = eFrame_Col(8,atom2)
551 uz_j(1) = eFrame_Col(3,atom2)
552 uz_j(2) = eFrame_Col(6,atom2)
553 uz_j(3) = eFrame_Col(9,atom2)
554 #else
555 ux_j(1) = eFrame(1,atom2)
556 ux_j(2) = eFrame(4,atom2)
557 ux_j(3) = eFrame(7,atom2)
558 uy_j(1) = eFrame(2,atom2)
559 uy_j(2) = eFrame(5,atom2)
560 uy_j(3) = eFrame(8,atom2)
561 uz_j(1) = eFrame(3,atom2)
562 uz_j(2) = eFrame(6,atom2)
563 uz_j(3) = eFrame(9,atom2)
564 #endif
565 cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
566 cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
567 cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
568 endif
569
570 !!$ switcher = 1.0d0
571 !!$ dswitcher = 0.0d0
572 !!$ ebalance = 0.0d0
573 !!$ ! weaken the dipole interaction at close range for TAP water
574 !!$ if (j_is_Tap .and. i_is_Tap) then
575 !!$ call calc_switch(rij, mu_i, switcher, dswitcher)
576 !!$ endif
577
578 epot = 0.0_dp
579 dudx = 0.0_dp
580 dudy = 0.0_dp
581 dudz = 0.0_dp
582
583 dudux_i = 0.0_dp
584 duduy_i = 0.0_dp
585 duduz_i = 0.0_dp
586
587 dudux_j = 0.0_dp
588 duduy_j = 0.0_dp
589 duduz_j = 0.0_dp
590
591 if (i_is_Charge) then
592
593 if (j_is_Charge) then
594
595 if (corrMethod .eq. 1) then
596 vterm = pre11 * q_i * q_j * (riji - rcuti)
597
598 vpair = vpair + vterm
599 epot = epot + sw * vterm
600
601 dudr = - sw * pre11 * q_i * q_j * (riji*riji*riji - rcuti2*rcuti)
602
603 dudx = dudx + dudr * d(1)
604 dudy = dudy + dudr * d(2)
605 dudz = dudz + dudr * d(3)
606
607 else
608 vterm = pre11 * q_i * q_j * riji
609
610 vpair = vpair + vterm
611 epot = epot + sw * vterm
612
613 dudr = - sw * vterm * riji
614
615 dudx = dudx + dudr * xhat
616 dudy = dudy + dudr * yhat
617 dudz = dudz + dudr * zhat
618
619 endif
620
621 endif
622
623 if (j_is_Dipole) then
624
625 pref = sw * pre12 * q_i * mu_j
626
627 if (corrMethod .eq. 1) then
628 ri2 = riji * riji
629 ri3 = ri2 * riji
630
631 vterm = - pref * ct_j * (ri2 - rcuti2)
632 vpair = vpair + swi*vterm
633 epot = epot + vterm
634
635 !! this has a + sign in the () because the rij vector is
636 !! r_j - r_i and the charge-dipole potential takes the origin
637 !! as the point dipole, which is atom j in this case.
638
639 dudx = dudx - pref * ( ri3*( uz_j(1) - 3.0d0*ct_j*xhat) &
640 - rcuti3*( uz_j(1) - 3.0d0*ct_j*d(1)*rcuti ) )
641 dudy = dudy - pref * ( ri3*( uz_j(2) - 3.0d0*ct_j*yhat) &
642 - rcuti3*( uz_j(2) - 3.0d0*ct_j*d(2)*rcuti ) )
643 dudz = dudz - pref * ( ri3*( uz_j(3) - 3.0d0*ct_j*zhat) &
644 - rcuti3*( uz_j(3) - 3.0d0*ct_j*d(3)*rcuti ) )
645
646 duduz_j(1) = duduz_j(1) - pref*( ri2*xhat - d(1)*rcuti3 )
647 duduz_j(2) = duduz_j(2) - pref*( ri2*yhat - d(2)*rcuti3 )
648 duduz_j(3) = duduz_j(3) - pref*( ri2*zhat - d(3)*rcuti3 )
649
650 else
651 if (j_is_SplitDipole) then
652 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
653 ri = 1.0_dp / BigR
654 scale = rij * ri
655 else
656 ri = riji
657 scale = 1.0_dp
658 endif
659
660 ri2 = ri * ri
661 ri3 = ri2 * ri
662 sc2 = scale * scale
663
664 vterm = - pref * ct_j * ri2 * scale
665 vpair = vpair + swi * vterm
666 epot = epot + vterm
667
668 !! this has a + sign in the () because the rij vector is
669 !! r_j - r_i and the charge-dipole potential takes the origin
670 !! as the point dipole, which is atom j in this case.
671
672 dudx = dudx - pref * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2)
673 dudy = dudy - pref * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2)
674 dudz = dudz - pref * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2)
675
676 duduz_j(1) = duduz_j(1) - pref * ri2 * xhat * scale
677 duduz_j(2) = duduz_j(2) - pref * ri2 * yhat * scale
678 duduz_j(3) = duduz_j(3) - pref * ri2 * zhat * scale
679
680 endif
681 endif
682
683 if (j_is_Quadrupole) then
684 ri2 = riji * riji
685 ri3 = ri2 * riji
686 ri4 = ri2 * ri2
687 cx2 = cx_j * cx_j
688 cy2 = cy_j * cy_j
689 cz2 = cz_j * cz_j
690
691
692 pref = sw * pre14 * q_i / 3.0_dp
693
694 if (corrMethod .eq. 1) then
695 vterm1 = pref * ri3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
696 qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
697 qzz_j * (3.0_dp*cz2 - 1.0_dp) )
698 vterm2 = pref * rcuti3*( qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
699 qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
700 qzz_j * (3.0_dp*cz2 - 1.0_dp) )
701 vpair = vpair + swi*( vterm1 - vterm2 )
702 epot = epot + ( vterm1 - vterm2 )
703
704 dudx = dudx - (5.0_dp * &
705 (vterm1*riji*xhat - vterm2*rcuti2*d(1))) + pref * ( &
706 (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(1)) - &
707 qxx_j*2.0_dp*(xhat - rcuti*d(1))) + &
708 (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(1)) - &
709 qyy_j*2.0_dp*(xhat - rcuti*d(1))) + &
710 (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(1)) - &
711 qzz_j*2.0_dp*(xhat - rcuti*d(1))) )
712 dudy = dudy - (5.0_dp * &
713 (vterm1*riji*yhat - vterm2*rcuti2*d(2))) + pref * ( &
714 (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(2)) - &
715 qxx_j*2.0_dp*(yhat - rcuti*d(2))) + &
716 (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(2)) - &
717 qyy_j*2.0_dp*(yhat - rcuti*d(2))) + &
718 (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(2)) - &
719 qzz_j*2.0_dp*(yhat - rcuti*d(2))) )
720 dudz = dudz - (5.0_dp * &
721 (vterm1*riji*zhat - vterm2*rcuti2*d(3))) + pref * ( &
722 (ri4 - rcuti4)*(qxx_j*(6.0_dp*cx_j*ux_j(3)) - &
723 qxx_j*2.0_dp*(zhat - rcuti*d(3))) + &
724 (ri4 - rcuti4)*(qyy_j*(6.0_dp*cy_j*uy_j(3)) - &
725 qyy_j*2.0_dp*(zhat - rcuti*d(3))) + &
726 (ri4 - rcuti4)*(qzz_j*(6.0_dp*cz_j*uz_j(3)) - &
727 qzz_j*2.0_dp*(zhat - rcuti*d(3))) )
728
729 dudux_j(1) = dudux_j(1) + pref * (ri3*(qxx_j*6.0_dp*cx_j*xhat) - &
730 rcuti4*(qxx_j*6.0_dp*cx_j*d(1)))
731 dudux_j(2) = dudux_j(2) + pref * (ri3*(qxx_j*6.0_dp*cx_j*yhat) - &
732 rcuti4*(qxx_j*6.0_dp*cx_j*d(2)))
733 dudux_j(3) = dudux_j(3) + pref * (ri3*(qxx_j*6.0_dp*cx_j*zhat) - &
734 rcuti4*(qxx_j*6.0_dp*cx_j*d(3)))
735
736 duduy_j(1) = duduy_j(1) + pref * (ri3*(qyy_j*6.0_dp*cy_j*xhat) - &
737 rcuti4*(qyy_j*6.0_dp*cx_j*d(1)))
738 duduy_j(2) = duduy_j(2) + pref * (ri3*(qyy_j*6.0_dp*cy_j*yhat) - &
739 rcuti4*(qyy_j*6.0_dp*cx_j*d(2)))
740 duduy_j(3) = duduy_j(3) + pref * (ri3*(qyy_j*6.0_dp*cy_j*zhat) - &
741 rcuti4*(qyy_j*6.0_dp*cx_j*d(3)))
742
743 duduz_j(1) = duduz_j(1) + pref * (ri3*(qzz_j*6.0_dp*cz_j*xhat) - &
744 rcuti4*(qzz_j*6.0_dp*cx_j*d(1)))
745 duduz_j(2) = duduz_j(2) + pref * (ri3*(qzz_j*6.0_dp*cz_j*yhat) - &
746 rcuti4*(qzz_j*6.0_dp*cx_j*d(2)))
747 duduz_j(3) = duduz_j(3) + pref * (ri3*(qzz_j*6.0_dp*cz_j*zhat) - &
748 rcuti4*(qzz_j*6.0_dp*cx_j*d(3)))
749
750 else
751 vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
752 qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
753 qzz_j * (3.0_dp*cz2 - 1.0_dp))
754 vpair = vpair + swi * vterm
755 epot = epot + vterm
756
757 dudx = dudx - 5.0_dp*vterm*riji*xhat + pref * ri4 * ( &
758 qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
759 qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
760 qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
761 dudy = dudy - 5.0_dp*vterm*riji*yhat + pref * ri4 * ( &
762 qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
763 qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
764 qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
765 dudz = dudz - 5.0_dp*vterm*riji*zhat + pref * ri4 * ( &
766 qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
767 qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
768 qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
769
770 dudux_j(1) = dudux_j(1) + pref * ri3*(qxx_j*6.0_dp*cx_j*xhat)
771 dudux_j(2) = dudux_j(2) + pref * ri3*(qxx_j*6.0_dp*cx_j*yhat)
772 dudux_j(3) = dudux_j(3) + pref * ri3*(qxx_j*6.0_dp*cx_j*zhat)
773
774 duduy_j(1) = duduy_j(1) + pref * ri3*(qyy_j*6.0_dp*cy_j*xhat)
775 duduy_j(2) = duduy_j(2) + pref * ri3*(qyy_j*6.0_dp*cy_j*yhat)
776 duduy_j(3) = duduy_j(3) + pref * ri3*(qyy_j*6.0_dp*cy_j*zhat)
777
778 duduz_j(1) = duduz_j(1) + pref * ri3*(qzz_j*6.0_dp*cz_j*xhat)
779 duduz_j(2) = duduz_j(2) + pref * ri3*(qzz_j*6.0_dp*cz_j*yhat)
780 duduz_j(3) = duduz_j(3) + pref * ri3*(qzz_j*6.0_dp*cz_j*zhat)
781
782 endif
783 endif
784 endif
785
786 if (i_is_Dipole) then
787
788 if (j_is_Charge) then
789
790 pref = sw * pre12 * q_j * mu_i
791
792 if (corrMethod .eq. 1) then
793 ri2 = riji * riji
794 ri3 = ri2 * riji
795
796 vterm = pref * ct_i * (ri2 - rcuti2)
797 vpair = vpair + swi * vterm
798 epot = epot + vterm
799
800 !! this has a + sign in the () because the rij vector is
801 !! r_j - r_i and the charge-dipole potential takes the origin
802 !! as the point dipole, which is atom j in this case.
803
804 dudx = dudx + pref * ( ri3*( uz_i(1) - 3.0d0*ct_i*xhat) &
805 - rcuti3*( uz_i(1) - 3.0d0*ct_i*d(1)*rcuti ) )
806 dudy = dudy + pref * ( ri3*( uz_i(2) - 3.0d0*ct_i*yhat) &
807 - rcuti3*( uz_i(2) - 3.0d0*ct_i*d(2)*rcuti ) )
808 dudz = dudz + pref * ( ri3*( uz_i(3) - 3.0d0*ct_i*zhat) &
809 - rcuti3*( uz_i(3) - 3.0d0*ct_i*d(3)*rcuti ) )
810
811 duduz_i(1) = duduz_i(1) - pref*( ri2*xhat - d(1)*rcuti3 )
812 duduz_i(2) = duduz_i(2) - pref*( ri2*yhat - d(2)*rcuti3 )
813 duduz_i(3) = duduz_i(3) - pref*( ri2*zhat - d(3)*rcuti3 )
814
815 else
816 if (i_is_SplitDipole) then
817 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
818 ri = 1.0_dp / BigR
819 scale = rij * ri
820 else
821 ri = riji
822 scale = 1.0_dp
823 endif
824
825 ri2 = ri * ri
826 ri3 = ri2 * ri
827 sc2 = scale * scale
828
829 vterm = pref * ct_i * ri2 * scale
830 vpair = vpair + swi * vterm
831 epot = epot + vterm
832
833 dudx = dudx + pref * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
834 dudy = dudy + pref * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
835 dudz = dudz + pref * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
836
837 duduz_i(1) = duduz_i(1) + pref * ri2 * xhat * scale
838 duduz_i(2) = duduz_i(2) + pref * ri2 * yhat * scale
839 duduz_i(3) = duduz_i(3) + pref * ri2 * zhat * scale
840 endif
841 endif
842
843 if (j_is_Dipole) then
844
845 pref = sw * pre22 * mu_i * mu_j
846
847 if (corrMethod .eq. 1) then
848 ri2 = riji * riji
849 ri3 = ri2 * riji
850 ri4 = ri2 * ri2
851
852 vterm = pref * (ri3 - rcuti3) * (ct_ij - 3.0d0 * ct_i * ct_j)
853 vpair = vpair + swi * vterm
854 epot = epot + vterm
855
856 a1 = 5.0d0 * ct_i * ct_j - ct_ij
857
858 dudx = dudx + pref*3.0d0*ri4 &
859 *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) - &
860 pref*3.0d0*rcuti4*(a1*rcuti*d(1)-ct_i*uz_j(1)-ct_j*uz_i(1))
861 dudy = dudy + pref*3.0d0*ri4 &
862 *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) - &
863 pref*3.0d0*rcuti4*(a1*rcuti*d(2)-ct_i*uz_j(2)-ct_j*uz_i(2))
864 dudz = dudz + pref*3.0d0*ri4 &
865 *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) - &
866 pref*3.0d0*rcuti4*(a1*rcuti*d(3)-ct_i*uz_j(3)-ct_j*uz_i(3))
867
868 duduz_i(1) = duduz_i(1) + pref*(ri3*(uz_j(1) - 3.0d0*ct_j*xhat) &
869 - rcuti3*(uz_j(1) - 3.0d0*ct_j*d(1)*rcuti))
870 duduz_i(2) = duduz_i(2) + pref*(ri3*(uz_j(2) - 3.0d0*ct_j*yhat) &
871 - rcuti3*(uz_j(2) - 3.0d0*ct_j*d(2)*rcuti))
872 duduz_i(3) = duduz_i(3) + pref*(ri3*(uz_j(3) - 3.0d0*ct_j*zhat) &
873 - rcuti3*(uz_j(3) - 3.0d0*ct_j*d(3)*rcuti))
874 duduz_j(1) = duduz_j(1) + pref*(ri3*(uz_i(1) - 3.0d0*ct_i*xhat) &
875 - rcuti3*(uz_i(1) - 3.0d0*ct_i*d(1)*rcuti))
876 duduz_j(2) = duduz_j(2) + pref*(ri3*(uz_i(2) - 3.0d0*ct_i*yhat) &
877 - rcuti3*(uz_i(2) - 3.0d0*ct_i*d(2)*rcuti))
878 duduz_j(3) = duduz_j(3) + pref*(ri3*(uz_i(3) - 3.0d0*ct_i*zhat) &
879 - rcuti3*(uz_i(3) - 3.0d0*ct_i*d(3)*rcuti))
880 else
881
882 if (i_is_SplitDipole) then
883 if (j_is_SplitDipole) then
884 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
885 else
886 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
887 endif
888 ri = 1.0_dp / BigR
889 scale = rij * ri
890 else
891 if (j_is_SplitDipole) then
892 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
893 ri = 1.0_dp / BigR
894 scale = rij * ri
895 else
896 ri = riji
897 scale = 1.0_dp
898 endif
899 endif
900
901 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
902
903 ri2 = ri * ri
904 ri3 = ri2 * ri
905 ri4 = ri2 * ri2
906 sc2 = scale * scale
907
908 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
909 vpair = vpair + swi * vterm
910 epot = epot + vterm
911
912 a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
913
914 dudx = dudx + pref*3.0d0*ri4*scale &
915 *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
916 dudy = dudy + pref*3.0d0*ri4*scale &
917 *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
918 dudz = dudz + pref*3.0d0*ri4*scale &
919 *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
920
921 duduz_i(1) = duduz_i(1) + pref*ri3 &
922 *(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
923 duduz_i(2) = duduz_i(2) + pref*ri3 &
924 *(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
925 duduz_i(3) = duduz_i(3) + pref*ri3 &
926 *(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
927
928 duduz_j(1) = duduz_j(1) + pref*ri3 &
929 *(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
930 duduz_j(2) = duduz_j(2) + pref*ri3 &
931 *(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
932 duduz_j(3) = duduz_j(3) + pref*ri3 &
933 *(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
934 endif
935 endif
936 endif
937
938 if (i_is_Quadrupole) then
939 if (j_is_Charge) then
940
941 ri2 = riji * riji
942 ri3 = ri2 * riji
943 ri4 = ri2 * ri2
944 cx2 = cx_i * cx_i
945 cy2 = cy_i * cy_i
946 cz2 = cz_i * cz_i
947
948 pref = sw * pre14 * q_j / 3.0_dp
949
950 if (corrMethod .eq. 1) then
951 vterm1 = pref * ri3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
952 qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
953 qzz_i * (3.0_dp*cz2 - 1.0_dp) )
954 vterm2 = pref * rcuti3*( qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
955 qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
956 qzz_i * (3.0_dp*cz2 - 1.0_dp) )
957 vpair = vpair + swi * ( vterm1 - vterm2 )
958 epot = epot + ( vterm1 - vterm2 )
959
960 dudx = dudx - (5.0_dp*(vterm1*riji*xhat - vterm2*rcuti2*d(1))) + &
961 pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(1)) - &
962 qxx_i*2.0_dp*(xhat - rcuti*d(1))) + &
963 (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(1)) - &
964 qyy_i*2.0_dp*(xhat - rcuti*d(1))) + &
965 (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(1)) - &
966 qzz_i*2.0_dp*(xhat - rcuti*d(1))) )
967 dudy = dudy - (5.0_dp*(vterm1*riji*yhat - vterm2*rcuti2*d(2))) + &
968 pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(2)) - &
969 qxx_i*2.0_dp*(yhat - rcuti*d(2))) + &
970 (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(2)) - &
971 qyy_i*2.0_dp*(yhat - rcuti*d(2))) + &
972 (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(2)) - &
973 qzz_i*2.0_dp*(yhat - rcuti*d(2))) )
974 dudz = dudz - (5.0_dp*(vterm1*riji*zhat - vterm2*rcuti2*d(3))) + &
975 pref * ( (ri4 - rcuti4)*(qxx_i*(6.0_dp*cx_i*ux_i(3)) - &
976 qxx_i*2.0_dp*(zhat - rcuti*d(3))) + &
977 (ri4 - rcuti4)*(qyy_i*(6.0_dp*cy_i*uy_i(3)) - &
978 qyy_i*2.0_dp*(zhat - rcuti*d(3))) + &
979 (ri4 - rcuti4)*(qzz_i*(6.0_dp*cz_i*uz_i(3)) - &
980 qzz_i*2.0_dp*(zhat - rcuti*d(3))) )
981
982 dudux_i(1) = dudux_i(1) + pref * (ri3*(qxx_i*6.0_dp*cx_i*xhat) - &
983 rcuti4*(qxx_i*6.0_dp*cx_i*d(1)))
984 dudux_i(2) = dudux_i(2) + pref * (ri3*(qxx_i*6.0_dp*cx_i*yhat) - &
985 rcuti4*(qxx_i*6.0_dp*cx_i*d(2)))
986 dudux_i(3) = dudux_i(3) + pref * (ri3*(qxx_i*6.0_dp*cx_i*zhat) - &
987 rcuti4*(qxx_i*6.0_dp*cx_i*d(3)))
988
989 duduy_i(1) = duduy_i(1) + pref * (ri3*(qyy_i*6.0_dp*cy_i*xhat) - &
990 rcuti4*(qyy_i*6.0_dp*cx_i*d(1)))
991 duduy_i(2) = duduy_i(2) + pref * (ri3*(qyy_i*6.0_dp*cy_i*yhat) - &
992 rcuti4*(qyy_i*6.0_dp*cx_i*d(2)))
993 duduy_i(3) = duduy_i(3) + pref * (ri3*(qyy_i*6.0_dp*cy_i*zhat) - &
994 rcuti4*(qyy_i*6.0_dp*cx_i*d(3)))
995
996 duduz_i(1) = duduz_i(1) + pref * (ri3*(qzz_i*6.0_dp*cz_i*xhat) - &
997 rcuti4*(qzz_i*6.0_dp*cx_i*d(1)))
998 duduz_i(2) = duduz_i(2) + pref * (ri3*(qzz_i*6.0_dp*cz_i*yhat) - &
999 rcuti4*(qzz_i*6.0_dp*cx_i*d(2)))
1000 duduz_i(3) = duduz_i(3) + pref * (ri3*(qzz_i*6.0_dp*cz_i*zhat) - &
1001 rcuti4*(qzz_i*6.0_dp*cx_i*d(3)))
1002
1003 else
1004 vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
1005 qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
1006 qzz_i * (3.0_dp*cz2 - 1.0_dp))
1007 vpair = vpair + swi * vterm
1008 epot = epot + vterm
1009
1010 dudx = dudx - 5.0_dp*vterm*riji*xhat + pref * ri4 * ( &
1011 qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
1012 qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
1013 qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
1014 dudy = dudy - 5.0_dp*vterm*riji*yhat + pref * ri4 * ( &
1015 qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
1016 qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
1017 qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
1018 dudz = dudz - 5.0_dp*vterm*riji*zhat + pref * ri4 * ( &
1019 qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
1020 qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
1021 qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
1022
1023 dudux_i(1) = dudux_i(1) + pref * ri3*(qxx_i*6.0_dp*cx_i*xhat)
1024 dudux_i(2) = dudux_i(2) + pref * ri3*(qxx_i*6.0_dp*cx_i*yhat)
1025 dudux_i(3) = dudux_i(3) + pref * ri3*(qxx_i*6.0_dp*cx_i*zhat)
1026
1027 duduy_i(1) = duduy_i(1) + pref * ri3*(qyy_i*6.0_dp*cy_i*xhat)
1028 duduy_i(2) = duduy_i(2) + pref * ri3*(qyy_i*6.0_dp*cy_i*yhat)
1029 duduy_i(3) = duduy_i(3) + pref * ri3*(qyy_i*6.0_dp*cy_i*zhat)
1030
1031 duduz_i(1) = duduz_i(1) + pref * ri3*(qzz_i*6.0_dp*cz_i*xhat)
1032 duduz_i(2) = duduz_i(2) + pref * ri3*(qzz_i*6.0_dp*cz_i*yhat)
1033 duduz_i(3) = duduz_i(3) + pref * ri3*(qzz_i*6.0_dp*cz_i*zhat)
1034 endif
1035 endif
1036 endif
1037
1038
1039 if (do_pot) then
1040 #ifdef IS_MPI
1041 pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
1042 pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
1043 #else
1044 pot = pot + epot
1045 #endif
1046 endif
1047
1048 #ifdef IS_MPI
1049 f_Row(1,atom1) = f_Row(1,atom1) + dudx
1050 f_Row(2,atom1) = f_Row(2,atom1) + dudy
1051 f_Row(3,atom1) = f_Row(3,atom1) + dudz
1052
1053 f_Col(1,atom2) = f_Col(1,atom2) - dudx
1054 f_Col(2,atom2) = f_Col(2,atom2) - dudy
1055 f_Col(3,atom2) = f_Col(3,atom2) - dudz
1056
1057 if (i_is_Dipole .or. i_is_Quadrupole) then
1058 t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1059 t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1060 t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1061 endif
1062 if (i_is_Quadrupole) then
1063 t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
1064 t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
1065 t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
1066
1067 t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
1068 t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
1069 t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
1070 endif
1071
1072 if (j_is_Dipole .or. j_is_Quadrupole) then
1073 t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1074 t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1075 t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1076 endif
1077 if (j_is_Quadrupole) then
1078 t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
1079 t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
1080 t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
1081
1082 t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
1083 t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
1084 t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
1085 endif
1086
1087 #else
1088 f(1,atom1) = f(1,atom1) + dudx
1089 f(2,atom1) = f(2,atom1) + dudy
1090 f(3,atom1) = f(3,atom1) + dudz
1091
1092 f(1,atom2) = f(1,atom2) - dudx
1093 f(2,atom2) = f(2,atom2) - dudy
1094 f(3,atom2) = f(3,atom2) - dudz
1095
1096 if (i_is_Dipole .or. i_is_Quadrupole) then
1097 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
1098 t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
1099 t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
1100 endif
1101 if (i_is_Quadrupole) then
1102 t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
1103 t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
1104 t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
1105
1106 t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
1107 t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
1108 t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
1109 endif
1110
1111 if (j_is_Dipole .or. j_is_Quadrupole) then
1112 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
1113 t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
1114 t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
1115 endif
1116 if (j_is_Quadrupole) then
1117 t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
1118 t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
1119 t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
1120
1121 t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
1122 t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
1123 t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
1124 endif
1125
1126 #endif
1127
1128 #ifdef IS_MPI
1129 id1 = AtomRowToGlobal(atom1)
1130 id2 = AtomColToGlobal(atom2)
1131 #else
1132 id1 = atom1
1133 id2 = atom2
1134 #endif
1135
1136 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
1137
1138 fpair(1) = fpair(1) + dudx
1139 fpair(2) = fpair(2) + dudy
1140 fpair(3) = fpair(3) + dudz
1141
1142 endif
1143
1144 return
1145 end subroutine doElectrostaticPair
1146
1147 !! calculates the switching functions and their derivatives for a given
1148 subroutine calc_switch(r, mu, scale, dscale)
1149
1150 real (kind=dp), intent(in) :: r, mu
1151 real (kind=dp), intent(inout) :: scale, dscale
1152 real (kind=dp) :: rl, ru, mulow, minRatio, temp, scaleVal
1153
1154 ! distances must be in angstroms
1155 rl = 2.75d0
1156 ru = 3.75d0
1157 mulow = 0.0d0 !3.3856d0 ! 1.84 * 1.84
1158 minRatio = mulow / (mu*mu)
1159 scaleVal = 1.0d0 - minRatio
1160
1161 if (r.lt.rl) then
1162 scale = minRatio
1163 dscale = 0.0d0
1164 elseif (r.gt.ru) then
1165 scale = 1.0d0
1166 dscale = 0.0d0
1167 else
1168 scale = 1.0d0 - scaleVal*((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) &
1169 / ((ru - rl)**3)
1170 dscale = -scaleVal * 6.0d0 * (r-ru)*(r-rl)/((ru - rl)**3)
1171 endif
1172
1173 return
1174 end subroutine calc_switch
1175
1176 subroutine destroyElectrostaticTypes()
1177
1178 if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap)
1179
1180 end subroutine destroyElectrostaticTypes
1181
1182 end module electrostatic_module