ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2129
Committed: Mon Mar 21 20:51:10 2005 UTC (19 years, 3 months ago) by chrisfen
File size: 29093 byte(s)
Log Message:
Make sure electrostatic_module provides data for 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 !! 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) :: drdxj, drdyj, drdzj
343 real (kind=dp) :: scale, sc2, bigR
344
345 if (.not.allocated(ElectrostaticMap)) then
346 call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
347 return
348 end if
349
350 #ifdef IS_MPI
351 me1 = atid_Row(atom1)
352 me2 = atid_Col(atom2)
353 #else
354 me1 = atid(atom1)
355 me2 = atid(atom2)
356 #endif
357
358 !! some variables we'll need independent of electrostatic type:
359
360 riji = 1.0d0 / rij
361
362 xhat = d(1) * riji
363 yhat = d(2) * riji
364 zhat = d(3) * riji
365
366 drdxj = xhat
367 drdyj = yhat
368 drdzj = zhat
369
370 !! logicals
371
372 i_is_Charge = ElectrostaticMap(me1)%is_Charge
373 i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
374 i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
375 i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
376
377 j_is_Charge = ElectrostaticMap(me2)%is_Charge
378 j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
379 j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
380 j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
381
382 if (i_is_Charge) then
383 q_i = ElectrostaticMap(me1)%charge
384 endif
385
386 if (i_is_Dipole) then
387 mu_i = ElectrostaticMap(me1)%dipole_moment
388 #ifdef IS_MPI
389 uz_i(1) = eFrame_Row(3,atom1)
390 uz_i(2) = eFrame_Row(6,atom1)
391 uz_i(3) = eFrame_Row(9,atom1)
392 #else
393 uz_i(1) = eFrame(3,atom1)
394 uz_i(2) = eFrame(6,atom1)
395 uz_i(3) = eFrame(9,atom1)
396 #endif
397 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
398
399 if (i_is_SplitDipole) then
400 d_i = ElectrostaticMap(me1)%split_dipole_distance
401 endif
402
403 endif
404
405 if (i_is_Quadrupole) then
406 qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
407 qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
408 qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
409 #ifdef IS_MPI
410 ux_i(1) = eFrame_Row(1,atom1)
411 ux_i(2) = eFrame_Row(4,atom1)
412 ux_i(3) = eFrame_Row(7,atom1)
413 uy_i(1) = eFrame_Row(2,atom1)
414 uy_i(2) = eFrame_Row(5,atom1)
415 uy_i(3) = eFrame_Row(8,atom1)
416 uz_i(1) = eFrame_Row(3,atom1)
417 uz_i(2) = eFrame_Row(6,atom1)
418 uz_i(3) = eFrame_Row(9,atom1)
419 #else
420 ux_i(1) = eFrame(1,atom1)
421 ux_i(2) = eFrame(4,atom1)
422 ux_i(3) = eFrame(7,atom1)
423 uy_i(1) = eFrame(2,atom1)
424 uy_i(2) = eFrame(5,atom1)
425 uy_i(3) = eFrame(8,atom1)
426 uz_i(1) = eFrame(3,atom1)
427 uz_i(2) = eFrame(6,atom1)
428 uz_i(3) = eFrame(9,atom1)
429 #endif
430 cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
431 cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
432 cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
433 endif
434
435
436 if (j_is_Charge) then
437 q_j = ElectrostaticMap(me2)%charge
438 endif
439
440 if (j_is_Dipole) then
441 mu_j = ElectrostaticMap(me2)%dipole_moment
442 #ifdef IS_MPI
443 uz_j(1) = eFrame_Col(3,atom2)
444 uz_j(2) = eFrame_Col(6,atom2)
445 uz_j(3) = eFrame_Col(9,atom2)
446 #else
447 uz_j(1) = eFrame(3,atom2)
448 uz_j(2) = eFrame(6,atom2)
449 uz_j(3) = eFrame(9,atom2)
450 #endif
451 ct_j = uz_j(1)*drdxj + uz_j(2)*drdyj + uz_j(3)*drdzj
452
453 if (j_is_SplitDipole) then
454 d_j = ElectrostaticMap(me2)%split_dipole_distance
455 endif
456 endif
457
458 if (j_is_Quadrupole) then
459 qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
460 qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
461 qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
462 #ifdef IS_MPI
463 ux_j(1) = eFrame_Col(1,atom2)
464 ux_j(2) = eFrame_Col(4,atom2)
465 ux_j(3) = eFrame_Col(7,atom2)
466 uy_j(1) = eFrame_Col(2,atom2)
467 uy_j(2) = eFrame_Col(5,atom2)
468 uy_j(3) = eFrame_Col(8,atom2)
469 uz_j(1) = eFrame_Col(3,atom2)
470 uz_j(2) = eFrame_Col(6,atom2)
471 uz_j(3) = eFrame_Col(9,atom2)
472 #else
473 ux_j(1) = eFrame(1,atom2)
474 ux_j(2) = eFrame(4,atom2)
475 ux_j(3) = eFrame(7,atom2)
476 uy_j(1) = eFrame(2,atom2)
477 uy_j(2) = eFrame(5,atom2)
478 uy_j(3) = eFrame(8,atom2)
479 uz_j(1) = eFrame(3,atom2)
480 uz_j(2) = eFrame(6,atom2)
481 uz_j(3) = eFrame(9,atom2)
482 #endif
483 cx_j = ux_j(1)*xhat + ux_j(2)*yhat + ux_j(3)*zhat
484 cy_j = uy_j(1)*xhat + uy_j(2)*yhat + uy_j(3)*zhat
485 cz_j = uz_j(1)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
486 endif
487
488 epot = 0.0_dp
489 dudx = 0.0_dp
490 dudy = 0.0_dp
491 dudz = 0.0_dp
492
493 dudux_i = 0.0_dp
494 duduy_i = 0.0_dp
495 duduz_i = 0.0_dp
496
497 dudux_j = 0.0_dp
498 duduy_j = 0.0_dp
499 duduz_j = 0.0_dp
500
501 if (i_is_Charge) then
502
503 if (j_is_Charge) then
504
505 vterm = pre11 * q_i * q_j * riji
506 vpair = vpair + vterm
507 epot = epot + sw*vterm
508
509 dudr = - sw * vterm * riji
510
511 dudx = dudx + dudr * drdxj
512 dudy = dudy + dudr * drdyj
513 dudz = dudz + dudr * drdzj
514
515 endif
516
517 if (j_is_Dipole) then
518
519 if (j_is_SplitDipole) then
520 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
521 ri = 1.0_dp / BigR
522 scale = rij * ri
523 else
524 ri = riji
525 scale = 1.0_dp
526 endif
527
528 ri2 = ri * ri
529 ri3 = ri2 * ri
530 sc2 = scale * scale
531
532 pref = pre12 * q_i * mu_j
533 vterm = pref * ct_j * ri2 * scale
534 vpair = vpair + vterm
535 epot = epot + sw * vterm
536
537 !! this has a + sign in the () because the rij vector is
538 !! r_j - r_i and the charge-dipole potential takes the origin
539 !! as the point dipole, which is atom j in this case.
540
541 dudx = dudx + pref * sw * ri3 * ( uz_j(1) + 3.0d0*ct_j*xhat*sc2)
542 dudy = dudy + pref * sw * ri3 * ( uz_j(2) + 3.0d0*ct_j*yhat*sc2)
543 dudz = dudz + pref * sw * ri3 * ( uz_j(3) + 3.0d0*ct_j*zhat*sc2)
544
545 duduz_j(1) = duduz_j(1) - pref * sw * ri2 * xhat * scale
546 duduz_j(2) = duduz_j(2) - pref * sw * ri2 * yhat * scale
547 duduz_j(3) = duduz_j(3) - pref * sw * ri2 * zhat * scale
548
549 endif
550
551 if (j_is_Quadrupole) then
552 ri2 = riji * riji
553 ri3 = ri2 * riji
554 ri4 = ri2 * ri2
555 cx2 = cx_j * cx_j
556 cy2 = cy_j * cy_j
557 cz2 = cz_j * cz_j
558
559
560 pref = pre14 * q_i / 6.0_dp
561 vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
562 qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
563 qzz_j * (3.0_dp*cz2 - 1.0_dp))
564 vpair = vpair + vterm
565 epot = epot + sw * vterm
566
567 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat - pref * sw * ri4 * ( &
568 qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
569 qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
570 qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
571 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat - pref * sw * ri4 * ( &
572 qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
573 qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
574 qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
575 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat - pref * sw * ri4 * ( &
576 qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
577 qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
578 qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
579
580 dudux_j(1) = dudux_j(1) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*xhat)
581 dudux_j(2) = dudux_j(2) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*yhat)
582 dudux_j(3) = dudux_j(3) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*zhat)
583
584 duduy_j(1) = duduy_j(1) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*xhat)
585 duduy_j(2) = duduy_j(2) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*yhat)
586 duduy_j(3) = duduy_j(3) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*zhat)
587
588 duduz_j(1) = duduz_j(1) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*xhat)
589 duduz_j(2) = duduz_j(2) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*yhat)
590 duduz_j(3) = duduz_j(3) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*zhat)
591 endif
592
593 endif
594
595 if (i_is_Dipole) then
596
597 if (j_is_Charge) then
598
599 if (i_is_SplitDipole) then
600 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
601 ri = 1.0_dp / BigR
602 scale = rij * ri
603 else
604 ri = riji
605 scale = 1.0_dp
606 endif
607
608 ri2 = ri * ri
609 ri3 = ri2 * ri
610 sc2 = scale * scale
611
612 pref = pre12 * q_j * mu_i
613 vterm = pref * ct_i * ri2 * scale
614 vpair = vpair + vterm
615 epot = epot + sw * vterm
616
617 dudx = dudx + pref * sw * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
618 dudy = dudy + pref * sw * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
619 dudz = dudz + pref * sw * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
620
621 duduz_i(1) = duduz_i(1) + pref * sw * ri2 * xhat * scale
622 duduz_i(2) = duduz_i(2) + pref * sw * ri2 * yhat * scale
623 duduz_i(3) = duduz_i(3) + pref * sw * ri2 * zhat * scale
624 endif
625
626 if (j_is_Dipole) then
627
628 if (i_is_SplitDipole) then
629 if (j_is_SplitDipole) then
630 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
631 else
632 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
633 endif
634 ri = 1.0_dp / BigR
635 scale = rij * ri
636 else
637 if (j_is_SplitDipole) then
638 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
639 ri = 1.0_dp / BigR
640 scale = rij * ri
641 else
642 ri = riji
643 scale = 1.0_dp
644 endif
645 endif
646
647 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
648
649 ri2 = ri * ri
650 ri3 = ri2 * ri
651 ri4 = ri2 * ri2
652 sc2 = scale * scale
653
654 pref = pre22 * mu_i * mu_j
655 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
656 vpair = vpair + vterm
657 epot = epot + sw * vterm
658
659 a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
660
661 dudx=dudx+pref*sw*3.0d0*ri4*scale*(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
662 dudy=dudy+pref*sw*3.0d0*ri4*scale*(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
663 dudz=dudz+pref*sw*3.0d0*ri4*scale*(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
664
665 duduz_i(1) = duduz_i(1) + pref*sw*ri3*(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
666 duduz_i(2) = duduz_i(2) + pref*sw*ri3*(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
667 duduz_i(3) = duduz_i(3) + pref*sw*ri3*(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
668
669 duduz_j(1) = duduz_j(1) + pref*sw*ri3*(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
670 duduz_j(2) = duduz_j(2) + pref*sw*ri3*(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
671 duduz_j(3) = duduz_j(3) + pref*sw*ri3*(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
672 endif
673
674 endif
675
676 if (i_is_Quadrupole) then
677 if (j_is_Charge) then
678
679 ri2 = riji * riji
680 ri3 = ri2 * riji
681 ri4 = ri2 * ri2
682 cx2 = cx_i * cx_i
683 cy2 = cy_i * cy_i
684 cz2 = cz_i * cz_i
685
686 pref = pre14 * q_j / 6.0_dp
687 vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
688 qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
689 qzz_i * (3.0_dp*cz2 - 1.0_dp))
690 vpair = vpair + vterm
691 epot = epot + sw * vterm
692
693 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat - pref * sw * ri4 * ( &
694 qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
695 qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
696 qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
697 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat - pref * sw * ri4 * ( &
698 qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
699 qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
700 qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
701 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat - pref * sw * ri4 * ( &
702 qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
703 qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
704 qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
705
706 dudux_i(1) = dudux_i(1) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*xhat)
707 dudux_i(2) = dudux_i(2) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*yhat)
708 dudux_i(3) = dudux_i(3) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*zhat)
709
710 duduy_i(1) = duduy_i(1) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*xhat)
711 duduy_i(2) = duduy_i(2) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*yhat)
712 duduy_i(3) = duduy_i(3) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*zhat)
713
714 duduz_i(1) = duduz_i(1) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*xhat)
715 duduz_i(2) = duduz_i(2) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*yhat)
716 duduz_i(3) = duduz_i(3) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*zhat)
717 endif
718 endif
719
720
721 if (do_pot) then
722 #ifdef IS_MPI
723 pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
724 pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
725 #else
726 pot = pot + epot
727 #endif
728 endif
729
730 #ifdef IS_MPI
731 f_Row(1,atom1) = f_Row(1,atom1) + dudx
732 f_Row(2,atom1) = f_Row(2,atom1) + dudy
733 f_Row(3,atom1) = f_Row(3,atom1) + dudz
734
735 f_Col(1,atom2) = f_Col(1,atom2) - dudx
736 f_Col(2,atom2) = f_Col(2,atom2) - dudy
737 f_Col(3,atom2) = f_Col(3,atom2) - dudz
738
739 if (i_is_Dipole .or. i_is_Quadrupole) then
740 t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
741 t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
742 t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
743 endif
744 if (i_is_Quadrupole) then
745 t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
746 t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
747 t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
748
749 t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
750 t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
751 t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
752 endif
753
754 if (j_is_Dipole .or. j_is_Quadrupole) then
755 t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
756 t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
757 t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
758 endif
759 if (j_is_Quadrupole) then
760 t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
761 t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
762 t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
763
764 t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
765 t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
766 t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
767 endif
768
769 #else
770 f(1,atom1) = f(1,atom1) + dudx
771 f(2,atom1) = f(2,atom1) + dudy
772 f(3,atom1) = f(3,atom1) + dudz
773
774 f(1,atom2) = f(1,atom2) - dudx
775 f(2,atom2) = f(2,atom2) - dudy
776 f(3,atom2) = f(3,atom2) - dudz
777
778 if (i_is_Dipole .or. i_is_Quadrupole) then
779 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
780 t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
781 t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
782 endif
783 if (i_is_Quadrupole) then
784 t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
785 t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
786 t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
787
788 t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
789 t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
790 t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
791 endif
792
793 if (j_is_Dipole .or. j_is_Quadrupole) then
794 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
795 t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
796 t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
797 endif
798 if (j_is_Quadrupole) then
799 t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
800 t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
801 t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
802
803 t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
804 t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
805 t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
806 endif
807
808 #endif
809
810 #ifdef IS_MPI
811 id1 = AtomRowToGlobal(atom1)
812 id2 = AtomColToGlobal(atom2)
813 #else
814 id1 = atom1
815 id2 = atom2
816 #endif
817
818 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
819
820 fpair(1) = fpair(1) + dudx
821 fpair(2) = fpair(2) + dudy
822 fpair(3) = fpair(3) + dudz
823
824 endif
825
826 return
827 end subroutine doElectrostaticPair
828
829 end module electrostatic_module