ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2306
Committed: Fri Sep 16 20:37:05 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 44968 byte(s)
Log Message:
fixing some summation method issues

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