ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2302
Committed: Fri Sep 16 16:07:39 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 44728 byte(s)
Log Message:
some fixes but even more breaking (cutting out the old way to do reaction field)

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