ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2162
Committed: Mon Apr 11 20:19:22 2005 UTC (19 years, 5 months ago) by chrisfen
File size: 28995 byte(s)
Log Message:
fixing of the quadrupoles.  look!  it's divide by 3 like stone says!

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
81 type :: Electrostatic
82 integer :: c_ident
83 logical :: is_Charge = .false.
84 logical :: is_Dipole = .false.
85 logical :: is_SplitDipole = .false.
86 logical :: is_Quadrupole = .false.
87 real(kind=DP) :: charge = 0.0_DP
88 real(kind=DP) :: dipole_moment = 0.0_DP
89 real(kind=DP) :: split_dipole_distance = 0.0_DP
90 real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP
91 end type Electrostatic
92
93 type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
94
95 contains
96
97 subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
98 is_SplitDipole, is_Quadrupole, status)
99
100 integer, intent(in) :: c_ident
101 logical, intent(in) :: is_Charge
102 logical, intent(in) :: is_Dipole
103 logical, intent(in) :: is_SplitDipole
104 logical, intent(in) :: is_Quadrupole
105 integer, intent(out) :: status
106 integer :: nAtypes, myATID, i, j
107
108 status = 0
109 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
110
111 !! Be simple-minded and assume that we need an ElectrostaticMap that
112 !! is the same size as the total number of atom types
113
114 if (.not.allocated(ElectrostaticMap)) then
115
116 nAtypes = getSize(atypes)
117
118 if (nAtypes == 0) then
119 status = -1
120 return
121 end if
122
123 if (.not. allocated(ElectrostaticMap)) then
124 allocate(ElectrostaticMap(nAtypes))
125 endif
126
127 end if
128
129 if (myATID .gt. size(ElectrostaticMap)) then
130 status = -1
131 return
132 endif
133
134 ! set the values for ElectrostaticMap for this atom type:
135
136 ElectrostaticMap(myATID)%c_ident = c_ident
137 ElectrostaticMap(myATID)%is_Charge = is_Charge
138 ElectrostaticMap(myATID)%is_Dipole = is_Dipole
139 ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole
140 ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole
141
142 end subroutine newElectrostaticType
143
144 subroutine setCharge(c_ident, charge, status)
145 integer, intent(in) :: c_ident
146 real(kind=dp), intent(in) :: charge
147 integer, intent(out) :: status
148 integer :: myATID
149
150 status = 0
151 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
152
153 if (.not.allocated(ElectrostaticMap)) then
154 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!")
155 status = -1
156 return
157 end if
158
159 if (myATID .gt. size(ElectrostaticMap)) then
160 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!")
161 status = -1
162 return
163 endif
164
165 if (.not.ElectrostaticMap(myATID)%is_Charge) then
166 call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!")
167 status = -1
168 return
169 endif
170
171 ElectrostaticMap(myATID)%charge = charge
172 end subroutine setCharge
173
174 subroutine setDipoleMoment(c_ident, dipole_moment, status)
175 integer, intent(in) :: c_ident
176 real(kind=dp), intent(in) :: dipole_moment
177 integer, intent(out) :: status
178 integer :: myATID
179
180 status = 0
181 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
182
183 if (.not.allocated(ElectrostaticMap)) then
184 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!")
185 status = -1
186 return
187 end if
188
189 if (myATID .gt. size(ElectrostaticMap)) then
190 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!")
191 status = -1
192 return
193 endif
194
195 if (.not.ElectrostaticMap(myATID)%is_Dipole) then
196 call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!")
197 status = -1
198 return
199 endif
200
201 ElectrostaticMap(myATID)%dipole_moment = dipole_moment
202 end subroutine setDipoleMoment
203
204 subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status)
205 integer, intent(in) :: c_ident
206 real(kind=dp), intent(in) :: split_dipole_distance
207 integer, intent(out) :: status
208 integer :: myATID
209
210 status = 0
211 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
212
213 if (.not.allocated(ElectrostaticMap)) then
214 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!")
215 status = -1
216 return
217 end if
218
219 if (myATID .gt. size(ElectrostaticMap)) then
220 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!")
221 status = -1
222 return
223 endif
224
225 if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then
226 call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!")
227 status = -1
228 return
229 endif
230
231 ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance
232 end subroutine setSplitDipoleDistance
233
234 subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status)
235 integer, intent(in) :: c_ident
236 real(kind=dp), intent(in), dimension(3) :: quadrupole_moments
237 integer, intent(out) :: status
238 integer :: myATID, i, j
239
240 status = 0
241 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
242
243 if (.not.allocated(ElectrostaticMap)) then
244 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!")
245 status = -1
246 return
247 end if
248
249 if (myATID .gt. size(ElectrostaticMap)) then
250 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!")
251 status = -1
252 return
253 endif
254
255 if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then
256 call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!")
257 status = -1
258 return
259 endif
260
261 do i = 1, 3
262 ElectrostaticMap(myATID)%quadrupole_moments(i) = &
263 quadrupole_moments(i)
264 enddo
265
266 end subroutine setQuadrupoleMoments
267
268
269 function getCharge(atid) result (c)
270 integer, intent(in) :: atid
271 integer :: localError
272 real(kind=dp) :: c
273
274 if (.not.allocated(ElectrostaticMap)) then
275 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
276 return
277 end if
278
279 if (.not.ElectrostaticMap(atid)%is_Charge) then
280 call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!")
281 return
282 endif
283
284 c = ElectrostaticMap(atid)%charge
285 end function getCharge
286
287 function getDipoleMoment(atid) result (dm)
288 integer, intent(in) :: atid
289 integer :: localError
290 real(kind=dp) :: dm
291
292 if (.not.allocated(ElectrostaticMap)) then
293 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
294 return
295 end if
296
297 if (.not.ElectrostaticMap(atid)%is_Dipole) then
298 call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!")
299 return
300 endif
301
302 dm = ElectrostaticMap(atid)%dipole_moment
303 end function getDipoleMoment
304
305 subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, &
306 vpair, fpair, pot, eFrame, f, t, do_pot)
307
308 logical, intent(in) :: do_pot
309
310 integer, intent(in) :: atom1, atom2
311 integer :: localError
312
313 real(kind=dp), intent(in) :: rij, r2, sw
314 real(kind=dp), intent(in), dimension(3) :: d
315 real(kind=dp), intent(inout) :: vpair
316 real(kind=dp), intent(inout), dimension(3) :: fpair
317
318 real( kind = dp ) :: pot
319 real( kind = dp ), dimension(9,nLocal) :: eFrame
320 real( kind = dp ), dimension(3,nLocal) :: f
321 real( kind = dp ), dimension(3,nLocal) :: t
322
323 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
324 real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
325 real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
326 real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
327
328 logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
329 logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
330 integer :: me1, me2, id1, id2
331 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
332 real (kind=dp) :: qxx_i, qyy_i, qzz_i
333 real (kind=dp) :: qxx_j, qyy_j, qzz_j
334 real (kind=dp) :: cx_i, cy_i, cz_i
335 real (kind=dp) :: cx_j, cy_j, cz_j
336 real (kind=dp) :: cx2, cy2, cz2
337 real (kind=dp) :: ct_i, ct_j, ct_ij, a1
338 real (kind=dp) :: riji, ri, ri2, ri3, ri4
339 real (kind=dp) :: pref, vterm, epot, dudr
340 real (kind=dp) :: xhat, yhat, zhat
341 real (kind=dp) :: dudx, dudy, dudz
342 real (kind=dp) :: scale, sc2, bigR
343
344 if (.not.allocated(ElectrostaticMap)) then
345 call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
346 return
347 end if
348
349 #ifdef IS_MPI
350 me1 = atid_Row(atom1)
351 me2 = atid_Col(atom2)
352 #else
353 me1 = atid(atom1)
354 me2 = atid(atom2)
355 #endif
356
357 !! some variables we'll need independent of electrostatic type:
358
359 riji = 1.0d0 / rij
360
361 xhat = d(1) * riji
362 yhat = d(2) * riji
363 zhat = d(3) * riji
364
365 !! logicals
366
367 i_is_Charge = ElectrostaticMap(me1)%is_Charge
368 i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
369 i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
370 i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
371
372 j_is_Charge = ElectrostaticMap(me2)%is_Charge
373 j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
374 j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
375 j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
376
377 if (i_is_Charge) then
378 q_i = ElectrostaticMap(me1)%charge
379 endif
380
381 if (i_is_Dipole) then
382 mu_i = ElectrostaticMap(me1)%dipole_moment
383 #ifdef IS_MPI
384 uz_i(1) = eFrame_Row(3,atom1)
385 uz_i(2) = eFrame_Row(6,atom1)
386 uz_i(3) = eFrame_Row(9,atom1)
387 #else
388 uz_i(1) = eFrame(3,atom1)
389 uz_i(2) = eFrame(6,atom1)
390 uz_i(3) = eFrame(9,atom1)
391 #endif
392 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
393
394 if (i_is_SplitDipole) then
395 d_i = ElectrostaticMap(me1)%split_dipole_distance
396 endif
397
398 endif
399
400 if (i_is_Quadrupole) then
401 qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
402 qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
403 qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
404 #ifdef IS_MPI
405 ux_i(1) = eFrame_Row(1,atom1)
406 ux_i(2) = eFrame_Row(4,atom1)
407 ux_i(3) = eFrame_Row(7,atom1)
408 uy_i(1) = eFrame_Row(2,atom1)
409 uy_i(2) = eFrame_Row(5,atom1)
410 uy_i(3) = eFrame_Row(8,atom1)
411 uz_i(1) = eFrame_Row(3,atom1)
412 uz_i(2) = eFrame_Row(6,atom1)
413 uz_i(3) = eFrame_Row(9,atom1)
414 #else
415 ux_i(1) = eFrame(1,atom1)
416 ux_i(2) = eFrame(4,atom1)
417 ux_i(3) = eFrame(7,atom1)
418 uy_i(1) = eFrame(2,atom1)
419 uy_i(2) = eFrame(5,atom1)
420 uy_i(3) = eFrame(8,atom1)
421 uz_i(1) = eFrame(3,atom1)
422 uz_i(2) = eFrame(6,atom1)
423 uz_i(3) = eFrame(9,atom1)
424 #endif
425 cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
426 cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
427 cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
428 endif
429
430
431 if (j_is_Charge) then
432 q_j = ElectrostaticMap(me2)%charge
433 endif
434
435 if (j_is_Dipole) then
436 mu_j = ElectrostaticMap(me2)%dipole_moment
437 #ifdef IS_MPI
438 uz_j(1) = eFrame_Col(3,atom2)
439 uz_j(2) = eFrame_Col(6,atom2)
440 uz_j(3) = eFrame_Col(9,atom2)
441 #else
442 uz_j(1) = eFrame(3,atom2)
443 uz_j(2) = eFrame(6,atom2)
444 uz_j(3) = eFrame(9,atom2)
445 #endif
446 ct_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
447
448 if (j_is_SplitDipole) then
449 d_j = ElectrostaticMap(me2)%split_dipole_distance
450 endif
451 endif
452
453 if (j_is_Quadrupole) then
454 qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
455 qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
456 qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
457 #ifdef IS_MPI
458 ux_j(1) = eFrame_Col(1,atom2)
459 ux_j(2) = eFrame_Col(4,atom2)
460 ux_j(3) = eFrame_Col(7,atom2)
461 uy_j(1) = eFrame_Col(2,atom2)
462 uy_j(2) = eFrame_Col(5,atom2)
463 uy_j(3) = eFrame_Col(8,atom2)
464 uz_j(1) = eFrame_Col(3,atom2)
465 uz_j(2) = eFrame_Col(6,atom2)
466 uz_j(3) = eFrame_Col(9,atom2)
467 #else
468 ux_j(1) = eFrame(1,atom2)
469 ux_j(2) = eFrame(4,atom2)
470 ux_j(3) = eFrame(7,atom2)
471 uy_j(1) = eFrame(2,atom2)
472 uy_j(2) = eFrame(5,atom2)
473 uy_j(3) = eFrame(8,atom2)
474 uz_j(1) = eFrame(3,atom2)
475 uz_j(2) = eFrame(6,atom2)
476 uz_j(3) = eFrame(9,atom2)
477 #endif
478 cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
479 cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
480 cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
481 endif
482
483 epot = 0.0_dp
484 dudx = 0.0_dp
485 dudy = 0.0_dp
486 dudz = 0.0_dp
487
488 dudux_i = 0.0_dp
489 duduy_i = 0.0_dp
490 duduz_i = 0.0_dp
491
492 dudux_j = 0.0_dp
493 duduy_j = 0.0_dp
494 duduz_j = 0.0_dp
495
496 if (i_is_Charge) then
497
498 if (j_is_Charge) then
499
500 vterm = pre11 * q_i * q_j * riji
501 vpair = vpair + vterm
502 epot = epot + sw*vterm
503
504 dudr = - sw * vterm * riji
505
506 dudx = dudx + dudr * xhat
507 dudy = dudy + dudr * yhat
508 dudz = dudz + dudr * zhat
509
510 endif
511
512 if (j_is_Dipole) then
513
514 if (j_is_SplitDipole) then
515 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
516 ri = 1.0_dp / BigR
517 scale = rij * ri
518 else
519 ri = riji
520 scale = 1.0_dp
521 endif
522
523 ri2 = ri * ri
524 ri3 = ri2 * ri
525 sc2 = scale * scale
526
527 pref = pre12 * q_i * mu_j
528 vterm = - pref * ct_j * ri2 * scale
529 vpair = vpair + vterm
530 epot = epot + sw * vterm
531
532 !! this has a + sign in the () because the rij vector is
533 !! r_j - r_i and the charge-dipole potential takes the origin
534 !! as the point dipole, which is atom j in this case.
535
536 dudx = dudx - pref * sw * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2)
537 dudy = dudy - pref * sw * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2)
538 dudz = dudz - pref * sw * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2)
539
540 duduz_j(1) = duduz_j(1) - pref * sw * ri2 * xhat * scale
541 duduz_j(2) = duduz_j(2) - pref * sw * ri2 * yhat * scale
542 duduz_j(3) = duduz_j(3) - pref * sw * ri2 * zhat * scale
543
544 endif
545
546 if (j_is_Quadrupole) then
547 ri2 = riji * riji
548 ri3 = ri2 * riji
549 ri4 = ri2 * ri2
550 cx2 = cx_j * cx_j
551 cy2 = cy_j * cy_j
552 cz2 = cz_j * cz_j
553
554
555 pref = pre14 * q_i / 3.0_dp
556 vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
557 qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
558 qzz_j * (3.0_dp*cz2 - 1.0_dp))
559 vpair = vpair + vterm
560 epot = epot + sw * vterm
561
562 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
563 qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
564 qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
565 qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
566 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
567 qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
568 qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
569 qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
570 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
571 qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
572 qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
573 qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
574
575 dudux_j(1) = dudux_j(1) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*xhat)
576 dudux_j(2) = dudux_j(2) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*yhat)
577 dudux_j(3) = dudux_j(3) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*zhat)
578
579 duduy_j(1) = duduy_j(1) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*xhat)
580 duduy_j(2) = duduy_j(2) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*yhat)
581 duduy_j(3) = duduy_j(3) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*zhat)
582
583 duduz_j(1) = duduz_j(1) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*xhat)
584 duduz_j(2) = duduz_j(2) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*yhat)
585 duduz_j(3) = duduz_j(3) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*zhat)
586 endif
587
588 endif
589
590 if (i_is_Dipole) then
591
592 if (j_is_Charge) then
593
594 if (i_is_SplitDipole) then
595 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
596 ri = 1.0_dp / BigR
597 scale = rij * ri
598 else
599 ri = riji
600 scale = 1.0_dp
601 endif
602
603 ri2 = ri * ri
604 ri3 = ri2 * ri
605 sc2 = scale * scale
606
607 pref = pre12 * q_j * mu_i
608 vterm = pref * ct_i * ri2 * scale
609 vpair = vpair + vterm
610 epot = epot + sw * vterm
611
612 dudx = dudx + pref * sw * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
613 dudy = dudy + pref * sw * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
614 dudz = dudz + pref * sw * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
615
616 duduz_i(1) = duduz_i(1) + pref * sw * ri2 * xhat * scale
617 duduz_i(2) = duduz_i(2) + pref * sw * ri2 * yhat * scale
618 duduz_i(3) = duduz_i(3) + pref * sw * ri2 * zhat * scale
619 endif
620
621 if (j_is_Dipole) then
622
623 if (i_is_SplitDipole) then
624 if (j_is_SplitDipole) then
625 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
626 else
627 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
628 endif
629 ri = 1.0_dp / BigR
630 scale = rij * ri
631 else
632 if (j_is_SplitDipole) then
633 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
634 ri = 1.0_dp / BigR
635 scale = rij * ri
636 else
637 ri = riji
638 scale = 1.0_dp
639 endif
640 endif
641
642 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
643
644 ri2 = ri * ri
645 ri3 = ri2 * ri
646 ri4 = ri2 * ri2
647 sc2 = scale * scale
648
649 pref = pre22 * mu_i * mu_j
650 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
651 vpair = vpair + vterm
652 epot = epot + sw * vterm
653
654 a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
655
656 dudx=dudx+pref*sw*3.0d0*ri4*scale*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
657 dudy=dudy+pref*sw*3.0d0*ri4*scale*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
658 dudz=dudz+pref*sw*3.0d0*ri4*scale*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
659
660 duduz_i(1) = duduz_i(1) + pref*sw*ri3*(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
661 duduz_i(2) = duduz_i(2) + pref*sw*ri3*(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
662 duduz_i(3) = duduz_i(3) + pref*sw*ri3*(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
663
664 duduz_j(1) = duduz_j(1) + pref*sw*ri3*(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
665 duduz_j(2) = duduz_j(2) + pref*sw*ri3*(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
666 duduz_j(3) = duduz_j(3) + pref*sw*ri3*(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
667 endif
668
669 endif
670
671 if (i_is_Quadrupole) then
672 if (j_is_Charge) then
673
674 ri2 = riji * riji
675 ri3 = ri2 * riji
676 ri4 = ri2 * ri2
677 cx2 = cx_i * cx_i
678 cy2 = cy_i * cy_i
679 cz2 = cz_i * cz_i
680
681 pref = pre14 * q_j / 3.0_dp
682 vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
683 qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
684 qzz_i * (3.0_dp*cz2 - 1.0_dp))
685 vpair = vpair + vterm
686 epot = epot + sw * vterm
687
688 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
689 qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
690 qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
691 qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
692 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
693 qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
694 qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
695 qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
696 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
697 qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
698 qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
699 qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
700
701 dudux_i(1) = dudux_i(1) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*xhat)
702 dudux_i(2) = dudux_i(2) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*yhat)
703 dudux_i(3) = dudux_i(3) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*zhat)
704
705 duduy_i(1) = duduy_i(1) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*xhat)
706 duduy_i(2) = duduy_i(2) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*yhat)
707 duduy_i(3) = duduy_i(3) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*zhat)
708
709 duduz_i(1) = duduz_i(1) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*xhat)
710 duduz_i(2) = duduz_i(2) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*yhat)
711 duduz_i(3) = duduz_i(3) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*zhat)
712 endif
713 endif
714
715
716 if (do_pot) then
717 #ifdef IS_MPI
718 pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
719 pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
720 #else
721 pot = pot + epot
722 #endif
723 endif
724
725 #ifdef IS_MPI
726 f_Row(1,atom1) = f_Row(1,atom1) + dudx
727 f_Row(2,atom1) = f_Row(2,atom1) + dudy
728 f_Row(3,atom1) = f_Row(3,atom1) + dudz
729
730 f_Col(1,atom2) = f_Col(1,atom2) - dudx
731 f_Col(2,atom2) = f_Col(2,atom2) - dudy
732 f_Col(3,atom2) = f_Col(3,atom2) - dudz
733
734 if (i_is_Dipole .or. i_is_Quadrupole) then
735 t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
736 t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
737 t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
738 endif
739 if (i_is_Quadrupole) then
740 t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
741 t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
742 t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
743
744 t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
745 t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
746 t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
747 endif
748
749 if (j_is_Dipole .or. j_is_Quadrupole) then
750 t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
751 t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
752 t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
753 endif
754 if (j_is_Quadrupole) then
755 t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
756 t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
757 t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
758
759 t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
760 t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
761 t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
762 endif
763
764 #else
765 f(1,atom1) = f(1,atom1) + dudx
766 f(2,atom1) = f(2,atom1) + dudy
767 f(3,atom1) = f(3,atom1) + dudz
768
769 f(1,atom2) = f(1,atom2) - dudx
770 f(2,atom2) = f(2,atom2) - dudy
771 f(3,atom2) = f(3,atom2) - dudz
772
773 if (i_is_Dipole .or. i_is_Quadrupole) then
774 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
775 t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
776 t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
777 endif
778 if (i_is_Quadrupole) then
779 t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
780 t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
781 t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
782
783 t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
784 t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
785 t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
786 endif
787
788 if (j_is_Dipole .or. j_is_Quadrupole) then
789 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
790 t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
791 t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
792 endif
793 if (j_is_Quadrupole) then
794 t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
795 t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
796 t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
797
798 t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
799 t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
800 t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
801 endif
802
803 #endif
804
805 #ifdef IS_MPI
806 id1 = AtomRowToGlobal(atom1)
807 id2 = AtomColToGlobal(atom2)
808 #else
809 id1 = atom1
810 id2 = atom2
811 #endif
812
813 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
814
815 fpair(1) = fpair(1) + dudx
816 fpair(2) = fpair(2) + dudy
817 fpair(3) = fpair(3) + dudz
818
819 endif
820
821 return
822 end subroutine doElectrostaticPair
823
824 end module electrostatic_module