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