ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 434
Committed: Fri Mar 11 15:53:18 2005 UTC (20 years, 4 months ago) by gezelter
File size: 20533 byte(s)
Log Message:
settled on a unit for quadrupoles

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
80 type :: Electrostatic
81 integer :: c_ident
82 logical :: is_Charge = .false.
83 logical :: is_Dipole = .false.
84 logical :: is_SplitDipole = .false.
85 logical :: is_Quadrupole = .false.
86 real(kind=DP) :: charge = 0.0_DP
87 real(kind=DP) :: dipole_moment = 0.0_DP
88 real(kind=DP) :: split_dipole_distance = 0.0_DP
89 real(kind=DP), dimension(3) :: quadrupole_moments = 0.0_DP
90 end type Electrostatic
91
92 type(Electrostatic), dimension(:), allocatable :: ElectrostaticMap
93
94 contains
95
96 subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
97 is_SplitDipole, is_Quadrupole, status)
98
99 integer, intent(in) :: c_ident
100 logical, intent(in) :: is_Charge
101 logical, intent(in) :: is_Dipole
102 logical, intent(in) :: is_SplitDipole
103 logical, intent(in) :: is_Quadrupole
104 integer, intent(out) :: status
105 integer :: nAtypes, myATID, i, j
106
107 status = 0
108 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
109
110 !! Be simple-minded and assume that we need an ElectrostaticMap that
111 !! is the same size as the total number of atom types
112
113 if (.not.allocated(ElectrostaticMap)) then
114
115 nAtypes = getSize(atypes)
116
117 if (nAtypes == 0) then
118 status = -1
119 return
120 end if
121
122 if (.not. allocated(ElectrostaticMap)) then
123 allocate(ElectrostaticMap(nAtypes))
124 endif
125
126 end if
127
128 if (myATID .gt. size(ElectrostaticMap)) then
129 status = -1
130 return
131 endif
132
133 ! set the values for ElectrostaticMap for this atom type:
134
135 ElectrostaticMap(myATID)%c_ident = c_ident
136 ElectrostaticMap(myATID)%is_Charge = is_Charge
137 ElectrostaticMap(myATID)%is_Dipole = is_Dipole
138 ElectrostaticMap(myATID)%is_SplitDipole = is_SplitDipole
139 ElectrostaticMap(myATID)%is_Quadrupole = is_Quadrupole
140
141 end subroutine newElectrostaticType
142
143 subroutine setCharge(c_ident, charge, status)
144 integer, intent(in) :: c_ident
145 real(kind=dp), intent(in) :: charge
146 integer, intent(out) :: status
147 integer :: myATID
148
149 status = 0
150 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
151
152 if (.not.allocated(ElectrostaticMap)) then
153 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setCharge!")
154 status = -1
155 return
156 end if
157
158 if (myATID .gt. size(ElectrostaticMap)) then
159 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setCharge!")
160 status = -1
161 return
162 endif
163
164 if (.not.ElectrostaticMap(myATID)%is_Charge) then
165 call handleError("electrostatic", "Attempt to setCharge of an atom type that is not a charge!")
166 status = -1
167 return
168 endif
169
170 ElectrostaticMap(myATID)%charge = charge
171 end subroutine setCharge
172
173 subroutine setDipoleMoment(c_ident, dipole_moment, status)
174 integer, intent(in) :: c_ident
175 real(kind=dp), intent(in) :: dipole_moment
176 integer, intent(out) :: status
177 integer :: myATID
178
179 status = 0
180 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
181
182 if (.not.allocated(ElectrostaticMap)) then
183 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setDipoleMoment!")
184 status = -1
185 return
186 end if
187
188 if (myATID .gt. size(ElectrostaticMap)) then
189 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setDipoleMoment!")
190 status = -1
191 return
192 endif
193
194 if (.not.ElectrostaticMap(myATID)%is_Dipole) then
195 call handleError("electrostatic", "Attempt to setDipoleMoment of an atom type that is not a dipole!")
196 status = -1
197 return
198 endif
199
200 ElectrostaticMap(myATID)%dipole_moment = dipole_moment
201 end subroutine setDipoleMoment
202
203 subroutine setSplitDipoleDistance(c_ident, split_dipole_distance, status)
204 integer, intent(in) :: c_ident
205 real(kind=dp), intent(in) :: split_dipole_distance
206 integer, intent(out) :: status
207 integer :: myATID
208
209 status = 0
210 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
211
212 if (.not.allocated(ElectrostaticMap)) then
213 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setSplitDipoleDistance!")
214 status = -1
215 return
216 end if
217
218 if (myATID .gt. size(ElectrostaticMap)) then
219 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setSplitDipoleDistance!")
220 status = -1
221 return
222 endif
223
224 if (.not.ElectrostaticMap(myATID)%is_SplitDipole) then
225 call handleError("electrostatic", "Attempt to setSplitDipoleDistance of an atom type that is not a splitDipole!")
226 status = -1
227 return
228 endif
229
230 ElectrostaticMap(myATID)%split_dipole_distance = split_dipole_distance
231 end subroutine setSplitDipoleDistance
232
233 subroutine setQuadrupoleMoments(c_ident, quadrupole_moments, status)
234 integer, intent(in) :: c_ident
235 real(kind=dp), intent(in), dimension(3) :: quadrupole_moments
236 integer, intent(out) :: status
237 integer :: myATID, i, j
238
239 status = 0
240 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
241
242 if (.not.allocated(ElectrostaticMap)) then
243 call handleError("electrostatic", "no ElectrostaticMap was present before first call of setQuadrupoleMoments!")
244 status = -1
245 return
246 end if
247
248 if (myATID .gt. size(ElectrostaticMap)) then
249 call handleError("electrostatic", "ElectrostaticMap was found to be too small during setQuadrupoleMoments!")
250 status = -1
251 return
252 endif
253
254 if (.not.ElectrostaticMap(myATID)%is_Quadrupole) then
255 call handleError("electrostatic", "Attempt to setQuadrupoleMoments of an atom type that is not a quadrupole!")
256 status = -1
257 return
258 endif
259
260 do i = 1, 3
261 ElectrostaticMap(myATID)%quadrupole_moments(i) = &
262 quadrupole_moments(i)
263 enddo
264
265 end subroutine setQuadrupoleMoments
266
267
268 function getCharge(atid) result (c)
269 integer, intent(in) :: atid
270 integer :: localError
271 real(kind=dp) :: c
272
273 if (.not.allocated(ElectrostaticMap)) then
274 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getCharge!")
275 return
276 end if
277
278 if (.not.ElectrostaticMap(atid)%is_Charge) then
279 call handleError("electrostatic", "getCharge was called for an atom type that isn't a charge!")
280 return
281 endif
282
283 c = ElectrostaticMap(atid)%charge
284 end function getCharge
285
286 function getDipoleMoment(atid) result (dm)
287 integer, intent(in) :: atid
288 integer :: localError
289 real(kind=dp) :: dm
290
291 if (.not.allocated(ElectrostaticMap)) then
292 call handleError("electrostatic", "no ElectrostaticMap was present before first call of getDipoleMoment!")
293 return
294 end if
295
296 if (.not.ElectrostaticMap(atid)%is_Dipole) then
297 call handleError("electrostatic", "getDipoleMoment was called for an atom type that isn't a dipole!")
298 return
299 endif
300
301 dm = ElectrostaticMap(atid)%dipole_moment
302 end function getDipoleMoment
303
304 subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, &
305 vpair, fpair, pot, eFrame, f, t, do_pot)
306
307 logical, intent(in) :: do_pot
308
309 integer, intent(in) :: atom1, atom2
310 integer :: localError
311
312 real(kind=dp), intent(in) :: rij, r2, sw
313 real(kind=dp), intent(in), dimension(3) :: d
314 real(kind=dp), intent(inout) :: vpair
315 real(kind=dp), intent(inout), dimension(3) :: fpair
316
317 real( kind = dp ) :: pot
318 real( kind = dp ), dimension(9,nLocal) :: eFrame
319 real( kind = dp ), dimension(3,nLocal) :: f
320 real( kind = dp ), dimension(3,nLocal) :: t
321
322 real (kind = dp), dimension(3) :: ul_i
323 real (kind = dp), dimension(3) :: ul_j
324
325 logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
326 logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
327 integer :: me1, me2, id1, id2
328 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
329 real (kind=dp) :: ct_i, ct_j, ct_ij, a1
330 real (kind=dp) :: riji, ri, ri2, ri3, ri4
331 real (kind=dp) :: pref, vterm, epot, dudr
332 real (kind=dp) :: xhat, yhat, zhat
333 real (kind=dp) :: dudx, dudy, dudz
334 real (kind=dp) :: drdxj, drdyj, drdzj
335 real (kind=dp) :: duduix, duduiy, duduiz, dudujx, dudujy, dudujz
336 real (kind=dp) :: scale, sc2, bigR
337
338 if (.not.allocated(ElectrostaticMap)) then
339 call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
340 return
341 end if
342
343 #ifdef IS_MPI
344 me1 = atid_Row(atom1)
345 me2 = atid_Col(atom2)
346 #else
347 me1 = atid(atom1)
348 me2 = atid(atom2)
349 #endif
350
351 !! some variables we'll need independent of electrostatic type:
352
353 riji = 1.0d0 / rij
354
355 xhat = d(1) * riji
356 yhat = d(2) * riji
357 zhat = d(3) * riji
358
359 drdxj = xhat
360 drdyj = yhat
361 drdzj = zhat
362
363 !! logicals
364
365 i_is_Charge = ElectrostaticMap(me1)%is_Charge
366 i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
367 i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
368 i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
369
370 j_is_Charge = ElectrostaticMap(me2)%is_Charge
371 j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
372 j_is_SplitDipole = ElectrostaticMap(me2)%is_SplitDipole
373 j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
374
375 if (i_is_Charge) then
376 q_i = ElectrostaticMap(me1)%charge
377 endif
378
379 if (i_is_Dipole) then
380 mu_i = ElectrostaticMap(me1)%dipole_moment
381 #ifdef IS_MPI
382 ul_i(1) = eFrame_Row(3,atom1)
383 ul_i(2) = eFrame_Row(6,atom1)
384 ul_i(3) = eFrame_Row(9,atom1)
385 #else
386 ul_i(1) = eFrame(3,atom1)
387 ul_i(2) = eFrame(6,atom1)
388 ul_i(3) = eFrame(9,atom1)
389 #endif
390 ct_i = ul_i(1)*drdxj + ul_i(2)*drdyj + ul_i(3)*drdzj
391
392 if (i_is_SplitDipole) then
393 d_i = ElectrostaticMap(me1)%split_dipole_distance
394 endif
395
396 endif
397
398 if (j_is_Charge) then
399 q_j = ElectrostaticMap(me2)%charge
400 endif
401
402 if (j_is_Dipole) then
403 mu_j = ElectrostaticMap(me2)%dipole_moment
404 #ifdef IS_MPI
405 ul_j(1) = eFrame_Col(3,atom2)
406 ul_j(2) = eFrame_Col(6,atom2)
407 ul_j(3) = eFrame_Col(9,atom2)
408 #else
409 ul_j(1) = eFrame(3,atom2)
410 ul_j(2) = eFrame(6,atom2)
411 ul_j(3) = eFrame(9,atom2)
412 #endif
413 ct_j = ul_j(1)*drdxj + ul_j(2)*drdyj + ul_j(3)*drdzj
414
415 if (j_is_SplitDipole) then
416 d_j = ElectrostaticMap(me2)%split_dipole_distance
417 endif
418 endif
419
420 epot = 0.0_dp
421 dudx = 0.0_dp
422 dudy = 0.0_dp
423 dudz = 0.0_dp
424
425 duduix = 0.0_dp
426 duduiy = 0.0_dp
427 duduiz = 0.0_dp
428
429 dudujx = 0.0_dp
430 dudujy = 0.0_dp
431 dudujz = 0.0_dp
432
433 if (i_is_Charge) then
434
435 if (j_is_Charge) then
436
437 vterm = pre11 * q_i * q_j * riji
438 vpair = vpair + vterm
439 epot = epot + sw*vterm
440
441 dudr = - sw * vterm * riji
442
443 dudx = dudx + dudr * drdxj
444 dudy = dudy + dudr * drdyj
445 dudz = dudz + dudr * drdzj
446
447 endif
448
449 if (j_is_Dipole) then
450
451 if (j_is_SplitDipole) then
452 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
453 ri = 1.0_dp / BigR
454 scale = rij * ri
455 else
456 ri = riji
457 scale = 1.0_dp
458 endif
459
460 ri2 = ri * ri
461 ri3 = ri2 * ri
462 sc2 = scale * scale
463
464 pref = pre12 * q_i * mu_j
465 vterm = pref * ct_j * ri2 * scale
466 vpair = vpair + vterm
467 epot = epot + sw * vterm
468
469 !! this has a + sign in the () because the rij vector is
470 !! r_j - r_i and the charge-dipole potential takes the origin
471 !! as the point dipole, which is atom j in this case.
472
473 dudx = dudx + pref * sw * ri3 * ( ul_j(1) + 3.0d0*ct_j*xhat*sc2)
474 dudy = dudy + pref * sw * ri3 * ( ul_j(2) + 3.0d0*ct_j*yhat*sc2)
475 dudz = dudz + pref * sw * ri3 * ( ul_j(3) + 3.0d0*ct_j*zhat*sc2)
476
477 dudujx = dudujx - pref * sw * ri2 * xhat * scale
478 dudujy = dudujy - pref * sw * ri2 * yhat * scale
479 dudujz = dudujz - pref * sw * ri2 * zhat * scale
480
481 endif
482
483 endif
484
485 if (i_is_Dipole) then
486
487 if (j_is_Charge) then
488
489 if (i_is_SplitDipole) then
490 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
491 ri = 1.0_dp / BigR
492 scale = rij * ri
493 else
494 ri = riji
495 scale = 1.0_dp
496 endif
497
498 ri2 = ri * ri
499 ri3 = ri2 * ri
500 sc2 = scale * scale
501
502 pref = pre12 * q_j * mu_i
503 vterm = pref * ct_i * ri2 * scale
504 vpair = vpair + vterm
505 epot = epot + sw * vterm
506
507 dudx = dudx + pref * sw * ri3 * ( ul_i(1) - 3.0d0 * ct_i * xhat*sc2)
508 dudy = dudy + pref * sw * ri3 * ( ul_i(2) - 3.0d0 * ct_i * yhat*sc2)
509 dudz = dudz + pref * sw * ri3 * ( ul_i(3) - 3.0d0 * ct_i * zhat*sc2)
510
511 duduix = duduix + pref * sw * ri2 * xhat * scale
512 duduiy = duduiy + pref * sw * ri2 * yhat * scale
513 duduiz = duduiz + pref * sw * ri2 * zhat * scale
514 endif
515
516 if (j_is_Dipole) then
517
518 if (i_is_SplitDipole) then
519 if (j_is_SplitDipole) then
520 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
521 else
522 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
523 endif
524 ri = 1.0_dp / BigR
525 scale = rij * ri
526 else
527 if (j_is_SplitDipole) then
528 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
529 ri = 1.0_dp / BigR
530 scale = rij * ri
531 else
532 ri = riji
533 scale = 1.0_dp
534 endif
535 endif
536
537 ct_ij = ul_i(1)*ul_j(1) + ul_i(2)*ul_j(2) + ul_i(3)*ul_j(3)
538
539 ri2 = ri * ri
540 ri3 = ri2 * ri
541 ri4 = ri2 * ri2
542 sc2 = scale * scale
543
544 pref = pre22 * mu_i * mu_j
545 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
546 vpair = vpair + vterm
547 epot = epot + sw * vterm
548
549 a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
550
551 dudx=dudx+pref*sw*3.0d0*ri4*scale*(a1*xhat-ct_i*ul_j(1)-ct_j*ul_i(1))
552 dudy=dudy+pref*sw*3.0d0*ri4*scale*(a1*yhat-ct_i*ul_j(2)-ct_j*ul_i(2))
553 dudz=dudz+pref*sw*3.0d0*ri4*scale*(a1*zhat-ct_i*ul_j(3)-ct_j*ul_i(3))
554
555 duduix = duduix + pref*sw*ri3*(ul_j(1) - 3.0d0*ct_j*xhat*sc2)
556 duduiy = duduiy + pref*sw*ri3*(ul_j(2) - 3.0d0*ct_j*yhat*sc2)
557 duduiz = duduiz + pref*sw*ri3*(ul_j(3) - 3.0d0*ct_j*zhat*sc2)
558
559 dudujx = dudujx + pref*sw*ri3*(ul_i(1) - 3.0d0*ct_i*xhat*sc2)
560 dudujy = dudujy + pref*sw*ri3*(ul_i(2) - 3.0d0*ct_i*yhat*sc2)
561 dudujz = dudujz + pref*sw*ri3*(ul_i(3) - 3.0d0*ct_i*zhat*sc2)
562 endif
563
564 endif
565
566 if (do_pot) then
567 #ifdef IS_MPI
568 pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
569 pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
570 #else
571 pot = pot + epot
572 #endif
573 endif
574
575 #ifdef IS_MPI
576 f_Row(1,atom1) = f_Row(1,atom1) + dudx
577 f_Row(2,atom1) = f_Row(2,atom1) + dudy
578 f_Row(3,atom1) = f_Row(3,atom1) + dudz
579
580 f_Col(1,atom2) = f_Col(1,atom2) - dudx
581 f_Col(2,atom2) = f_Col(2,atom2) - dudy
582 f_Col(3,atom2) = f_Col(3,atom2) - dudz
583
584 if (i_is_Dipole .or. i_is_Quadrupole) then
585 t_Row(1,atom1) = t_Row(1,atom1) - ul_i(2)*duduiz + ul_i(3)*duduiy
586 t_Row(2,atom1) = t_Row(2,atom1) - ul_i(3)*duduix + ul_i(1)*duduiz
587 t_Row(3,atom1) = t_Row(3,atom1) - ul_i(1)*duduiy + ul_i(2)*duduix
588 endif
589
590 if (j_is_Dipole .or. j_is_Quadrupole) then
591 t_Col(1,atom2) = t_Col(1,atom2) - ul_j(2)*dudujz + ul_j(3)*dudujy
592 t_Col(2,atom2) = t_Col(2,atom2) - ul_j(3)*dudujx + ul_j(1)*dudujz
593 t_Col(3,atom2) = t_Col(3,atom2) - ul_j(1)*dudujy + ul_j(2)*dudujx
594 endif
595
596 #else
597 f(1,atom1) = f(1,atom1) + dudx
598 f(2,atom1) = f(2,atom1) + dudy
599 f(3,atom1) = f(3,atom1) + dudz
600
601 f(1,atom2) = f(1,atom2) - dudx
602 f(2,atom2) = f(2,atom2) - dudy
603 f(3,atom2) = f(3,atom2) - dudz
604
605 if (i_is_Dipole .or. i_is_Quadrupole) then
606 t(1,atom1) = t(1,atom1) - ul_i(2)*duduiz + ul_i(3)*duduiy
607 t(2,atom1) = t(2,atom1) - ul_i(3)*duduix + ul_i(1)*duduiz
608 t(3,atom1) = t(3,atom1) - ul_i(1)*duduiy + ul_i(2)*duduix
609 endif
610
611 if (j_is_Dipole .or. j_is_Quadrupole) then
612 t(1,atom2) = t(1,atom2) - ul_j(2)*dudujz + ul_j(3)*dudujy
613 t(2,atom2) = t(2,atom2) - ul_j(3)*dudujx + ul_j(1)*dudujz
614 t(3,atom2) = t(3,atom2) - ul_j(1)*dudujy + ul_j(2)*dudujx
615 endif
616 #endif
617
618 #ifdef IS_MPI
619 id1 = AtomRowToGlobal(atom1)
620 id2 = AtomColToGlobal(atom2)
621 #else
622 id1 = atom1
623 id2 = atom2
624 #endif
625
626 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
627
628 fpair(1) = fpair(1) + dudx
629 fpair(2) = fpair(2) + dudy
630 fpair(3) = fpair(3) + dudz
631
632 endif
633
634 return
635 end subroutine doElectrostaticPair
636
637 end module electrostatic_module